xref: /petsc/src/dm/impls/plex/plexfem.c (revision 1c41a8ca76d8605dc7d579627879ffaa2f37c0cf)
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 
33955f2e967SMatthew G. Knepley #undef __FUNCT__
34055f2e967SMatthew G. Knepley #define __FUNCT__ "DMPlexInsertBoundaryValues"
34155f2e967SMatthew G. Knepley PetscErrorCode DMPlexInsertBoundaryValues(DM dm, Vec localX)
34255f2e967SMatthew G. Knepley {
34355f2e967SMatthew G. Knepley   void        (**funcs)(const PetscReal x[], PetscScalar *u, void *ctx);
34455f2e967SMatthew G. Knepley   void         **ctxs;
34555f2e967SMatthew G. Knepley   PetscFE       *fe;
34655f2e967SMatthew G. Knepley   PetscInt       numFields, f, numBd, b;
34755f2e967SMatthew G. Knepley   PetscErrorCode ierr;
34855f2e967SMatthew G. Knepley 
34955f2e967SMatthew G. Knepley   PetscFunctionBegin;
35055f2e967SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
35155f2e967SMatthew G. Knepley   PetscValidHeaderSpecific(localX, VEC_CLASSID, 2);
35255f2e967SMatthew G. Knepley   ierr = DMGetNumFields(dm, &numFields);CHKERRQ(ierr);
35355f2e967SMatthew G. Knepley   ierr = PetscMalloc3(numFields,&fe,numFields,&funcs,numFields,&ctxs);CHKERRQ(ierr);
35455f2e967SMatthew G. Knepley   for (f = 0; f < numFields; ++f) {ierr = DMGetField(dm, f, (PetscObject *) &fe[f]);CHKERRQ(ierr);}
35555f2e967SMatthew G. Knepley   /* OPT: Could attempt to do multiple BCs at once */
35655f2e967SMatthew G. Knepley   ierr = DMPlexGetNumBoundary(dm, &numBd);CHKERRQ(ierr);
35755f2e967SMatthew G. Knepley   for (b = 0; b < numBd; ++b) {
358a18a7fb9SMatthew G. Knepley     DMLabel         label;
35955f2e967SMatthew G. Knepley     const PetscInt *ids;
360a18a7fb9SMatthew G. Knepley     const char     *name;
36155f2e967SMatthew G. Knepley     PetscInt        numids, field;
36255f2e967SMatthew G. Knepley     PetscBool       isEssential;
36355f2e967SMatthew G. Knepley     void          (*func)();
36455f2e967SMatthew G. Knepley     void           *ctx;
36555f2e967SMatthew G. Knepley 
36655f2e967SMatthew G. Knepley     /* TODO: We need to set only the part indicated by the ids */
367a18a7fb9SMatthew G. Knepley     ierr = DMPlexGetBoundary(dm, b, &isEssential, &name, &field, &func, &numids, &ids, &ctx);CHKERRQ(ierr);
368a18a7fb9SMatthew G. Knepley     ierr = DMPlexGetLabel(dm, name, &label);CHKERRQ(ierr);
36955f2e967SMatthew G. Knepley     for (f = 0; f < numFields; ++f) {
37055f2e967SMatthew G. Knepley       funcs[f] = field == f ? (void (*)(const PetscReal[], PetscScalar *, void *)) func : NULL;
37155f2e967SMatthew G. Knepley       ctxs[f]  = field == f ? ctx : NULL;
37255f2e967SMatthew G. Knepley     }
373a18a7fb9SMatthew G. Knepley     ierr = DMPlexProjectFunctionLabelLocal(dm, label, numids, ids, fe, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr);
37455f2e967SMatthew G. Knepley   }
37555f2e967SMatthew G. Knepley   ierr = PetscFree3(fe,funcs,ctxs);CHKERRQ(ierr);
37655f2e967SMatthew G. Knepley   PetscFunctionReturn(0);
37755f2e967SMatthew G. Knepley }
37855f2e967SMatthew G. Knepley 
379cb1e1211SMatthew G Knepley #undef __FUNCT__
380cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexComputeL2Diff"
381cb1e1211SMatthew G Knepley /*@C
382cb1e1211SMatthew G Knepley   DMPlexComputeL2Diff - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h.
383cb1e1211SMatthew G Knepley 
384cb1e1211SMatthew G Knepley   Input Parameters:
385cb1e1211SMatthew G Knepley + dm    - The DM
386c5bbbd5bSMatthew G. Knepley . fe    - The PetscFE object for each field
387cb1e1211SMatthew G Knepley . funcs - The functions to evaluate for each field component
38851259fa3SMatthew G. Knepley . ctxs  - Optional array of contexts to pass to each function, or NULL.
389cb1e1211SMatthew G Knepley - X     - The coefficient vector u_h
390cb1e1211SMatthew G Knepley 
391cb1e1211SMatthew G Knepley   Output Parameter:
392cb1e1211SMatthew G Knepley . diff - The diff ||u - u_h||_2
393cb1e1211SMatthew G Knepley 
394cb1e1211SMatthew G Knepley   Level: developer
395cb1e1211SMatthew G Knepley 
39623d86601SMatthew G. Knepley .seealso: DMPlexProjectFunction(), DMPlexComputeL2GradientDiff()
397878cb397SSatish Balay @*/
398c110b1eeSGeoffrey Irving PetscErrorCode DMPlexComputeL2Diff(DM dm, PetscFE fe[], void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, PetscReal *diff)
399cb1e1211SMatthew G Knepley {
400cb1e1211SMatthew G Knepley   const PetscInt  debug = 0;
401cb1e1211SMatthew G Knepley   PetscSection    section;
402c5bbbd5bSMatthew G. Knepley   PetscQuadrature quad;
403cb1e1211SMatthew G Knepley   Vec             localX;
40472f94c41SMatthew G. Knepley   PetscScalar    *funcVal;
405cb1e1211SMatthew G Knepley   PetscReal      *coords, *v0, *J, *invJ, detJ;
406cb1e1211SMatthew G Knepley   PetscReal       localDiff = 0.0;
407cb1e1211SMatthew G Knepley   PetscInt        dim, numFields, numComponents = 0, cStart, cEnd, c, field, fieldOffset, comp;
408cb1e1211SMatthew G Knepley   PetscErrorCode  ierr;
409cb1e1211SMatthew G Knepley 
410cb1e1211SMatthew G Knepley   PetscFunctionBegin;
411cb1e1211SMatthew G Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
412cb1e1211SMatthew G Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
413cb1e1211SMatthew G Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
414cb1e1211SMatthew G Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
415cb1e1211SMatthew G Knepley   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
416cb1e1211SMatthew G Knepley   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
417cb1e1211SMatthew G Knepley   for (field = 0; field < numFields; ++field) {
418c5bbbd5bSMatthew G. Knepley     PetscInt Nc;
419c5bbbd5bSMatthew G. Knepley 
420c5bbbd5bSMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe[field], &Nc);CHKERRQ(ierr);
421c5bbbd5bSMatthew G. Knepley     numComponents += Nc;
422cb1e1211SMatthew G Knepley   }
423c110b1eeSGeoffrey Irving   ierr = DMPlexProjectFunctionLocal(dm, fe, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr);
424dcca6d9dSJed Brown   ierr = PetscMalloc5(numComponents,&funcVal,dim,&coords,dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr);
425cb1e1211SMatthew G Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
426c5bbbd5bSMatthew G. Knepley   ierr = PetscFEGetQuadrature(fe[0], &quad);CHKERRQ(ierr);
427cb1e1211SMatthew G Knepley   for (c = cStart; c < cEnd; ++c) {
428a1e44745SMatthew G. Knepley     PetscScalar *x = NULL;
429cb1e1211SMatthew G Knepley     PetscReal    elemDiff = 0.0;
430cb1e1211SMatthew G Knepley 
431cb1e1211SMatthew G Knepley     ierr = DMPlexComputeCellGeometry(dm, c, v0, J, invJ, &detJ);CHKERRQ(ierr);
432cb1e1211SMatthew G Knepley     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
433cb1e1211SMatthew G Knepley     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
434cb1e1211SMatthew G Knepley 
435cb1e1211SMatthew G Knepley     for (field = 0, comp = 0, fieldOffset = 0; field < numFields; ++field) {
436c110b1eeSGeoffrey Irving       void * const     ctx = ctxs ? ctxs[field] : NULL;
43721454ff5SMatthew G. Knepley       const PetscReal *quadPoints, *quadWeights;
438c5bbbd5bSMatthew G. Knepley       PetscReal       *basis;
43921454ff5SMatthew G. Knepley       PetscInt         numQuadPoints, numBasisFuncs, numBasisComps, q, d, e, fc, f;
440cb1e1211SMatthew G Knepley 
44121454ff5SMatthew G. Knepley       ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr);
442c5bbbd5bSMatthew G. Knepley       ierr = PetscFEGetDimension(fe[field], &numBasisFuncs);CHKERRQ(ierr);
443c5bbbd5bSMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe[field], &numBasisComps);CHKERRQ(ierr);
444c5bbbd5bSMatthew G. Knepley       ierr = PetscFEGetDefaultTabulation(fe[field], &basis, NULL, NULL);CHKERRQ(ierr);
445cb1e1211SMatthew G Knepley       if (debug) {
446cb1e1211SMatthew G Knepley         char title[1024];
447cb1e1211SMatthew G Knepley         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
448cb1e1211SMatthew G Knepley         ierr = DMPrintCellVector(c, title, numBasisFuncs*numBasisComps, &x[fieldOffset]);CHKERRQ(ierr);
449cb1e1211SMatthew G Knepley       }
450cb1e1211SMatthew G Knepley       for (q = 0; q < numQuadPoints; ++q) {
451cb1e1211SMatthew G Knepley         for (d = 0; d < dim; d++) {
452cb1e1211SMatthew G Knepley           coords[d] = v0[d];
453cb1e1211SMatthew G Knepley           for (e = 0; e < dim; e++) {
454cb1e1211SMatthew G Knepley             coords[d] += J[d*dim+e]*(quadPoints[q*dim+e] + 1.0);
455cb1e1211SMatthew G Knepley           }
456cb1e1211SMatthew G Knepley         }
457c110b1eeSGeoffrey Irving         (*funcs[field])(coords, funcVal, ctx);
458cb1e1211SMatthew G Knepley         for (fc = 0; fc < numBasisComps; ++fc) {
459a1d24da5SMatthew G. Knepley           PetscScalar interpolant = 0.0;
460a1d24da5SMatthew G. Knepley 
461cb1e1211SMatthew G Knepley           for (f = 0; f < numBasisFuncs; ++f) {
462cb1e1211SMatthew G Knepley             const PetscInt fidx = f*numBasisComps+fc;
463a1d24da5SMatthew G. Knepley             interpolant += x[fieldOffset+fidx]*basis[q*numBasisFuncs*numBasisComps+fidx];
464cb1e1211SMatthew G Knepley           }
46572f94c41SMatthew 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);}
46672f94c41SMatthew G. Knepley           elemDiff += PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ;
467cb1e1211SMatthew G Knepley         }
468cb1e1211SMatthew G Knepley       }
469cb1e1211SMatthew G Knepley       comp        += numBasisComps;
470cb1e1211SMatthew G Knepley       fieldOffset += numBasisFuncs*numBasisComps;
471cb1e1211SMatthew G Knepley     }
472cb1e1211SMatthew G Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
473cb1e1211SMatthew G Knepley     if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);}
474cb1e1211SMatthew G Knepley     localDiff += elemDiff;
475cb1e1211SMatthew G Knepley   }
47672f94c41SMatthew G. Knepley   ierr  = PetscFree5(funcVal,coords,v0,J,invJ);CHKERRQ(ierr);
477cb1e1211SMatthew G Knepley   ierr  = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
47886a74ee0SMatthew G. Knepley   ierr  = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
479cb1e1211SMatthew G Knepley   *diff = PetscSqrtReal(*diff);
480cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
481cb1e1211SMatthew G Knepley }
482cb1e1211SMatthew G Knepley 
483cb1e1211SMatthew G Knepley #undef __FUNCT__
48440e14135SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeL2GradientDiff"
48540e14135SMatthew G. Knepley /*@C
48640e14135SMatthew 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.
48740e14135SMatthew G. Knepley 
48840e14135SMatthew G. Knepley   Input Parameters:
48940e14135SMatthew G. Knepley + dm    - The DM
49040e14135SMatthew G. Knepley . fe    - The PetscFE object for each field
49140e14135SMatthew G. Knepley . funcs - The gradient functions to evaluate for each field component
49251259fa3SMatthew G. Knepley . ctxs  - Optional array of contexts to pass to each function, or NULL.
49340e14135SMatthew G. Knepley . X     - The coefficient vector u_h
49440e14135SMatthew G. Knepley - n     - The vector to project along
49540e14135SMatthew G. Knepley 
49640e14135SMatthew G. Knepley   Output Parameter:
49740e14135SMatthew G. Knepley . diff - The diff ||(grad u - grad u_h) . n||_2
49840e14135SMatthew G. Knepley 
49940e14135SMatthew G. Knepley   Level: developer
50040e14135SMatthew G. Knepley 
50140e14135SMatthew G. Knepley .seealso: DMPlexProjectFunction(), DMPlexComputeL2Diff()
50240e14135SMatthew G. Knepley @*/
50351259fa3SMatthew 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)
504cb1e1211SMatthew G Knepley {
50540e14135SMatthew G. Knepley   const PetscInt  debug = 0;
506cb1e1211SMatthew G Knepley   PetscSection    section;
50740e14135SMatthew G. Knepley   PetscQuadrature quad;
50840e14135SMatthew G. Knepley   Vec             localX;
50940e14135SMatthew G. Knepley   PetscScalar    *funcVal, *interpolantVec;
51040e14135SMatthew G. Knepley   PetscReal      *coords, *realSpaceDer, *v0, *J, *invJ, detJ;
51140e14135SMatthew G. Knepley   PetscReal       localDiff = 0.0;
51240e14135SMatthew G. Knepley   PetscInt        dim, numFields, numComponents = 0, cStart, cEnd, c, field, fieldOffset, comp;
513cb1e1211SMatthew G Knepley   PetscErrorCode  ierr;
514cb1e1211SMatthew G Knepley 
515cb1e1211SMatthew G Knepley   PetscFunctionBegin;
51640e14135SMatthew G. Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
51740e14135SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
51840e14135SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
51940e14135SMatthew G. Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
52040e14135SMatthew G. Knepley   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
52140e14135SMatthew G. Knepley   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
522652b88e8SMatthew G. Knepley   for (field = 0; field < numFields; ++field) {
52340e14135SMatthew G. Knepley     PetscInt Nc;
524652b88e8SMatthew G. Knepley 
52540e14135SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe[field], &Nc);CHKERRQ(ierr);
52640e14135SMatthew G. Knepley     numComponents += Nc;
527652b88e8SMatthew G. Knepley   }
52840e14135SMatthew G. Knepley   /* ierr = DMPlexProjectFunctionLocal(dm, fe, funcs, INSERT_BC_VALUES, localX);CHKERRQ(ierr); */
52940e14135SMatthew G. Knepley   ierr = PetscMalloc7(numComponents,&funcVal,dim,&coords,dim,&realSpaceDer,dim,&v0,dim*dim,&J,dim*dim,&invJ,dim,&interpolantVec);CHKERRQ(ierr);
53040e14135SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
53140e14135SMatthew G. Knepley   ierr = PetscFEGetQuadrature(fe[0], &quad);CHKERRQ(ierr);
53240e14135SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
53340e14135SMatthew G. Knepley     PetscScalar *x = NULL;
53440e14135SMatthew G. Knepley     PetscReal    elemDiff = 0.0;
535652b88e8SMatthew G. Knepley 
53640e14135SMatthew G. Knepley     ierr = DMPlexComputeCellGeometry(dm, c, v0, J, invJ, &detJ);CHKERRQ(ierr);
53740e14135SMatthew G. Knepley     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
53840e14135SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
53940e14135SMatthew G. Knepley 
54040e14135SMatthew G. Knepley     for (field = 0, comp = 0, fieldOffset = 0; field < numFields; ++field) {
54151259fa3SMatthew G. Knepley       void * const     ctx = ctxs ? ctxs[field] : NULL;
54221454ff5SMatthew G. Knepley       const PetscReal *quadPoints, *quadWeights;
54340e14135SMatthew G. Knepley       PetscReal       *basisDer;
54421454ff5SMatthew G. Knepley       PetscInt         numQuadPoints, Nb, Ncomp, q, d, e, fc, f, g;
54540e14135SMatthew G. Knepley 
54621454ff5SMatthew G. Knepley       ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr);
54740e14135SMatthew G. Knepley       ierr = PetscFEGetDimension(fe[field], &Nb);CHKERRQ(ierr);
54840e14135SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe[field], &Ncomp);CHKERRQ(ierr);
54940e14135SMatthew G. Knepley       ierr = PetscFEGetDefaultTabulation(fe[field], NULL, &basisDer, NULL);CHKERRQ(ierr);
55040e14135SMatthew G. Knepley       if (debug) {
55140e14135SMatthew G. Knepley         char title[1024];
55240e14135SMatthew G. Knepley         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
55340e14135SMatthew G. Knepley         ierr = DMPrintCellVector(c, title, Nb*Ncomp, &x[fieldOffset]);CHKERRQ(ierr);
554652b88e8SMatthew G. Knepley       }
55540e14135SMatthew G. Knepley       for (q = 0; q < numQuadPoints; ++q) {
55640e14135SMatthew G. Knepley         for (d = 0; d < dim; d++) {
55740e14135SMatthew G. Knepley           coords[d] = v0[d];
55840e14135SMatthew G. Knepley           for (e = 0; e < dim; e++) {
55940e14135SMatthew G. Knepley             coords[d] += J[d*dim+e]*(quadPoints[q*dim+e] + 1.0);
560652b88e8SMatthew G. Knepley           }
56140e14135SMatthew G. Knepley         }
56251259fa3SMatthew G. Knepley         (*funcs[field])(coords, n, funcVal, ctx);
56340e14135SMatthew G. Knepley         for (fc = 0; fc < Ncomp; ++fc) {
56440e14135SMatthew G. Knepley           PetscScalar interpolant = 0.0;
56540e14135SMatthew G. Knepley 
56640e14135SMatthew G. Knepley           for (d = 0; d < dim; ++d) interpolantVec[d] = 0.0;
56740e14135SMatthew G. Knepley           for (f = 0; f < Nb; ++f) {
56840e14135SMatthew G. Knepley             const PetscInt fidx = f*Ncomp+fc;
56940e14135SMatthew G. Knepley 
57040e14135SMatthew G. Knepley             for (d = 0; d < dim; ++d) {
57140e14135SMatthew G. Knepley               realSpaceDer[d] = 0.0;
57240e14135SMatthew G. Knepley               for (g = 0; g < dim; ++g) {
57340e14135SMatthew G. Knepley                 realSpaceDer[d] += invJ[g*dim+d]*basisDer[(q*Nb*Ncomp+fidx)*dim+g];
57440e14135SMatthew G. Knepley               }
57540e14135SMatthew G. Knepley               interpolantVec[d] += x[fieldOffset+fidx]*realSpaceDer[d];
57640e14135SMatthew G. Knepley             }
57740e14135SMatthew G. Knepley           }
57840e14135SMatthew G. Knepley           for (d = 0; d < dim; ++d) interpolant += interpolantVec[d]*n[d];
57940e14135SMatthew 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);}
58040e14135SMatthew G. Knepley           elemDiff += PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ;
58140e14135SMatthew G. Knepley         }
58240e14135SMatthew G. Knepley       }
58340e14135SMatthew G. Knepley       comp        += Ncomp;
58440e14135SMatthew G. Knepley       fieldOffset += Nb*Ncomp;
58540e14135SMatthew G. Knepley     }
58640e14135SMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
58740e14135SMatthew G. Knepley     if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);}
58840e14135SMatthew G. Knepley     localDiff += elemDiff;
58940e14135SMatthew G. Knepley   }
59040e14135SMatthew G. Knepley   ierr  = PetscFree7(funcVal,coords,realSpaceDer,v0,J,invJ,interpolantVec);CHKERRQ(ierr);
59140e14135SMatthew G. Knepley   ierr  = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
59240e14135SMatthew G. Knepley   ierr  = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
59340e14135SMatthew G. Knepley   *diff = PetscSqrtReal(*diff);
594cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
595cb1e1211SMatthew G Knepley }
596cb1e1211SMatthew G Knepley 
597a0845e3aSMatthew G. Knepley #undef __FUNCT__
598a0845e3aSMatthew G. Knepley #define __FUNCT__ "DMPlexComputeResidualFEM"
599a0845e3aSMatthew G. Knepley /*@
600a0845e3aSMatthew G. Knepley   DMPlexComputeResidualFEM - Form the local residual F from the local input X using pointwise functions specified by the user
601a0845e3aSMatthew G. Knepley 
602a0845e3aSMatthew G. Knepley   Input Parameters:
603a0845e3aSMatthew G. Knepley + dm - The mesh
604a0845e3aSMatthew G. Knepley . X  - Local input vector
605a0845e3aSMatthew G. Knepley - user - The user context
606a0845e3aSMatthew G. Knepley 
607a0845e3aSMatthew G. Knepley   Output Parameter:
608a0845e3aSMatthew G. Knepley . F  - Local output vector
609a0845e3aSMatthew G. Knepley 
610a0845e3aSMatthew G. Knepley   Note:
6118026896cSMatthew G. Knepley   The first member of the user context must be an FEMContext.
612a0845e3aSMatthew G. Knepley 
613a0845e3aSMatthew G. Knepley   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
614a0845e3aSMatthew G. Knepley   like a GPU, or vectorize on a multicore machine.
615a0845e3aSMatthew G. Knepley 
616a0845e3aSMatthew G. Knepley   Level: developer
617a0845e3aSMatthew G. Knepley 
618a0845e3aSMatthew G. Knepley .seealso: DMPlexComputeJacobianActionFEM()
619a0845e3aSMatthew G. Knepley @*/
620a0845e3aSMatthew G. Knepley PetscErrorCode DMPlexComputeResidualFEM(DM dm, Vec X, Vec F, void *user)
621a0845e3aSMatthew G. Knepley {
622a0845e3aSMatthew G. Knepley   DM_Plex          *mesh  = (DM_Plex *) dm->data;
6239a559087SMatthew G. Knepley   PetscFEM         *fem   = (PetscFEM *) user;
624a0845e3aSMatthew G. Knepley   PetscFE          *fe    = fem->fe;
6259a559087SMatthew G. Knepley   PetscFE          *feAux = fem->feAux;
626f1ea0e2fSMatthew G. Knepley   PetscFE          *feBd  = fem->feBd;
627a0845e3aSMatthew G. Knepley   const char       *name  = "Residual";
6289a559087SMatthew G. Knepley   DM                dmAux;
6299a559087SMatthew G. Knepley   Vec               A;
630a0845e3aSMatthew G. Knepley   PetscQuadrature   q;
631a0845e3aSMatthew G. Knepley   PetscCellGeometry geom;
6329a559087SMatthew G. Knepley   PetscSection      section, sectionAux;
633a0845e3aSMatthew G. Knepley   PetscReal        *v0, *J, *invJ, *detJ;
63401599b20SMatthew G. Knepley   PetscScalar      *elemVec, *u, *a = NULL;
6359a559087SMatthew G. Knepley   PetscInt          dim, Nf, NfAux = 0, f, numCells, cStart, cEnd, c;
636a0845e3aSMatthew G. Knepley   PetscInt          cellDof = 0, numComponents = 0;
6379a559087SMatthew G. Knepley   PetscInt          cellDofAux = 0, numComponentsAux = 0;
638a0845e3aSMatthew G. Knepley   PetscErrorCode    ierr;
639a0845e3aSMatthew G. Knepley 
640a0845e3aSMatthew G. Knepley   PetscFunctionBegin;
641a0845e3aSMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_ResidualFEM,dm,0,0,0);CHKERRQ(ierr);
642a0845e3aSMatthew G. Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
643a0845e3aSMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
6449a559087SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
645a0845e3aSMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
646a0845e3aSMatthew G. Knepley   numCells = cEnd - cStart;
6479a559087SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
648a0845e3aSMatthew G. Knepley     PetscInt Nb, Nc;
649a0845e3aSMatthew G. Knepley 
650a0845e3aSMatthew G. Knepley     ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr);
651a0845e3aSMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe[f], &Nc);CHKERRQ(ierr);
652a0845e3aSMatthew G. Knepley     cellDof       += Nb*Nc;
653a0845e3aSMatthew G. Knepley     numComponents += Nc;
654a0845e3aSMatthew G. Knepley   }
6559a559087SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
6569a559087SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr);
6579a559087SMatthew G. Knepley   if (dmAux) {
6589a559087SMatthew G. Knepley     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
6599a559087SMatthew G. Knepley     ierr = PetscSectionGetNumFields(sectionAux, &NfAux);CHKERRQ(ierr);
6609a559087SMatthew G. Knepley   }
6619a559087SMatthew G. Knepley   for (f = 0; f < NfAux; ++f) {
6629a559087SMatthew G. Knepley     PetscInt Nb, Nc;
6639a559087SMatthew G. Knepley 
6649a559087SMatthew G. Knepley     ierr = PetscFEGetDimension(feAux[f], &Nb);CHKERRQ(ierr);
6659a559087SMatthew G. Knepley     ierr = PetscFEGetNumComponents(feAux[f], &Nc);CHKERRQ(ierr);
6669a559087SMatthew G. Knepley     cellDofAux       += Nb*Nc;
6679a559087SMatthew G. Knepley     numComponentsAux += Nc;
6689a559087SMatthew G. Knepley   }
669bd5ce809SMatthew G. Knepley   ierr = DMPlexInsertBoundaryValues(dm, X);CHKERRQ(ierr);
670a0845e3aSMatthew G. Knepley   ierr = VecSet(F, 0.0);CHKERRQ(ierr);
671dcca6d9dSJed Brown   ierr = PetscMalloc6(numCells*cellDof,&u,numCells*dim,&v0,numCells*dim*dim,&J,numCells*dim*dim,&invJ,numCells,&detJ,numCells*cellDof,&elemVec);CHKERRQ(ierr);
672785e854fSJed Brown   if (dmAux) {ierr = PetscMalloc1(numCells*cellDofAux, &a);CHKERRQ(ierr);}
673a0845e3aSMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
674a0845e3aSMatthew G. Knepley     PetscScalar *x = NULL;
675a0845e3aSMatthew G. Knepley     PetscInt     i;
676a0845e3aSMatthew G. Knepley 
677a0845e3aSMatthew G. Knepley     ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr);
678a0845e3aSMatthew G. Knepley     if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c);
679a0845e3aSMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr);
680a0845e3aSMatthew G. Knepley     for (i = 0; i < cellDof; ++i) u[c*cellDof+i] = x[i];
681a0845e3aSMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr);
6829a559087SMatthew G. Knepley     if (dmAux) {
6839a559087SMatthew G. Knepley       ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
6849a559087SMatthew G. Knepley       for (i = 0; i < cellDofAux; ++i) a[c*cellDofAux+i] = x[i];
6859a559087SMatthew G. Knepley       ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
686a0845e3aSMatthew G. Knepley     }
6879a559087SMatthew G. Knepley   }
6889a559087SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
689c012ea0aSMatthew G. Knepley     void   (*f0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->f0Funcs[f];
690c012ea0aSMatthew G. Knepley     void   (*f1)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->f1Funcs[f];
69121454ff5SMatthew G. Knepley     PetscInt numQuadPoints, Nb;
692a0845e3aSMatthew G. Knepley     /* Conforming batches */
693f30c5766SMatthew G. Knepley     PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
694a0845e3aSMatthew G. Knepley     /* Remainder */
695a0845e3aSMatthew G. Knepley     PetscInt Nr, offset;
696a0845e3aSMatthew G. Knepley 
697a0845e3aSMatthew G. Knepley     ierr = PetscFEGetQuadrature(fe[f], &q);CHKERRQ(ierr);
698a0845e3aSMatthew G. Knepley     ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr);
699f30c5766SMatthew G. Knepley     ierr = PetscFEGetTileSizes(fe[f], NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
70021454ff5SMatthew G. Knepley     ierr = PetscQuadratureGetData(q, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr);
70121454ff5SMatthew G. Knepley     blockSize = Nb*numQuadPoints;
702a0845e3aSMatthew G. Knepley     batchSize = numBlocks * blockSize;
703f30c5766SMatthew G. Knepley     ierr =  PetscFESetTileSizes(fe[f], blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
704a0845e3aSMatthew G. Knepley     numChunks = numCells / (numBatches*batchSize);
705a0845e3aSMatthew G. Knepley     Ne        = numChunks*numBatches*batchSize;
706a0845e3aSMatthew G. Knepley     Nr        = numCells % (numBatches*batchSize);
707a0845e3aSMatthew G. Knepley     offset    = numCells - Nr;
708a0845e3aSMatthew G. Knepley     geom.v0   = v0;
709a0845e3aSMatthew G. Knepley     geom.J    = J;
710a0845e3aSMatthew G. Knepley     geom.invJ = invJ;
711a0845e3aSMatthew G. Knepley     geom.detJ = detJ;
7129a559087SMatthew G. Knepley     ierr = PetscFEIntegrateResidual(fe[f], Ne, Nf, fe, f, geom, u, NfAux, feAux, a, f0, f1, elemVec);CHKERRQ(ierr);
713a0845e3aSMatthew G. Knepley     geom.v0   = &v0[offset*dim];
714a0845e3aSMatthew G. Knepley     geom.J    = &J[offset*dim*dim];
715a0845e3aSMatthew G. Knepley     geom.invJ = &invJ[offset*dim*dim];
716a0845e3aSMatthew G. Knepley     geom.detJ = &detJ[offset];
7179a559087SMatthew 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);
718a0845e3aSMatthew G. Knepley   }
719a0845e3aSMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
720a0845e3aSMatthew G. Knepley     if (mesh->printFEM > 1) {ierr = DMPrintCellVector(c, name, cellDof, &elemVec[c*cellDof]);CHKERRQ(ierr);}
721a0845e3aSMatthew G. Knepley     ierr = DMPlexVecSetClosure(dm, section, F, c, &elemVec[c*cellDof], ADD_VALUES);CHKERRQ(ierr);
722a0845e3aSMatthew G. Knepley   }
723a0845e3aSMatthew G. Knepley   ierr = PetscFree6(u,v0,J,invJ,detJ,elemVec);CHKERRQ(ierr);
7249a559087SMatthew G. Knepley   if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);}
725f1ea0e2fSMatthew G. Knepley   if (feBd) {
726075da914SMatthew G. Knepley     DMLabel         label, depth;
727f1ea0e2fSMatthew G. Knepley     IS              pointIS;
728f1ea0e2fSMatthew G. Knepley     const PetscInt *points;
729075da914SMatthew G. Knepley     PetscInt        dep, numPoints, p, numFaces;
730f1ea0e2fSMatthew G. Knepley     PetscReal      *n;
731f1ea0e2fSMatthew G. Knepley 
732f1ea0e2fSMatthew G. Knepley     ierr = DMPlexGetLabel(dm, "boundary", &label);CHKERRQ(ierr);
733075da914SMatthew G. Knepley     ierr = DMPlexGetDepthLabel(dm, &depth);CHKERRQ(ierr);
734f1ea0e2fSMatthew G. Knepley     ierr = DMLabelGetStratumSize(label, 1, &numPoints);CHKERRQ(ierr);
735f1ea0e2fSMatthew G. Knepley     ierr = DMLabelGetStratumIS(label, 1, &pointIS);CHKERRQ(ierr);
736f1ea0e2fSMatthew G. Knepley     ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
737f1ea0e2fSMatthew G. Knepley     for (f = 0, cellDof = 0, numComponents = 0; f < Nf; ++f) {
738f1ea0e2fSMatthew G. Knepley       PetscInt Nb, Nc;
739f1ea0e2fSMatthew G. Knepley 
740f1ea0e2fSMatthew G. Knepley       ierr = PetscFEGetDimension(feBd[f], &Nb);CHKERRQ(ierr);
741f1ea0e2fSMatthew G. Knepley       ierr = PetscFEGetNumComponents(feBd[f], &Nc);CHKERRQ(ierr);
742f1ea0e2fSMatthew G. Knepley       cellDof       += Nb*Nc;
743f1ea0e2fSMatthew G. Knepley       numComponents += Nc;
744f1ea0e2fSMatthew G. Knepley     }
745075da914SMatthew G. Knepley     for (p = 0, numFaces = 0; p < numPoints; ++p) {
746075da914SMatthew G. Knepley       ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr);
747075da914SMatthew G. Knepley       if (dep == dim-1) ++numFaces;
748075da914SMatthew G. Knepley     }
749dcca6d9dSJed 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);
750075da914SMatthew G. Knepley     for (p = 0, f = 0; p < numPoints; ++p) {
751f1ea0e2fSMatthew G. Knepley       const PetscInt point = points[p];
752f1ea0e2fSMatthew G. Knepley       PetscScalar   *x     = NULL;
753f1ea0e2fSMatthew G. Knepley       PetscInt       i;
754f1ea0e2fSMatthew G. Knepley 
755075da914SMatthew G. Knepley       ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr);
756075da914SMatthew G. Knepley       if (dep != dim-1) continue;
757075da914SMatthew G. Knepley       ierr = DMPlexComputeCellGeometry(dm, point, &v0[f*dim], &J[f*dim*dim], &invJ[f*dim*dim], &detJ[f]);CHKERRQ(ierr);
758a8007bbfSMatthew G. Knepley       ierr = DMPlexComputeCellGeometryFVM(dm, point, NULL, NULL, &n[f*dim]);
759075da914SMatthew G. Knepley       if (detJ[f] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for face %d", detJ[f], point);
760f1ea0e2fSMatthew G. Knepley       ierr = DMPlexVecGetClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr);
761075da914SMatthew G. Knepley       for (i = 0; i < cellDof; ++i) u[f*cellDof+i] = x[i];
762f1ea0e2fSMatthew G. Knepley       ierr = DMPlexVecRestoreClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr);
763075da914SMatthew G. Knepley       ++f;
764f1ea0e2fSMatthew G. Knepley     }
765f1ea0e2fSMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
766f1ea0e2fSMatthew G. Knepley       void   (*f0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->f0BdFuncs[f];
767f1ea0e2fSMatthew G. Knepley       void   (*f1)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->f1BdFuncs[f];
76821454ff5SMatthew G. Knepley       PetscInt numQuadPoints, Nb;
769f1ea0e2fSMatthew G. Knepley       /* Conforming batches */
770f1ea0e2fSMatthew G. Knepley       PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
771f1ea0e2fSMatthew G. Knepley       /* Remainder */
772f1ea0e2fSMatthew G. Knepley       PetscInt Nr, offset;
773f1ea0e2fSMatthew G. Knepley 
774f1ea0e2fSMatthew G. Knepley       ierr = PetscFEGetQuadrature(feBd[f], &q);CHKERRQ(ierr);
775f1ea0e2fSMatthew G. Knepley       ierr = PetscFEGetDimension(feBd[f], &Nb);CHKERRQ(ierr);
776f1ea0e2fSMatthew G. Knepley       ierr = PetscFEGetTileSizes(feBd[f], NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
77721454ff5SMatthew G. Knepley       ierr = PetscQuadratureGetData(q, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr);
77821454ff5SMatthew G. Knepley       blockSize = Nb*numQuadPoints;
779f1ea0e2fSMatthew G. Knepley       batchSize = numBlocks * blockSize;
780f1ea0e2fSMatthew G. Knepley       ierr =  PetscFESetTileSizes(feBd[f], blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
781075da914SMatthew G. Knepley       numChunks = numFaces / (numBatches*batchSize);
782f1ea0e2fSMatthew G. Knepley       Ne        = numChunks*numBatches*batchSize;
783075da914SMatthew G. Knepley       Nr        = numFaces % (numBatches*batchSize);
784075da914SMatthew G. Knepley       offset    = numFaces - Nr;
785f1ea0e2fSMatthew G. Knepley       geom.v0   = v0;
786f1ea0e2fSMatthew G. Knepley       geom.n    = n;
787f1ea0e2fSMatthew G. Knepley       geom.J    = J;
788f1ea0e2fSMatthew G. Knepley       geom.invJ = invJ;
789f1ea0e2fSMatthew G. Knepley       geom.detJ = detJ;
790f1ea0e2fSMatthew G. Knepley       ierr = PetscFEIntegrateBdResidual(feBd[f], Ne, Nf, feBd, f, geom, u, 0, NULL, NULL, f0, f1, elemVec);CHKERRQ(ierr);
791f1ea0e2fSMatthew G. Knepley       geom.v0   = &v0[offset*dim];
792f1ea0e2fSMatthew G. Knepley       geom.n    = &n[offset*dim];
793f1ea0e2fSMatthew G. Knepley       geom.J    = &J[offset*dim*dim];
794f1ea0e2fSMatthew G. Knepley       geom.invJ = &invJ[offset*dim*dim];
795f1ea0e2fSMatthew G. Knepley       geom.detJ = &detJ[offset];
796f1ea0e2fSMatthew G. Knepley       ierr = PetscFEIntegrateBdResidual(feBd[f], Nr, Nf, feBd, f, geom, &u[offset*cellDof], 0, NULL, NULL, f0, f1, &elemVec[offset*cellDof]);CHKERRQ(ierr);
797f1ea0e2fSMatthew G. Knepley     }
798075da914SMatthew G. Knepley     for (p = 0, f = 0; p < numPoints; ++p) {
799f1ea0e2fSMatthew G. Knepley       const PetscInt point = points[p];
800f1ea0e2fSMatthew G. Knepley 
801075da914SMatthew G. Knepley       ierr = DMLabelGetValue(depth, point, &dep);CHKERRQ(ierr);
802075da914SMatthew G. Knepley       if (dep != dim-1) continue;
803075da914SMatthew G. Knepley       if (mesh->printFEM > 1) {ierr = DMPrintCellVector(point, "BdResidual", cellDof, &elemVec[f*cellDof]);CHKERRQ(ierr);}
804075da914SMatthew G. Knepley       ierr = DMPlexVecSetClosure(dm, NULL, F, point, &elemVec[f*cellDof], ADD_VALUES);CHKERRQ(ierr);
805075da914SMatthew G. Knepley       ++f;
806f1ea0e2fSMatthew G. Knepley     }
807f1ea0e2fSMatthew G. Knepley     ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
808f1ea0e2fSMatthew G. Knepley     ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
809f1ea0e2fSMatthew G. Knepley     ierr = PetscFree7(u,v0,n,J,invJ,detJ,elemVec);CHKERRQ(ierr);
810f1ea0e2fSMatthew G. Knepley   }
8116113b454SMatthew G. Knepley   if (mesh->printFEM) {ierr = DMPrintLocalVec(dm, name, mesh->printTol, F);CHKERRQ(ierr);}
812a0845e3aSMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_ResidualFEM,dm,0,0,0);CHKERRQ(ierr);
813a0845e3aSMatthew G. Knepley   PetscFunctionReturn(0);
814a0845e3aSMatthew G. Knepley }
815a0845e3aSMatthew G. Knepley 
816cb1e1211SMatthew G Knepley #undef __FUNCT__
817af1eca97SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeIFunctionFEM"
818af1eca97SMatthew G. Knepley /*@
819af1eca97SMatthew G. Knepley   DMPlexComputeIFunctionFEM - Form the local implicit function F from the local input X, X_t using pointwise functions specified by the user
820af1eca97SMatthew G. Knepley 
821af1eca97SMatthew G. Knepley   Input Parameters:
822af1eca97SMatthew G. Knepley + dm - The mesh
823*1c41a8caSMatthew G. Knepley . time - The current time
824af1eca97SMatthew G. Knepley . X  - Local input vector
825af1eca97SMatthew G. Knepley . X_t  - Time derivative of the local input vector
826af1eca97SMatthew G. Knepley - user - The user context
827af1eca97SMatthew G. Knepley 
828af1eca97SMatthew G. Knepley   Output Parameter:
829af1eca97SMatthew G. Knepley . F  - Local output vector
830af1eca97SMatthew G. Knepley 
831af1eca97SMatthew G. Knepley   Note:
832af1eca97SMatthew G. Knepley   The first member of the user context must be an FEMContext.
833af1eca97SMatthew G. Knepley 
834af1eca97SMatthew G. Knepley   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
835af1eca97SMatthew G. Knepley   like a GPU, or vectorize on a multicore machine.
836af1eca97SMatthew G. Knepley 
837af1eca97SMatthew G. Knepley   Level: developer
838af1eca97SMatthew G. Knepley 
839af1eca97SMatthew G. Knepley .seealso: DMPlexComputeResidualFEM()
840af1eca97SMatthew G. Knepley @*/
841*1c41a8caSMatthew G. Knepley PetscErrorCode DMPlexComputeIFunctionFEM(DM dm, PetscReal time, Vec X, Vec X_t, Vec F, void *user)
842af1eca97SMatthew G. Knepley {
843af1eca97SMatthew G. Knepley   DM_Plex          *mesh  = (DM_Plex *) dm->data;
844af1eca97SMatthew G. Knepley   PetscFEM         *fem   = (PetscFEM *) user;
845af1eca97SMatthew G. Knepley   PetscFE          *fe    = fem->fe;
846af1eca97SMatthew G. Knepley   PetscFE          *feAux = fem->feAux;
847af1eca97SMatthew G. Knepley   PetscFE          *feBd  = fem->feBd;
848af1eca97SMatthew G. Knepley   const char       *name  = "Residual";
849af1eca97SMatthew G. Knepley   DM                dmAux;
850af1eca97SMatthew G. Knepley   Vec               A;
851af1eca97SMatthew G. Knepley   PetscQuadrature   q;
852af1eca97SMatthew G. Knepley   PetscCellGeometry geom;
853af1eca97SMatthew G. Knepley   PetscSection      section, sectionAux;
854af1eca97SMatthew G. Knepley   PetscReal        *v0, *J, *invJ, *detJ;
855af1eca97SMatthew G. Knepley   PetscScalar      *elemVec, *u, *u_t, *a = NULL;
856af1eca97SMatthew G. Knepley   PetscInt          dim, Nf, NfAux = 0, f, numCells, cStart, cEnd, c;
857af1eca97SMatthew G. Knepley   PetscInt          cellDof = 0, numComponents = 0;
858af1eca97SMatthew G. Knepley   PetscInt          cellDofAux = 0, numComponentsAux = 0;
859af1eca97SMatthew G. Knepley   PetscErrorCode    ierr;
860af1eca97SMatthew G. Knepley 
861af1eca97SMatthew G. Knepley   PetscFunctionBegin;
862af1eca97SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_ResidualFEM,dm,0,0,0);CHKERRQ(ierr);
863af1eca97SMatthew G. Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
864af1eca97SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
865af1eca97SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
866af1eca97SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
867af1eca97SMatthew G. Knepley   numCells = cEnd - cStart;
868af1eca97SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
869af1eca97SMatthew G. Knepley     PetscInt Nb, Nc;
870af1eca97SMatthew G. Knepley 
871af1eca97SMatthew G. Knepley     ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr);
872af1eca97SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe[f], &Nc);CHKERRQ(ierr);
873af1eca97SMatthew G. Knepley     cellDof       += Nb*Nc;
874af1eca97SMatthew G. Knepley     numComponents += Nc;
875af1eca97SMatthew G. Knepley   }
876af1eca97SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
877af1eca97SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr);
878af1eca97SMatthew G. Knepley   if (dmAux) {
879af1eca97SMatthew G. Knepley     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
880af1eca97SMatthew G. Knepley     ierr = PetscSectionGetNumFields(sectionAux, &NfAux);CHKERRQ(ierr);
881af1eca97SMatthew G. Knepley   }
882af1eca97SMatthew G. Knepley   for (f = 0; f < NfAux; ++f) {
883af1eca97SMatthew G. Knepley     PetscInt Nb, Nc;
884af1eca97SMatthew G. Knepley 
885af1eca97SMatthew G. Knepley     ierr = PetscFEGetDimension(feAux[f], &Nb);CHKERRQ(ierr);
886af1eca97SMatthew G. Knepley     ierr = PetscFEGetNumComponents(feAux[f], &Nc);CHKERRQ(ierr);
887af1eca97SMatthew G. Knepley     cellDofAux       += Nb*Nc;
888af1eca97SMatthew G. Knepley     numComponentsAux += Nc;
889af1eca97SMatthew G. Knepley   }
890bd5ce809SMatthew G. Knepley   ierr = DMPlexInsertBoundaryValues(dm, X);CHKERRQ(ierr);
891af1eca97SMatthew G. Knepley   ierr = VecSet(F, 0.0);CHKERRQ(ierr);
892af1eca97SMatthew 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);
893af1eca97SMatthew G. Knepley   if (dmAux) {ierr = PetscMalloc1(numCells*cellDofAux, &a);CHKERRQ(ierr);}
894af1eca97SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
895af1eca97SMatthew G. Knepley     PetscScalar *x = NULL, *x_t = NULL;
896af1eca97SMatthew G. Knepley     PetscInt     i;
897af1eca97SMatthew G. Knepley 
898af1eca97SMatthew G. Knepley     ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr);
899af1eca97SMatthew G. Knepley     if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c);
900af1eca97SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr);
901af1eca97SMatthew G. Knepley     for (i = 0; i < cellDof; ++i) u[c*cellDof+i] = x[i];
902af1eca97SMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr);
903af1eca97SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, section, X_t, c, NULL, &x_t);CHKERRQ(ierr);
904af1eca97SMatthew G. Knepley     for (i = 0; i < cellDof; ++i) u_t[c*cellDof+i] = x_t[i];
905af1eca97SMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, section, X_t, c, NULL, &x_t);CHKERRQ(ierr);
906af1eca97SMatthew G. Knepley     if (dmAux) {
907af1eca97SMatthew G. Knepley       PetscScalar *x_a = NULL;
908af1eca97SMatthew G. Knepley       ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x_a);CHKERRQ(ierr);
909af1eca97SMatthew G. Knepley       for (i = 0; i < cellDofAux; ++i) a[c*cellDofAux+i] = x_a[i];
910af1eca97SMatthew G. Knepley       ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x_a);CHKERRQ(ierr);
911af1eca97SMatthew G. Knepley     }
912af1eca97SMatthew G. Knepley   }
913af1eca97SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
914af1eca97SMatthew G. Knepley     void   (*f0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->f0IFuncs[f];
915af1eca97SMatthew G. Knepley     void   (*f1)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->f1IFuncs[f];
916af1eca97SMatthew G. Knepley     PetscInt numQuadPoints, Nb;
917af1eca97SMatthew G. Knepley     /* Conforming batches */
918af1eca97SMatthew G. Knepley     PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
919af1eca97SMatthew G. Knepley     /* Remainder */
920af1eca97SMatthew G. Knepley     PetscInt Nr, offset;
921af1eca97SMatthew G. Knepley 
922af1eca97SMatthew G. Knepley     ierr = PetscFEGetQuadrature(fe[f], &q);CHKERRQ(ierr);
923af1eca97SMatthew G. Knepley     ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr);
924af1eca97SMatthew G. Knepley     ierr = PetscFEGetTileSizes(fe[f], NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
925af1eca97SMatthew G. Knepley     ierr = PetscQuadratureGetData(q, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr);
926af1eca97SMatthew G. Knepley     blockSize = Nb*numQuadPoints;
927af1eca97SMatthew G. Knepley     batchSize = numBlocks * blockSize;
928af1eca97SMatthew G. Knepley     ierr =  PetscFESetTileSizes(fe[f], blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
929af1eca97SMatthew G. Knepley     numChunks = numCells / (numBatches*batchSize);
930af1eca97SMatthew G. Knepley     Ne        = numChunks*numBatches*batchSize;
931af1eca97SMatthew G. Knepley     Nr        = numCells % (numBatches*batchSize);
932af1eca97SMatthew G. Knepley     offset    = numCells - Nr;
933af1eca97SMatthew G. Knepley     geom.v0   = v0;
934af1eca97SMatthew G. Knepley     geom.J    = J;
935af1eca97SMatthew G. Knepley     geom.invJ = invJ;
936af1eca97SMatthew G. Knepley     geom.detJ = detJ;
937af1eca97SMatthew G. Knepley     ierr = PetscFEIntegrateIFunction(fe[f], Ne, Nf, fe, f, geom, u, u_t, NfAux, feAux, a, f0, f1, elemVec);CHKERRQ(ierr);
938af1eca97SMatthew G. Knepley     geom.v0   = &v0[offset*dim];
939af1eca97SMatthew G. Knepley     geom.J    = &J[offset*dim*dim];
940af1eca97SMatthew G. Knepley     geom.invJ = &invJ[offset*dim*dim];
941af1eca97SMatthew G. Knepley     geom.detJ = &detJ[offset];
942af1eca97SMatthew 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);
943af1eca97SMatthew G. Knepley   }
944af1eca97SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
945af1eca97SMatthew G. Knepley     if (mesh->printFEM > 1) {ierr = DMPrintCellVector(c, name, cellDof, &elemVec[c*cellDof]);CHKERRQ(ierr);}
946af1eca97SMatthew G. Knepley     ierr = DMPlexVecSetClosure(dm, section, F, c, &elemVec[c*cellDof], ADD_VALUES);CHKERRQ(ierr);
947af1eca97SMatthew G. Knepley   }
948af1eca97SMatthew G. Knepley   ierr = PetscFree7(u,u_t,v0,J,invJ,detJ,elemVec);CHKERRQ(ierr);
949af1eca97SMatthew G. Knepley   if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);}
950af1eca97SMatthew G. Knepley   if (feBd) {
951af1eca97SMatthew G. Knepley     DMLabel         label, depth;
952af1eca97SMatthew G. Knepley     IS              pointIS;
953af1eca97SMatthew G. Knepley     const PetscInt *points;
954af1eca97SMatthew G. Knepley     PetscInt        dep, numPoints, p, numFaces;
955af1eca97SMatthew G. Knepley     PetscReal      *n;
956af1eca97SMatthew G. Knepley 
957af1eca97SMatthew G. Knepley     ierr = DMPlexGetLabel(dm, "boundary", &label);CHKERRQ(ierr);
958af1eca97SMatthew G. Knepley     ierr = DMPlexGetDepthLabel(dm, &depth);CHKERRQ(ierr);
959af1eca97SMatthew G. Knepley     ierr = DMLabelGetStratumSize(label, 1, &numPoints);CHKERRQ(ierr);
960af1eca97SMatthew G. Knepley     ierr = DMLabelGetStratumIS(label, 1, &pointIS);CHKERRQ(ierr);
961af1eca97SMatthew G. Knepley     ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
962af1eca97SMatthew G. Knepley     for (f = 0, cellDof = 0, numComponents = 0; f < Nf; ++f) {
963af1eca97SMatthew G. Knepley       PetscInt Nb, Nc;
964af1eca97SMatthew G. Knepley 
965af1eca97SMatthew G. Knepley       ierr = PetscFEGetDimension(feBd[f], &Nb);CHKERRQ(ierr);
966af1eca97SMatthew G. Knepley       ierr = PetscFEGetNumComponents(feBd[f], &Nc);CHKERRQ(ierr);
967af1eca97SMatthew G. Knepley       cellDof       += Nb*Nc;
968af1eca97SMatthew G. Knepley       numComponents += Nc;
969af1eca97SMatthew G. Knepley     }
970af1eca97SMatthew G. Knepley     for (p = 0, numFaces = 0; p < numPoints; ++p) {
971af1eca97SMatthew G. Knepley       ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr);
972af1eca97SMatthew G. Knepley       if (dep == dim-1) ++numFaces;
973af1eca97SMatthew G. Knepley     }
974af1eca97SMatthew 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);
975af1eca97SMatthew G. Knepley     ierr = PetscMalloc1(numFaces*cellDof,&u_t);CHKERRQ(ierr);
976af1eca97SMatthew G. Knepley     for (p = 0, f = 0; p < numPoints; ++p) {
977af1eca97SMatthew G. Knepley       const PetscInt point = points[p];
978af1eca97SMatthew G. Knepley       PetscScalar   *x     = NULL;
979af1eca97SMatthew G. Knepley       PetscInt       i;
980af1eca97SMatthew G. Knepley 
981af1eca97SMatthew G. Knepley       ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr);
982af1eca97SMatthew G. Knepley       if (dep != dim-1) continue;
983af1eca97SMatthew G. Knepley       ierr = DMPlexComputeCellGeometry(dm, point, &v0[f*dim], &J[f*dim*dim], &invJ[f*dim*dim], &detJ[f]);CHKERRQ(ierr);
984af1eca97SMatthew G. Knepley       ierr = DMPlexComputeCellGeometryFVM(dm, point, NULL, NULL, &n[f*dim]);
985af1eca97SMatthew G. Knepley       if (detJ[f] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for face %d", detJ[f], point);
986af1eca97SMatthew G. Knepley       ierr = DMPlexVecGetClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr);
987af1eca97SMatthew G. Knepley       for (i = 0; i < cellDof; ++i) u[f*cellDof+i] = x[i];
988af1eca97SMatthew G. Knepley       ierr = DMPlexVecRestoreClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr);
989af1eca97SMatthew G. Knepley       ierr = DMPlexVecGetClosure(dm, section, X_t, point, NULL, &x);CHKERRQ(ierr);
990af1eca97SMatthew G. Knepley       for (i = 0; i < cellDof; ++i) u_t[f*cellDof+i] = x[i];
991af1eca97SMatthew G. Knepley       ierr = DMPlexVecRestoreClosure(dm, section, X_t, point, NULL, &x);CHKERRQ(ierr);
992af1eca97SMatthew G. Knepley       ++f;
993af1eca97SMatthew G. Knepley     }
994af1eca97SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
995af1eca97SMatthew G. Knepley       void   (*f0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->f0BdIFuncs[f];
996af1eca97SMatthew G. Knepley       void   (*f1)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->f1BdIFuncs[f];
997af1eca97SMatthew G. Knepley       PetscInt numQuadPoints, Nb;
998af1eca97SMatthew G. Knepley       /* Conforming batches */
999af1eca97SMatthew G. Knepley       PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
1000af1eca97SMatthew G. Knepley       /* Remainder */
1001af1eca97SMatthew G. Knepley       PetscInt Nr, offset;
1002af1eca97SMatthew G. Knepley 
1003af1eca97SMatthew G. Knepley       ierr = PetscFEGetQuadrature(feBd[f], &q);CHKERRQ(ierr);
1004af1eca97SMatthew G. Knepley       ierr = PetscFEGetDimension(feBd[f], &Nb);CHKERRQ(ierr);
1005af1eca97SMatthew G. Knepley       ierr = PetscFEGetTileSizes(feBd[f], NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
1006af1eca97SMatthew G. Knepley       ierr = PetscQuadratureGetData(q, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr);
1007af1eca97SMatthew G. Knepley       blockSize = Nb*numQuadPoints;
1008af1eca97SMatthew G. Knepley       batchSize = numBlocks * blockSize;
1009af1eca97SMatthew G. Knepley       ierr =  PetscFESetTileSizes(feBd[f], blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
1010af1eca97SMatthew G. Knepley       numChunks = numFaces / (numBatches*batchSize);
1011af1eca97SMatthew G. Knepley       Ne        = numChunks*numBatches*batchSize;
1012af1eca97SMatthew G. Knepley       Nr        = numFaces % (numBatches*batchSize);
1013af1eca97SMatthew G. Knepley       offset    = numFaces - Nr;
1014af1eca97SMatthew G. Knepley       geom.v0   = v0;
1015af1eca97SMatthew G. Knepley       geom.n    = n;
1016af1eca97SMatthew G. Knepley       geom.J    = J;
1017af1eca97SMatthew G. Knepley       geom.invJ = invJ;
1018af1eca97SMatthew G. Knepley       geom.detJ = detJ;
1019af1eca97SMatthew G. Knepley       ierr = PetscFEIntegrateBdIFunction(feBd[f], Ne, Nf, feBd, f, geom, u, u_t, 0, NULL, NULL, f0, f1, elemVec);CHKERRQ(ierr);
1020af1eca97SMatthew G. Knepley       geom.v0   = &v0[offset*dim];
1021af1eca97SMatthew G. Knepley       geom.n    = &n[offset*dim];
1022af1eca97SMatthew G. Knepley       geom.J    = &J[offset*dim*dim];
1023af1eca97SMatthew G. Knepley       geom.invJ = &invJ[offset*dim*dim];
1024af1eca97SMatthew G. Knepley       geom.detJ = &detJ[offset];
1025af1eca97SMatthew 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);
1026af1eca97SMatthew G. Knepley     }
1027af1eca97SMatthew G. Knepley     for (p = 0, f = 0; p < numPoints; ++p) {
1028af1eca97SMatthew G. Knepley       const PetscInt point = points[p];
1029af1eca97SMatthew G. Knepley 
1030af1eca97SMatthew G. Knepley       ierr = DMLabelGetValue(depth, point, &dep);CHKERRQ(ierr);
1031af1eca97SMatthew G. Knepley       if (dep != dim-1) continue;
1032af1eca97SMatthew G. Knepley       if (mesh->printFEM > 1) {ierr = DMPrintCellVector(point, "BdResidual", cellDof, &elemVec[f*cellDof]);CHKERRQ(ierr);}
1033af1eca97SMatthew G. Knepley       ierr = DMPlexVecSetClosure(dm, NULL, F, point, &elemVec[f*cellDof], ADD_VALUES);CHKERRQ(ierr);
1034af1eca97SMatthew G. Knepley       ++f;
1035af1eca97SMatthew G. Knepley     }
1036af1eca97SMatthew G. Knepley     ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
1037af1eca97SMatthew G. Knepley     ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
1038af1eca97SMatthew G. Knepley     ierr = PetscFree7(u,v0,n,J,invJ,detJ,elemVec);CHKERRQ(ierr);
1039af1eca97SMatthew G. Knepley     ierr = PetscFree(u_t);CHKERRQ(ierr);
1040af1eca97SMatthew G. Knepley   }
1041af1eca97SMatthew G. Knepley   if (mesh->printFEM) {ierr = DMPrintLocalVec(dm, name, mesh->printTol, F);CHKERRQ(ierr);}
1042af1eca97SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_ResidualFEM,dm,0,0,0);CHKERRQ(ierr);
1043af1eca97SMatthew G. Knepley   PetscFunctionReturn(0);
1044af1eca97SMatthew G. Knepley }
1045af1eca97SMatthew G. Knepley 
1046af1eca97SMatthew G. Knepley #undef __FUNCT__
1047cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexComputeJacobianActionFEM"
1048cb1e1211SMatthew G Knepley /*@C
1049cb1e1211SMatthew G Knepley   DMPlexComputeJacobianActionFEM - Form the local action of Jacobian J(u) on the local input X using pointwise functions specified by the user
1050cb1e1211SMatthew G Knepley 
1051cb1e1211SMatthew G Knepley   Input Parameters:
1052cb1e1211SMatthew G Knepley + dm - The mesh
1053cb1e1211SMatthew G Knepley . J  - The Jacobian shell matrix
1054cb1e1211SMatthew G Knepley . X  - Local input vector
1055cb1e1211SMatthew G Knepley - user - The user context
1056cb1e1211SMatthew G Knepley 
1057cb1e1211SMatthew G Knepley   Output Parameter:
1058cb1e1211SMatthew G Knepley . F  - Local output vector
1059cb1e1211SMatthew G Knepley 
1060cb1e1211SMatthew G Knepley   Note:
10618026896cSMatthew G. Knepley   The first member of the user context must be an FEMContext.
1062cb1e1211SMatthew G Knepley 
1063cb1e1211SMatthew G Knepley   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
1064cb1e1211SMatthew G Knepley   like a GPU, or vectorize on a multicore machine.
1065cb1e1211SMatthew G Knepley 
10660059ad2aSSatish Balay   Level: developer
10670059ad2aSSatish Balay 
1068cb1e1211SMatthew G Knepley .seealso: DMPlexComputeResidualFEM()
1069878cb397SSatish Balay @*/
1070cb1e1211SMatthew G Knepley PetscErrorCode DMPlexComputeJacobianActionFEM(DM dm, Mat Jac, Vec X, Vec F, void *user)
1071cb1e1211SMatthew G Knepley {
1072cb1e1211SMatthew G Knepley   DM_Plex          *mesh = (DM_Plex *) dm->data;
10739a559087SMatthew G. Knepley   PetscFEM         *fem  = (PetscFEM *) user;
10740483ade4SMatthew G. Knepley   PetscFE          *fe   = fem->fe;
10750483ade4SMatthew G. Knepley   PetscQuadrature   quad;
10760483ade4SMatthew G. Knepley   PetscCellGeometry geom;
1077cb1e1211SMatthew G Knepley   PetscSection      section;
1078cb1e1211SMatthew G Knepley   JacActionCtx     *jctx;
1079cb1e1211SMatthew G Knepley   PetscReal        *v0, *J, *invJ, *detJ;
1080cb1e1211SMatthew G Knepley   PetscScalar      *elemVec, *u, *a;
10810483ade4SMatthew G. Knepley   PetscInt          dim, numFields, field, numCells, cStart, cEnd, c;
1082cb1e1211SMatthew G Knepley   PetscInt          cellDof = 0;
1083cb1e1211SMatthew G Knepley   PetscErrorCode    ierr;
1084cb1e1211SMatthew G Knepley 
1085cb1e1211SMatthew G Knepley   PetscFunctionBegin;
10860483ade4SMatthew G. Knepley   /* ierr = PetscLogEventBegin(DMPLEX_JacobianActionFEM,dm,0,0,0);CHKERRQ(ierr); */
1087cb1e1211SMatthew G Knepley   ierr = MatShellGetContext(Jac, &jctx);CHKERRQ(ierr);
1088cb1e1211SMatthew G Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
1089cb1e1211SMatthew G Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
1090cb1e1211SMatthew G Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
1091cb1e1211SMatthew G Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
1092cb1e1211SMatthew G Knepley   numCells = cEnd - cStart;
1093cb1e1211SMatthew G Knepley   for (field = 0; field < numFields; ++field) {
10940483ade4SMatthew G. Knepley     PetscInt Nb, Nc;
10950483ade4SMatthew G. Knepley 
10960483ade4SMatthew G. Knepley     ierr = PetscFEGetDimension(fe[field], &Nb);CHKERRQ(ierr);
10970483ade4SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe[field], &Nc);CHKERRQ(ierr);
10980483ade4SMatthew G. Knepley     cellDof += Nb*Nc;
1099cb1e1211SMatthew G Knepley   }
1100cb1e1211SMatthew G Knepley   ierr = VecSet(F, 0.0);CHKERRQ(ierr);
1101dcca6d9dSJed 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);
1102cb1e1211SMatthew G Knepley   for (c = cStart; c < cEnd; ++c) {
1103a1e44745SMatthew G. Knepley     PetscScalar *x = NULL;
1104cb1e1211SMatthew G Knepley     PetscInt     i;
1105cb1e1211SMatthew G Knepley 
1106cb1e1211SMatthew G Knepley     ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr);
1107cb1e1211SMatthew G Knepley     if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c);
1108cb1e1211SMatthew G Knepley     ierr = DMPlexVecGetClosure(dm, NULL, jctx->u, c, NULL, &x);CHKERRQ(ierr);
1109cb1e1211SMatthew G Knepley     for (i = 0; i < cellDof; ++i) u[c*cellDof+i] = x[i];
1110cb1e1211SMatthew G Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, jctx->u, c, NULL, &x);CHKERRQ(ierr);
1111cb1e1211SMatthew G Knepley     ierr = DMPlexVecGetClosure(dm, NULL, X, c, NULL, &x);CHKERRQ(ierr);
1112cb1e1211SMatthew G Knepley     for (i = 0; i < cellDof; ++i) a[c*cellDof+i] = x[i];
1113cb1e1211SMatthew G Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, X, c, NULL, &x);CHKERRQ(ierr);
1114cb1e1211SMatthew G Knepley   }
1115cb1e1211SMatthew G Knepley   for (field = 0; field < numFields; ++field) {
111621454ff5SMatthew G. Knepley     PetscInt numQuadPoints, Nb;
1117cb1e1211SMatthew G Knepley     /* Conforming batches */
1118cb1e1211SMatthew G Knepley     PetscInt numBlocks  = 1;
11190483ade4SMatthew G. Knepley     PetscInt numBatches = 1;
11200483ade4SMatthew G. Knepley     PetscInt numChunks, Ne, blockSize, batchSize;
1121cb1e1211SMatthew G Knepley     /* Remainder */
11220483ade4SMatthew G. Knepley     PetscInt Nr, offset;
1123cb1e1211SMatthew G Knepley 
11240483ade4SMatthew G. Knepley     ierr = PetscFEGetQuadrature(fe[field], &quad);CHKERRQ(ierr);
11250483ade4SMatthew G. Knepley     ierr = PetscFEGetDimension(fe[field], &Nb);CHKERRQ(ierr);
112621454ff5SMatthew G. Knepley     ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr);
112721454ff5SMatthew G. Knepley     blockSize = Nb*numQuadPoints;
11280483ade4SMatthew G. Knepley     batchSize = numBlocks * blockSize;
11290483ade4SMatthew G. Knepley     numChunks = numCells / (numBatches*batchSize);
11300483ade4SMatthew G. Knepley     Ne        = numChunks*numBatches*batchSize;
11310483ade4SMatthew G. Knepley     Nr        = numCells % (numBatches*batchSize);
11320483ade4SMatthew G. Knepley     offset    = numCells - Nr;
11330483ade4SMatthew G. Knepley     geom.v0   = v0;
11340483ade4SMatthew G. Knepley     geom.J    = J;
11350483ade4SMatthew G. Knepley     geom.invJ = invJ;
11360483ade4SMatthew G. Knepley     geom.detJ = detJ;
11370483ade4SMatthew G. Knepley     ierr = PetscFEIntegrateJacobianAction(fe[field], Ne, numFields, fe, field, geom, u, a, fem->g0Funcs, fem->g1Funcs, fem->g2Funcs, fem->g3Funcs, elemVec);CHKERRQ(ierr);
11380483ade4SMatthew G. Knepley     geom.v0   = &v0[offset*dim];
11390483ade4SMatthew G. Knepley     geom.J    = &J[offset*dim*dim];
11400483ade4SMatthew G. Knepley     geom.invJ = &invJ[offset*dim*dim];
11410483ade4SMatthew G. Knepley     geom.detJ = &detJ[offset];
11420483ade4SMatthew G. Knepley     ierr = PetscFEIntegrateJacobianAction(fe[field], Nr, numFields, fe, field, geom, &u[offset*cellDof], &a[offset*cellDof],
1143cb1e1211SMatthew G Knepley                                           fem->g0Funcs, fem->g1Funcs, fem->g2Funcs, fem->g3Funcs, &elemVec[offset*cellDof]);CHKERRQ(ierr);
1144cb1e1211SMatthew G Knepley   }
1145cb1e1211SMatthew G Knepley   for (c = cStart; c < cEnd; ++c) {
1146cb1e1211SMatthew G Knepley     if (mesh->printFEM > 1) {ierr = DMPrintCellVector(c, "Jacobian Action", cellDof, &elemVec[c*cellDof]);CHKERRQ(ierr);}
1147cb1e1211SMatthew G Knepley     ierr = DMPlexVecSetClosure(dm, NULL, F, c, &elemVec[c*cellDof], ADD_VALUES);CHKERRQ(ierr);
1148cb1e1211SMatthew G Knepley   }
1149cb1e1211SMatthew G Knepley   ierr = PetscFree7(u,a,v0,J,invJ,detJ,elemVec);CHKERRQ(ierr);
1150cb1e1211SMatthew G Knepley   if (mesh->printFEM) {
1151cb1e1211SMatthew G Knepley     PetscMPIInt rank, numProcs;
1152cb1e1211SMatthew G Knepley     PetscInt    p;
1153cb1e1211SMatthew G Knepley 
1154cb1e1211SMatthew G Knepley     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr);
1155cb1e1211SMatthew G Knepley     ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm), &numProcs);CHKERRQ(ierr);
115686a74ee0SMatthew G. Knepley     ierr = PetscPrintf(PetscObjectComm((PetscObject)dm), "Jacobian Action:\n");CHKERRQ(ierr);
1157cb1e1211SMatthew G Knepley     for (p = 0; p < numProcs; ++p) {
1158cb1e1211SMatthew G Knepley       if (p == rank) {ierr = VecView(F, PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);}
1159cb1e1211SMatthew G Knepley       ierr = PetscBarrier((PetscObject) dm);CHKERRQ(ierr);
1160cb1e1211SMatthew G Knepley     }
1161cb1e1211SMatthew G Knepley   }
11620483ade4SMatthew G. Knepley   /* ierr = PetscLogEventEnd(DMPLEX_JacobianActionFEM,dm,0,0,0);CHKERRQ(ierr); */
1163cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
1164cb1e1211SMatthew G Knepley }
1165cb1e1211SMatthew G Knepley 
1166cb1e1211SMatthew G Knepley #undef __FUNCT__
1167cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexComputeJacobianFEM"
1168cb1e1211SMatthew G Knepley /*@
1169cb1e1211SMatthew G Knepley   DMPlexComputeJacobianFEM - Form the local portion of the Jacobian matrix J at the local solution X using pointwise functions specified by the user.
1170cb1e1211SMatthew G Knepley 
1171cb1e1211SMatthew G Knepley   Input Parameters:
1172cb1e1211SMatthew G Knepley + dm - The mesh
1173cb1e1211SMatthew G Knepley . X  - Local input vector
1174cb1e1211SMatthew G Knepley - user - The user context
1175cb1e1211SMatthew G Knepley 
1176cb1e1211SMatthew G Knepley   Output Parameter:
1177cb1e1211SMatthew G Knepley . Jac  - Jacobian matrix
1178cb1e1211SMatthew G Knepley 
1179cb1e1211SMatthew G Knepley   Note:
11808026896cSMatthew G. Knepley   The first member of the user context must be an FEMContext.
1181cb1e1211SMatthew G Knepley 
1182cb1e1211SMatthew G Knepley   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
1183cb1e1211SMatthew G Knepley   like a GPU, or vectorize on a multicore machine.
1184cb1e1211SMatthew G Knepley 
11850059ad2aSSatish Balay   Level: developer
11860059ad2aSSatish Balay 
1187cb1e1211SMatthew G Knepley .seealso: FormFunctionLocal()
1188878cb397SSatish Balay @*/
11890405ed22SMatthew G. Knepley PetscErrorCode DMPlexComputeJacobianFEM(DM dm, Vec X, Mat Jac, Mat JacP,void *user)
1190cb1e1211SMatthew G Knepley {
1191cb1e1211SMatthew G Knepley   DM_Plex          *mesh  = (DM_Plex *) dm->data;
11929a559087SMatthew G. Knepley   PetscFEM         *fem   = (PetscFEM *) user;
1193a319912fSMatthew G. Knepley   PetscFE          *fe    = fem->fe;
1194754551f4SMatthew G. Knepley   PetscFE          *feAux = fem->feAux;
1195a319912fSMatthew G. Knepley   const char       *name  = "Jacobian";
1196754551f4SMatthew G. Knepley   DM                dmAux;
1197754551f4SMatthew G. Knepley   Vec               A;
1198a319912fSMatthew G. Knepley   PetscQuadrature   quad;
1199a319912fSMatthew G. Knepley   PetscCellGeometry geom;
1200754551f4SMatthew G. Knepley   PetscSection      section, globalSection, sectionAux;
1201cb1e1211SMatthew G Knepley   PetscReal        *v0, *J, *invJ, *detJ;
1202754551f4SMatthew G. Knepley   PetscScalar      *elemMat, *u, *a;
1203754551f4SMatthew G. Knepley   PetscInt          dim, Nf, NfAux = 0, f, fieldI, fieldJ, numCells, cStart, cEnd, c;
1204cb1e1211SMatthew G Knepley   PetscInt          cellDof = 0, numComponents = 0;
1205754551f4SMatthew G. Knepley   PetscInt          cellDofAux = 0, numComponentsAux = 0;
1206cb1e1211SMatthew G Knepley   PetscBool         isShell;
1207cb1e1211SMatthew G Knepley   PetscErrorCode    ierr;
1208cb1e1211SMatthew G Knepley 
1209cb1e1211SMatthew G Knepley   PetscFunctionBegin;
1210a319912fSMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_JacobianFEM,dm,0,0,0);CHKERRQ(ierr);
1211cb1e1211SMatthew G Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
1212cb1e1211SMatthew G Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
1213a319912fSMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dm, &globalSection);CHKERRQ(ierr);
1214754551f4SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
1215cb1e1211SMatthew G Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
1216cb1e1211SMatthew G Knepley   numCells = cEnd - cStart;
1217754551f4SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
1218a319912fSMatthew G. Knepley     PetscInt Nb, Nc;
1219a319912fSMatthew G. Knepley 
1220a319912fSMatthew G. Knepley     ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr);
1221a319912fSMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe[f], &Nc);CHKERRQ(ierr);
1222a319912fSMatthew G. Knepley     cellDof       += Nb*Nc;
1223a319912fSMatthew G. Knepley     numComponents += Nc;
1224cb1e1211SMatthew G Knepley   }
1225754551f4SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
1226754551f4SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr);
1227754551f4SMatthew G. Knepley   if (dmAux) {
1228754551f4SMatthew G. Knepley     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
1229754551f4SMatthew G. Knepley     ierr = PetscSectionGetNumFields(sectionAux, &NfAux);CHKERRQ(ierr);
1230754551f4SMatthew G. Knepley   }
1231754551f4SMatthew G. Knepley   for (f = 0; f < NfAux; ++f) {
1232754551f4SMatthew G. Knepley     PetscInt Nb, Nc;
1233754551f4SMatthew G. Knepley 
1234754551f4SMatthew G. Knepley     ierr = PetscFEGetDimension(feAux[f], &Nb);CHKERRQ(ierr);
1235754551f4SMatthew G. Knepley     ierr = PetscFEGetNumComponents(feAux[f], &Nc);CHKERRQ(ierr);
1236754551f4SMatthew G. Knepley     cellDofAux       += Nb*Nc;
1237754551f4SMatthew G. Knepley     numComponentsAux += Nc;
1238754551f4SMatthew G. Knepley   }
1239bd5ce809SMatthew G. Knepley   ierr = DMPlexInsertBoundaryValues(dm, X);CHKERRQ(ierr);
1240cb1e1211SMatthew G Knepley   ierr = MatZeroEntries(JacP);CHKERRQ(ierr);
1241dcca6d9dSJed 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);
1242785e854fSJed Brown   if (dmAux) {ierr = PetscMalloc1(numCells*cellDofAux, &a);CHKERRQ(ierr);}
1243cb1e1211SMatthew G Knepley   for (c = cStart; c < cEnd; ++c) {
1244a1e44745SMatthew G. Knepley     PetscScalar *x = NULL;
1245cb1e1211SMatthew G Knepley     PetscInt     i;
1246cb1e1211SMatthew G Knepley 
1247cb1e1211SMatthew G Knepley     ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr);
1248cb1e1211SMatthew G Knepley     if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c);
1249a319912fSMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr);
1250cb1e1211SMatthew G Knepley     for (i = 0; i < cellDof; ++i) u[c*cellDof+i] = x[i];
1251a319912fSMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr);
1252754551f4SMatthew G. Knepley     if (dmAux) {
1253754551f4SMatthew G. Knepley       ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
1254754551f4SMatthew G. Knepley       for (i = 0; i < cellDofAux; ++i) a[c*cellDofAux+i] = x[i];
1255754551f4SMatthew G. Knepley       ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
1256754551f4SMatthew G. Knepley     }
1257cb1e1211SMatthew G Knepley   }
1258cb1e1211SMatthew G Knepley   ierr = PetscMemzero(elemMat, numCells*cellDof*cellDof * sizeof(PetscScalar));CHKERRQ(ierr);
1259754551f4SMatthew G. Knepley   for (fieldI = 0; fieldI < Nf; ++fieldI) {
126021454ff5SMatthew G. Knepley     PetscInt numQuadPoints, Nb;
1261cb1e1211SMatthew G Knepley     /* Conforming batches */
1262754551f4SMatthew G. Knepley     PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
1263cb1e1211SMatthew G Knepley     /* Remainder */
1264a319912fSMatthew G. Knepley     PetscInt Nr, offset;
1265cb1e1211SMatthew G Knepley 
1266754551f4SMatthew G. Knepley     ierr = PetscFEGetQuadrature(fe[fieldI], &quad);CHKERRQ(ierr);
1267754551f4SMatthew G. Knepley     ierr = PetscFEGetDimension(fe[fieldI], &Nb);CHKERRQ(ierr);
1268754551f4SMatthew G. Knepley     ierr = PetscFEGetTileSizes(fe[fieldI], NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
126921454ff5SMatthew G. Knepley     ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr);
127021454ff5SMatthew G. Knepley     blockSize = Nb*numQuadPoints;
1271a319912fSMatthew G. Knepley     batchSize = numBlocks * blockSize;
1272754551f4SMatthew G. Knepley     ierr = PetscFESetTileSizes(fe[fieldI], blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
1273a319912fSMatthew G. Knepley     numChunks = numCells / (numBatches*batchSize);
1274a319912fSMatthew G. Knepley     Ne        = numChunks*numBatches*batchSize;
1275a319912fSMatthew G. Knepley     Nr        = numCells % (numBatches*batchSize);
1276a319912fSMatthew G. Knepley     offset    = numCells - Nr;
1277754551f4SMatthew G. Knepley     for (fieldJ = 0; fieldJ < Nf; ++fieldJ) {
1278754551f4SMatthew G. Knepley       void   (*g0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g0Funcs[fieldI*Nf+fieldJ];
1279754551f4SMatthew G. Knepley       void   (*g1)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g1Funcs[fieldI*Nf+fieldJ];
1280754551f4SMatthew G. Knepley       void   (*g2)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g2Funcs[fieldI*Nf+fieldJ];
1281754551f4SMatthew G. Knepley       void   (*g3)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g3Funcs[fieldI*Nf+fieldJ];
1282754551f4SMatthew G. Knepley 
1283a319912fSMatthew G. Knepley       geom.v0   = v0;
1284a319912fSMatthew G. Knepley       geom.J    = J;
1285a319912fSMatthew G. Knepley       geom.invJ = invJ;
1286a319912fSMatthew G. Knepley       geom.detJ = detJ;
1287754551f4SMatthew G. Knepley       ierr = PetscFEIntegrateJacobian(fe[fieldI], Ne, Nf, fe, fieldI, fieldJ, geom, u, NfAux, feAux, a, g0, g1, g2, g3, elemMat);CHKERRQ(ierr);
1288a319912fSMatthew G. Knepley       geom.v0   = &v0[offset*dim];
1289a319912fSMatthew G. Knepley       geom.J    = &J[offset*dim*dim];
1290a319912fSMatthew G. Knepley       geom.invJ = &invJ[offset*dim*dim];
1291a319912fSMatthew G. Knepley       geom.detJ = &detJ[offset];
1292754551f4SMatthew 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);
1293cb1e1211SMatthew G Knepley     }
1294cb1e1211SMatthew G Knepley   }
1295cb1e1211SMatthew G Knepley   for (c = cStart; c < cEnd; ++c) {
1296a319912fSMatthew G. Knepley     if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(c, name, cellDof, cellDof, &elemMat[c*cellDof*cellDof]);CHKERRQ(ierr);}
1297a319912fSMatthew G. Knepley     ierr = DMPlexMatSetClosure(dm, section, globalSection, JacP, c, &elemMat[c*cellDof*cellDof], ADD_VALUES);CHKERRQ(ierr);
1298cb1e1211SMatthew G Knepley   }
1299cb1e1211SMatthew G Knepley   ierr = PetscFree6(u,v0,J,invJ,detJ,elemMat);CHKERRQ(ierr);
1300754551f4SMatthew G. Knepley   if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);}
1301cb1e1211SMatthew G Knepley   ierr = MatAssemblyBegin(JacP, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1302cb1e1211SMatthew G Knepley   ierr = MatAssemblyEnd(JacP, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1303cb1e1211SMatthew G Knepley   if (mesh->printFEM) {
1304a319912fSMatthew G. Knepley     ierr = PetscPrintf(PETSC_COMM_WORLD, "%s:\n", name);CHKERRQ(ierr);
1305cb1e1211SMatthew G Knepley     ierr = MatChop(JacP, 1.0e-10);CHKERRQ(ierr);
1306cb1e1211SMatthew G Knepley     ierr = MatView(JacP, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1307cb1e1211SMatthew G Knepley   }
1308a319912fSMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_JacobianFEM,dm,0,0,0);CHKERRQ(ierr);
1309cb1e1211SMatthew G Knepley   ierr = PetscObjectTypeCompare((PetscObject) Jac, MATSHELL, &isShell);CHKERRQ(ierr);
1310cb1e1211SMatthew G Knepley   if (isShell) {
1311cb1e1211SMatthew G Knepley     JacActionCtx *jctx;
1312cb1e1211SMatthew G Knepley 
1313cb1e1211SMatthew G Knepley     ierr = MatShellGetContext(Jac, &jctx);CHKERRQ(ierr);
1314cb1e1211SMatthew G Knepley     ierr = VecCopy(X, jctx->u);CHKERRQ(ierr);
1315cb1e1211SMatthew G Knepley   }
1316cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
1317cb1e1211SMatthew G Knepley }
1318bceba477SMatthew G. Knepley 
1319d69c5d34SMatthew G. Knepley #if 0
1320d69c5d34SMatthew G. Knepley 
132132356553SMatthew 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[])
132232356553SMatthew G. Knepley {
132332356553SMatthew G. Knepley   g3[0] = 1.0;
132432356553SMatthew G. Knepley }
132532356553SMatthew G. Knepley 
132632356553SMatthew 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[])
132732356553SMatthew G. Knepley {
132832356553SMatthew G. Knepley   g3[0] = g3[3] = 1.0;
132932356553SMatthew G. Knepley }
133032356553SMatthew G. Knepley 
133132356553SMatthew 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[])
133232356553SMatthew G. Knepley {
133332356553SMatthew G. Knepley   g3[0] = g3[4] = g3[8] = 1.0;
133432356553SMatthew G. Knepley }
133532356553SMatthew G. Knepley 
1336bceba477SMatthew G. Knepley #undef __FUNCT__
1337d69c5d34SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeInterpolatorFEMBroken"
1338d69c5d34SMatthew G. Knepley PetscErrorCode DMPlexComputeInterpolatorFEMBroken(DM dmc, DM dmf, Mat I, void *user)
1339bceba477SMatthew G. Knepley {
1340bceba477SMatthew G. Knepley   DM_Plex          *mesh  = (DM_Plex *) dmc->data;
1341bceba477SMatthew G. Knepley   PetscFEM         *fem   = (PetscFEM *) user;
1342bceba477SMatthew G. Knepley   PetscFE          *fe    = fem->fe;
1343bceba477SMatthew G. Knepley   const char       *name  = "Interpolator";
134432356553SMatthew G. Knepley   PetscFE          *feRef;
1345bceba477SMatthew G. Knepley   PetscQuadrature   quad, quadOld;
1346bceba477SMatthew G. Knepley   PetscCellGeometry geom;
134732356553SMatthew G. Knepley   PetscSection      fsection, fglobalSection;
134832356553SMatthew G. Knepley   PetscSection      csection, cglobalSection;
1349bceba477SMatthew G. Knepley   PetscReal        *v0, *J, *invJ, *detJ;
1350bceba477SMatthew G. Knepley   PetscScalar      *elemMat;
1351bceba477SMatthew G. Knepley   PetscInt          dim, Nf, f, fieldI, fieldJ, numCells, cStart, cEnd, c;
135232356553SMatthew G. Knepley   PetscInt          rCellDof = 0, cCellDof = 0, numComponents = 0;
1353bceba477SMatthew G. Knepley   PetscErrorCode    ierr;
1354bceba477SMatthew G. Knepley 
1355bceba477SMatthew G. Knepley   PetscFunctionBegin;
1356bceba477SMatthew G. Knepley #if 0
1357bceba477SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1358bceba477SMatthew G. Knepley #endif
1359bceba477SMatthew G. Knepley   ierr = DMPlexGetDimension(dmf, &dim);CHKERRQ(ierr);
136032356553SMatthew G. Knepley   ierr = DMGetDefaultSection(dmf, &fsection);CHKERRQ(ierr);
136132356553SMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dmf, &fglobalSection);CHKERRQ(ierr);
136232356553SMatthew G. Knepley   ierr = DMGetDefaultSection(dmc, &csection);CHKERRQ(ierr);
136332356553SMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dmc, &cglobalSection);CHKERRQ(ierr);
136432356553SMatthew G. Knepley   ierr = PetscSectionGetNumFields(fsection, &Nf);CHKERRQ(ierr);
1365bceba477SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dmc, 0, &cStart, &cEnd);CHKERRQ(ierr);
1366bceba477SMatthew G. Knepley   numCells = cEnd - cStart;
136732356553SMatthew G. Knepley   ierr = PetscMalloc1(Nf,&feRef);CHKERRQ(ierr);
136832356553SMatthew G. Knepley   for (fieldI = 0; fieldI < Nf; ++fieldI) {
136932356553SMatthew G. Knepley     ierr = PetscFERefine(fe[fieldI], &feRef[fieldI]);CHKERRQ(ierr);
137032356553SMatthew G. Knepley   }
1371bceba477SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
137232356553SMatthew G. Knepley     PetscInt rNb, cNb, Nc;
1373bceba477SMatthew G. Knepley 
137432356553SMatthew G. Knepley     ierr = PetscFEGetDimension(feRef[f], &rNb);CHKERRQ(ierr);
137532356553SMatthew G. Knepley     ierr = PetscFEGetDimension(fe[f], &cNb);CHKERRQ(ierr);
1376bceba477SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe[f], &Nc);CHKERRQ(ierr);
1377bceba477SMatthew G. Knepley     numComponents += Nc;
137832356553SMatthew G. Knepley     rCellDof += rNb*Nc;
137932356553SMatthew G. Knepley     cCellDof += cNb*Nc;
1380bceba477SMatthew G. Knepley   }
1381bceba477SMatthew G. Knepley   ierr = MatZeroEntries(I);CHKERRQ(ierr);
138232356553SMatthew G. Knepley   ierr = PetscMalloc5(numCells*dim,&v0,numCells*dim*dim,&J,numCells*dim*dim,&invJ,numCells,&detJ,numCells*rCellDof*cCellDof,&elemMat);CHKERRQ(ierr);
1383bceba477SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
1384bceba477SMatthew G. Knepley     ierr = DMPlexComputeCellGeometry(dmc, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr);
1385bceba477SMatthew G. Knepley     if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c);
1386bceba477SMatthew G. Knepley   }
138732356553SMatthew G. Knepley   ierr = PetscMemzero(elemMat, numCells*rCellDof*cCellDof * sizeof(PetscScalar));CHKERRQ(ierr);
1388bceba477SMatthew G. Knepley   for (fieldI = 0; fieldI < Nf; ++fieldI) {
1389bceba477SMatthew G. Knepley     PetscInt numQuadPoints, Nb;
1390bceba477SMatthew G. Knepley     /* Conforming batches */
1391bceba477SMatthew G. Knepley     PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
1392bceba477SMatthew G. Knepley     /* Remainder */
1393bceba477SMatthew G. Knepley     PetscInt Nr, offset;
1394bceba477SMatthew G. Knepley 
1395bceba477SMatthew G. Knepley     /* Make new fine FE which refines the ref cell and the quadrature rule */
139632356553SMatthew G. Knepley     ierr = PetscFEGetQuadrature(feRef[fieldI], &quad);CHKERRQ(ierr);
139732356553SMatthew G. Knepley     ierr = PetscFEGetDimension(feRef[fieldI], &Nb);CHKERRQ(ierr);
139832356553SMatthew G. Knepley     ierr = PetscFEGetTileSizes(feRef[fieldI], NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
1399bceba477SMatthew G. Knepley     ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr);
1400bceba477SMatthew G. Knepley     blockSize = Nb*numQuadPoints;
1401bceba477SMatthew G. Knepley     batchSize = numBlocks * blockSize;
140232356553SMatthew G. Knepley     ierr = PetscFESetTileSizes(feRef[fieldI], blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
1403bceba477SMatthew G. Knepley     numChunks = numCells / (numBatches*batchSize);
1404bceba477SMatthew G. Knepley     Ne        = numChunks*numBatches*batchSize;
1405bceba477SMatthew G. Knepley     Nr        = numCells % (numBatches*batchSize);
1406bceba477SMatthew G. Knepley     offset    = numCells - Nr;
1407d69c5d34SMatthew G. Knepley 
1408bceba477SMatthew G. Knepley     for (fieldJ = 0; fieldJ < Nf; ++fieldJ) {
140932356553SMatthew G. Knepley       /* void   (*g0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g0Funcs[fieldI*Nf+fieldJ]; */
141032356553SMatthew G. Knepley       void   (*g0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = g0_identity_2d_static;
1411bceba477SMatthew G. Knepley 
1412bceba477SMatthew G. Knepley       /* Replace quadrature in coarse FE with refined quadrature */
1413bceba477SMatthew G. Knepley       ierr = PetscFEGetQuadrature(fe[fieldJ], &quadOld);CHKERRQ(ierr);
1414bceba477SMatthew G. Knepley       ierr = PetscObjectReference((PetscObject) quadOld);CHKERRQ(ierr);
1415bceba477SMatthew G. Knepley       ierr = PetscFESetQuadrature(fe[fieldJ], quad);CHKERRQ(ierr);
1416bceba477SMatthew G. Knepley       geom.v0   = v0;
1417bceba477SMatthew G. Knepley       geom.J    = J;
1418bceba477SMatthew G. Knepley       geom.invJ = invJ;
1419bceba477SMatthew G. Knepley       geom.detJ = detJ;
142032356553SMatthew G. Knepley       ierr = PetscFEIntegrateInterpolator_Basic(feRef[fieldI], Ne, Nf, feRef, fieldI, fe, fieldJ, geom, g0, elemMat);CHKERRQ(ierr);
1421bceba477SMatthew G. Knepley       geom.v0   = &v0[offset*dim];
1422bceba477SMatthew G. Knepley       geom.J    = &J[offset*dim*dim];
1423bceba477SMatthew G. Knepley       geom.invJ = &invJ[offset*dim*dim];
1424bceba477SMatthew G. Knepley       geom.detJ = &detJ[offset];
142532356553SMatthew G. Knepley       ierr = PetscFEIntegrateInterpolator_Basic(feRef[fieldI], Nr, Nf, feRef, fieldI, fe, fieldJ, geom, g0, &elemMat[offset*rCellDof*cCellDof]);CHKERRQ(ierr);
1426bceba477SMatthew G. Knepley       ierr = PetscFESetQuadrature(fe[fieldJ], quadOld);CHKERRQ(ierr);
1427bceba477SMatthew G. Knepley     }
1428bceba477SMatthew G. Knepley   }
1429bceba477SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
143032356553SMatthew G. Knepley     if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(c, name, rCellDof, cCellDof, &elemMat[c*rCellDof*cCellDof]);CHKERRQ(ierr);}
143132356553SMatthew G. Knepley     ierr = DMPlexMatSetClosureRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, I, c, &elemMat[c*rCellDof*cCellDof], ADD_VALUES);CHKERRQ(ierr);
1432bceba477SMatthew G. Knepley   }
1433bceba477SMatthew G. Knepley   ierr = PetscFree5(v0,J,invJ,detJ,elemMat);CHKERRQ(ierr);
143432356553SMatthew G. Knepley   ierr = PetscFree(feRef);CHKERRQ(ierr);
1435bceba477SMatthew G. Knepley   ierr = MatAssemblyBegin(I, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1436bceba477SMatthew G. Knepley   ierr = MatAssemblyEnd(I, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1437bceba477SMatthew G. Knepley   if (mesh->printFEM) {
1438bceba477SMatthew G. Knepley     ierr = PetscPrintf(PETSC_COMM_WORLD, "%s:\n", name);CHKERRQ(ierr);
1439bceba477SMatthew G. Knepley     ierr = MatChop(I, 1.0e-10);CHKERRQ(ierr);
1440bceba477SMatthew G. Knepley     ierr = MatView(I, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1441bceba477SMatthew G. Knepley   }
1442bceba477SMatthew G. Knepley #if 0
1443bceba477SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1444bceba477SMatthew G. Knepley #endif
1445bceba477SMatthew G. Knepley   PetscFunctionReturn(0);
1446bceba477SMatthew G. Knepley }
1447d69c5d34SMatthew G. Knepley #endif
1448d69c5d34SMatthew G. Knepley 
1449d69c5d34SMatthew G. Knepley #undef __FUNCT__
1450d69c5d34SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeInterpolatorFEM"
1451d69c5d34SMatthew G. Knepley /*@
1452d69c5d34SMatthew G. Knepley   DMPlexComputeInterpolatorFEM - Form the local portion of the interpolation matrix I from the coarse DM to the uniformly refined DM.
1453d69c5d34SMatthew G. Knepley 
1454d69c5d34SMatthew G. Knepley   Input Parameters:
1455d69c5d34SMatthew G. Knepley + dmf  - The fine mesh
1456d69c5d34SMatthew G. Knepley . dmc  - The coarse mesh
1457d69c5d34SMatthew G. Knepley - user - The user context
1458d69c5d34SMatthew G. Knepley 
1459d69c5d34SMatthew G. Knepley   Output Parameter:
1460934789fcSMatthew G. Knepley . In  - The interpolation matrix
1461d69c5d34SMatthew G. Knepley 
1462d69c5d34SMatthew G. Knepley   Note:
1463d69c5d34SMatthew G. Knepley   The first member of the user context must be an FEMContext.
1464d69c5d34SMatthew G. Knepley 
1465d69c5d34SMatthew G. Knepley   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
1466d69c5d34SMatthew G. Knepley   like a GPU, or vectorize on a multicore machine.
1467d69c5d34SMatthew G. Knepley 
1468d69c5d34SMatthew G. Knepley   Level: developer
1469d69c5d34SMatthew G. Knepley 
1470d69c5d34SMatthew G. Knepley .seealso: DMPlexComputeJacobianFEM()
1471d69c5d34SMatthew G. Knepley @*/
1472934789fcSMatthew G. Knepley PetscErrorCode DMPlexComputeInterpolatorFEM(DM dmc, DM dmf, Mat In, void *user)
1473d69c5d34SMatthew G. Knepley {
1474d69c5d34SMatthew G. Knepley   DM_Plex          *mesh  = (DM_Plex *) dmc->data;
1475d69c5d34SMatthew G. Knepley   PetscFEM         *fem   = (PetscFEM *) user;
1476d69c5d34SMatthew G. Knepley   PetscFE          *fe    = fem->fe;
1477d69c5d34SMatthew G. Knepley   const char       *name  = "Interpolator";
1478d69c5d34SMatthew G. Knepley   PetscFE          *feRef;
1479d69c5d34SMatthew G. Knepley   PetscSection      fsection, fglobalSection;
1480d69c5d34SMatthew G. Knepley   PetscSection      csection, cglobalSection;
1481d69c5d34SMatthew G. Knepley   PetscScalar      *elemMat;
1482942a7a06SMatthew G. Knepley   PetscInt          dim, Nf, f, fieldI, fieldJ, offsetI, offsetJ, cStart, cEnd, c;
1483d69c5d34SMatthew G. Knepley   PetscInt          rCellDof = 0, cCellDof = 0;
1484d69c5d34SMatthew G. Knepley   PetscErrorCode    ierr;
1485d69c5d34SMatthew G. Knepley 
1486d69c5d34SMatthew G. Knepley   PetscFunctionBegin;
1487d69c5d34SMatthew G. Knepley #if 0
1488d69c5d34SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1489d69c5d34SMatthew G. Knepley #endif
1490d69c5d34SMatthew G. Knepley   ierr = DMPlexGetDimension(dmf, &dim);CHKERRQ(ierr);
1491d69c5d34SMatthew G. Knepley   ierr = DMGetDefaultSection(dmf, &fsection);CHKERRQ(ierr);
1492d69c5d34SMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dmf, &fglobalSection);CHKERRQ(ierr);
1493d69c5d34SMatthew G. Knepley   ierr = DMGetDefaultSection(dmc, &csection);CHKERRQ(ierr);
1494d69c5d34SMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dmc, &cglobalSection);CHKERRQ(ierr);
1495d69c5d34SMatthew G. Knepley   ierr = PetscSectionGetNumFields(fsection, &Nf);CHKERRQ(ierr);
1496d69c5d34SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dmc, 0, &cStart, &cEnd);CHKERRQ(ierr);
1497d69c5d34SMatthew G. Knepley   ierr = PetscMalloc1(Nf,&feRef);CHKERRQ(ierr);
1498d69c5d34SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
1499d69c5d34SMatthew G. Knepley     PetscInt rNb, cNb, Nc;
1500d69c5d34SMatthew G. Knepley 
1501d69c5d34SMatthew G. Knepley     ierr = PetscFERefine(fe[f], &feRef[f]);CHKERRQ(ierr);
1502d69c5d34SMatthew G. Knepley     ierr = PetscFEGetDimension(feRef[f], &rNb);CHKERRQ(ierr);
1503d69c5d34SMatthew G. Knepley     ierr = PetscFEGetDimension(fe[f], &cNb);CHKERRQ(ierr);
1504d69c5d34SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe[f], &Nc);CHKERRQ(ierr);
1505d69c5d34SMatthew G. Knepley     rCellDof += rNb*Nc;
1506d69c5d34SMatthew G. Knepley     cCellDof += cNb*Nc;
1507d69c5d34SMatthew G. Knepley   }
1508934789fcSMatthew G. Knepley   ierr = MatZeroEntries(In);CHKERRQ(ierr);
1509d69c5d34SMatthew G. Knepley   ierr = PetscMalloc1(rCellDof*cCellDof,&elemMat);CHKERRQ(ierr);
1510d69c5d34SMatthew G. Knepley   ierr = PetscMemzero(elemMat, rCellDof*cCellDof * sizeof(PetscScalar));CHKERRQ(ierr);
1511d69c5d34SMatthew G. Knepley   for (fieldI = 0, offsetI = 0; fieldI < Nf; ++fieldI) {
1512d69c5d34SMatthew G. Knepley     PetscDualSpace   Qref;
1513d69c5d34SMatthew G. Knepley     PetscQuadrature  f;
1514d69c5d34SMatthew G. Knepley     const PetscReal *qpoints, *qweights;
1515d69c5d34SMatthew G. Knepley     PetscReal       *points;
1516d69c5d34SMatthew G. Knepley     PetscInt         npoints = 0, Nc, Np, fpdim, i, k, p, d;
1517d69c5d34SMatthew G. Knepley 
1518d69c5d34SMatthew G. Knepley     /* Compose points from all dual basis functionals */
1519d69c5d34SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe[fieldI], &Nc);CHKERRQ(ierr);
1520d69c5d34SMatthew G. Knepley     ierr = PetscFEGetDualSpace(feRef[fieldI], &Qref);CHKERRQ(ierr);
1521d69c5d34SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(Qref, &fpdim);CHKERRQ(ierr);
1522d69c5d34SMatthew G. Knepley     for (i = 0; i < fpdim; ++i) {
1523d69c5d34SMatthew G. Knepley       ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1524d69c5d34SMatthew G. Knepley       ierr = PetscQuadratureGetData(f, NULL, &Np, NULL, NULL);CHKERRQ(ierr);
1525d69c5d34SMatthew G. Knepley       npoints += Np;
1526d69c5d34SMatthew G. Knepley     }
1527d69c5d34SMatthew G. Knepley     ierr = PetscMalloc1(npoints*dim,&points);CHKERRQ(ierr);
1528d69c5d34SMatthew G. Knepley     for (i = 0, k = 0; i < fpdim; ++i) {
1529d69c5d34SMatthew G. Knepley       ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1530d69c5d34SMatthew G. Knepley       ierr = PetscQuadratureGetData(f, NULL, &Np, &qpoints, NULL);CHKERRQ(ierr);
1531d69c5d34SMatthew G. Knepley       for (p = 0; p < Np; ++p, ++k) for (d = 0; d < dim; ++d) points[k*dim+d] = qpoints[p*dim+d];
1532d69c5d34SMatthew G. Knepley     }
1533d69c5d34SMatthew G. Knepley 
1534d69c5d34SMatthew G. Knepley     for (fieldJ = 0, offsetJ = 0; fieldJ < Nf; ++fieldJ) {
1535d69c5d34SMatthew G. Knepley       PetscReal *B;
153636a6d9c0SMatthew G. Knepley       PetscInt   NcJ, cpdim, j;
1537d69c5d34SMatthew G. Knepley 
1538d69c5d34SMatthew G. Knepley       /* Evaluate basis at points */
153936a6d9c0SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe[fieldJ], &NcJ);CHKERRQ(ierr);
154036a6d9c0SMatthew 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);
1541d69c5d34SMatthew G. Knepley       ierr = PetscFEGetDimension(fe[fieldJ], &cpdim);CHKERRQ(ierr);
1542ffe73a53SMatthew G. Knepley       /* For now, fields only interpolate themselves */
1543ffe73a53SMatthew G. Knepley       if (fieldI == fieldJ) {
1544d69c5d34SMatthew G. Knepley         ierr = PetscFEGetTabulation(fe[fieldJ], npoints, points, &B, NULL, NULL);CHKERRQ(ierr);
1545d69c5d34SMatthew G. Knepley         for (i = 0, k = 0; i < fpdim; ++i) {
1546d69c5d34SMatthew G. Knepley           ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1547d69c5d34SMatthew G. Knepley           ierr = PetscQuadratureGetData(f, NULL, &Np, NULL, &qweights);CHKERRQ(ierr);
1548d69c5d34SMatthew G. Knepley           for (p = 0; p < Np; ++p, ++k) {
154936a6d9c0SMatthew G. Knepley             for (j = 0; j < cpdim; ++j) {
155036a6d9c0SMatthew 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];
155136a6d9c0SMatthew G. Knepley             }
1552d69c5d34SMatthew G. Knepley           }
1553d69c5d34SMatthew G. Knepley         }
1554d69c5d34SMatthew G. Knepley         ierr = PetscFERestoreTabulation(fe[fieldJ], npoints, points, &B, NULL, NULL);CHKERRQ(ierr);CHKERRQ(ierr);
1555ffe73a53SMatthew G. Knepley       }
155636a6d9c0SMatthew G. Knepley       offsetJ += cpdim*NcJ;
1557d69c5d34SMatthew G. Knepley     }
1558d69c5d34SMatthew G. Knepley     offsetI += fpdim*Nc;
1559549a8adaSMatthew G. Knepley     ierr = PetscFree(points);CHKERRQ(ierr);
1560d69c5d34SMatthew G. Knepley   }
1561d69c5d34SMatthew G. Knepley   if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(0, name, rCellDof, cCellDof, elemMat);CHKERRQ(ierr);}
1562d69c5d34SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
1563934789fcSMatthew G. Knepley     ierr = DMPlexMatSetClosureRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, In, c, elemMat, INSERT_VALUES);CHKERRQ(ierr);
1564d69c5d34SMatthew G. Knepley   }
1565549a8adaSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {ierr = PetscFEDestroy(&feRef[f]);CHKERRQ(ierr);}
1566d69c5d34SMatthew G. Knepley   ierr = PetscFree(feRef);CHKERRQ(ierr);
1567549a8adaSMatthew G. Knepley   ierr = PetscFree(elemMat);CHKERRQ(ierr);
1568934789fcSMatthew G. Knepley   ierr = MatAssemblyBegin(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1569934789fcSMatthew G. Knepley   ierr = MatAssemblyEnd(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1570d69c5d34SMatthew G. Knepley   if (mesh->printFEM) {
1571d69c5d34SMatthew G. Knepley     ierr = PetscPrintf(PETSC_COMM_WORLD, "%s:\n", name);CHKERRQ(ierr);
1572934789fcSMatthew G. Knepley     ierr = MatChop(In, 1.0e-10);CHKERRQ(ierr);
1573934789fcSMatthew G. Knepley     ierr = MatView(In, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1574d69c5d34SMatthew G. Knepley   }
1575d69c5d34SMatthew G. Knepley #if 0
1576d69c5d34SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1577d69c5d34SMatthew G. Knepley #endif
1578d69c5d34SMatthew G. Knepley   PetscFunctionReturn(0);
1579d69c5d34SMatthew G. Knepley }
15806c73c22cSMatthew G. Knepley 
15816c73c22cSMatthew G. Knepley #undef __FUNCT__
15826c73c22cSMatthew G. Knepley #define __FUNCT__ "DMPlexAddBoundary"
15836c73c22cSMatthew G. Knepley /* The ids can be overridden by the command line option -bc_<boundary name> */
15846c73c22cSMatthew G. Knepley PetscErrorCode DMPlexAddBoundary(DM dm, PetscBool isEssential, const char name[], PetscInt field, void (*bcFunc)(), PetscInt numids, const PetscInt *ids, void *ctx)
15856c73c22cSMatthew G. Knepley {
15866c73c22cSMatthew G. Knepley   DM_Plex       *mesh = (DM_Plex *) dm->data;
15876c73c22cSMatthew G. Knepley   DMBoundary     b;
15886c73c22cSMatthew G. Knepley   PetscErrorCode ierr;
15896c73c22cSMatthew G. Knepley 
15906c73c22cSMatthew G. Knepley   PetscFunctionBegin;
15916c73c22cSMatthew G. Knepley   ierr = PetscNew(&b);CHKERRQ(ierr);
15926c73c22cSMatthew G. Knepley   ierr = PetscStrallocpy(name, (char **) &b->name);CHKERRQ(ierr);
15936c73c22cSMatthew G. Knepley   ierr = PetscMalloc1(numids, &b->ids);CHKERRQ(ierr);
15946c73c22cSMatthew G. Knepley   ierr = PetscMemcpy(b->ids, ids, numids*sizeof(PetscInt));CHKERRQ(ierr);
15956c73c22cSMatthew G. Knepley   b->essential   = isEssential;
15966c73c22cSMatthew G. Knepley   b->field       = field;
15976c73c22cSMatthew G. Knepley   b->func        = bcFunc;
15986c73c22cSMatthew G. Knepley   b->numids      = numids;
15996c73c22cSMatthew G. Knepley   b->ctx         = ctx;
16006c73c22cSMatthew G. Knepley   b->next        = mesh->boundary;
16016c73c22cSMatthew G. Knepley   mesh->boundary = b;
16026c73c22cSMatthew G. Knepley   PetscFunctionReturn(0);
16036c73c22cSMatthew G. Knepley }
16046c73c22cSMatthew G. Knepley 
16056c73c22cSMatthew G. Knepley #undef __FUNCT__
16066c73c22cSMatthew G. Knepley #define __FUNCT__ "DMPlexGetNumBoundary"
16076c73c22cSMatthew G. Knepley PetscErrorCode DMPlexGetNumBoundary(DM dm, PetscInt *numBd)
16086c73c22cSMatthew G. Knepley {
16096c73c22cSMatthew G. Knepley   DM_Plex   *mesh = (DM_Plex *) dm->data;
16106c73c22cSMatthew G. Knepley   DMBoundary b    = mesh->boundary;
16116c73c22cSMatthew G. Knepley 
16126c73c22cSMatthew G. Knepley   PetscFunctionBegin;
16136c73c22cSMatthew G. Knepley   *numBd = 0;
16146c73c22cSMatthew G. Knepley   while (b) {++(*numBd); b = b->next;}
16156c73c22cSMatthew G. Knepley   PetscFunctionReturn(0);
16166c73c22cSMatthew G. Knepley }
16176c73c22cSMatthew G. Knepley 
16186c73c22cSMatthew G. Knepley #undef __FUNCT__
16196c73c22cSMatthew G. Knepley #define __FUNCT__ "DMPlexGetBoundary"
16206c73c22cSMatthew G. Knepley PetscErrorCode DMPlexGetBoundary(DM dm, PetscInt bd, PetscBool *isEssential, const char **name, PetscInt *field, void (**func)(), PetscInt *numids, const PetscInt **ids, void **ctx)
16216c73c22cSMatthew G. Knepley {
16226c73c22cSMatthew G. Knepley   DM_Plex   *mesh = (DM_Plex *) dm->data;
16236c73c22cSMatthew G. Knepley   DMBoundary b    = mesh->boundary;
16246c73c22cSMatthew G. Knepley   PetscInt   n    = 0;
16256c73c22cSMatthew G. Knepley 
16266c73c22cSMatthew G. Knepley   PetscFunctionBegin;
16276c73c22cSMatthew G. Knepley   while (b) {
16286c73c22cSMatthew G. Knepley     if (n == bd) break;
16296c73c22cSMatthew G. Knepley     b = b->next;
16306c73c22cSMatthew G. Knepley     ++n;
16316c73c22cSMatthew G. Knepley   }
16326c73c22cSMatthew G. Knepley   if (n != bd) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Boundary %d is not in [0, %d)", bd, n);
16336c73c22cSMatthew G. Knepley   if (isEssential) {
16346c73c22cSMatthew G. Knepley     PetscValidPointer(isEssential, 3);
16356c73c22cSMatthew G. Knepley     *isEssential = b->essential;
16366c73c22cSMatthew G. Knepley   }
16376c73c22cSMatthew G. Knepley   if (name) {
16386c73c22cSMatthew G. Knepley     PetscValidPointer(name, 4);
16396c73c22cSMatthew G. Knepley     *name = b->name;
16406c73c22cSMatthew G. Knepley   }
16416c73c22cSMatthew G. Knepley   if (field) {
16426c73c22cSMatthew G. Knepley     PetscValidPointer(field, 5);
16436c73c22cSMatthew G. Knepley     *field = b->field;
16446c73c22cSMatthew G. Knepley   }
16456c73c22cSMatthew G. Knepley   if (func) {
16466c73c22cSMatthew G. Knepley     PetscValidPointer(func, 6);
16476c73c22cSMatthew G. Knepley     *func = b->func;
16486c73c22cSMatthew G. Knepley   }
16496c73c22cSMatthew G. Knepley   if (numids) {
16506c73c22cSMatthew G. Knepley     PetscValidPointer(numids, 7);
16516c73c22cSMatthew G. Knepley     *numids = b->numids;
16526c73c22cSMatthew G. Knepley   }
16536c73c22cSMatthew G. Knepley   if (ids) {
16546c73c22cSMatthew G. Knepley     PetscValidPointer(ids, 8);
16556c73c22cSMatthew G. Knepley     *ids = b->ids;
16566c73c22cSMatthew G. Knepley   }
16576c73c22cSMatthew G. Knepley   if (ctx) {
16586c73c22cSMatthew G. Knepley     PetscValidPointer(ctx, 9);
16596c73c22cSMatthew G. Knepley     *ctx = b->ctx;
16606c73c22cSMatthew G. Knepley   }
16616c73c22cSMatthew G. Knepley   PetscFunctionReturn(0);
16626c73c22cSMatthew G. Knepley }
1663