xref: /petsc/src/dm/impls/plex/plexfem.c (revision af1eca97be555adbd75c517f986ed82780294fad)
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__
184cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexProjectFunctionLocal"
185c110b1eeSGeoffrey Irving PetscErrorCode DMPlexProjectFunctionLocal(DM dm, PetscFE fe[], void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX)
186cb1e1211SMatthew G Knepley {
18772f94c41SMatthew G. Knepley   PetscDualSpace *sp;
18872f94c41SMatthew G. Knepley   PetscSection    section;
18972f94c41SMatthew G. Knepley   PetscScalar    *values;
19072f94c41SMatthew G. Knepley   PetscReal      *v0, *J, detJ;
191120386c5SMatthew G. Knepley   PetscInt        numFields, numComp, dim, spDim, totDim = 0, numValues, cStart, cEnd, c, f, d, v, comp;
192cb1e1211SMatthew G Knepley   PetscErrorCode  ierr;
193cb1e1211SMatthew G Knepley 
194cb1e1211SMatthew G Knepley   PetscFunctionBegin;
195cb1e1211SMatthew G Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
19672f94c41SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
197785e854fSJed Brown   ierr = PetscMalloc1(numFields, &sp);CHKERRQ(ierr);
19872f94c41SMatthew G. Knepley   for (f = 0; f < numFields; ++f) {
19972f94c41SMatthew G. Knepley     ierr = PetscFEGetDualSpace(fe[f], &sp[f]);CHKERRQ(ierr);
20072f94c41SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe[f], &numComp);CHKERRQ(ierr);
20172f94c41SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
20272f94c41SMatthew G. Knepley     totDim += spDim*numComp;
203cb1e1211SMatthew G Knepley   }
20472f94c41SMatthew G. Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
20572f94c41SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
20672f94c41SMatthew G. Knepley   ierr = DMPlexVecGetClosure(dm, section, localX, cStart, &numValues, NULL);CHKERRQ(ierr);
20772f94c41SMatthew 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);
20872f94c41SMatthew G. Knepley   ierr = DMGetWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
209dcca6d9dSJed Brown   ierr = PetscMalloc2(dim,&v0,dim*dim,&J);CHKERRQ(ierr);
21072f94c41SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
21172f94c41SMatthew G. Knepley     PetscCellGeometry geom;
212cb1e1211SMatthew G Knepley 
213cb1e1211SMatthew G Knepley     ierr = DMPlexComputeCellGeometry(dm, c, v0, J, NULL, &detJ);CHKERRQ(ierr);
21472f94c41SMatthew G. Knepley     geom.v0   = v0;
21572f94c41SMatthew G. Knepley     geom.J    = J;
21672f94c41SMatthew G. Knepley     geom.detJ = &detJ;
21772f94c41SMatthew G. Knepley     for (f = 0, v = 0; f < numFields; ++f) {
218c110b1eeSGeoffrey Irving       void * const ctx = ctxs ? ctxs[f] : NULL;
21972f94c41SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe[f], &numComp);CHKERRQ(ierr);
22072f94c41SMatthew G. Knepley       ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
22172f94c41SMatthew G. Knepley       for (d = 0; d < spDim; ++d) {
222120386c5SMatthew G. Knepley         if (funcs[f]) {
223c110b1eeSGeoffrey Irving           ierr = PetscDualSpaceApply(sp[f], d, geom, numComp, funcs[f], ctx, &values[v]);CHKERRQ(ierr);
224120386c5SMatthew G. Knepley         } else {
225120386c5SMatthew G. Knepley           for (comp = 0; comp < numComp; ++comp) values[v+comp] = 0.0;
226120386c5SMatthew G. Knepley         }
22772f94c41SMatthew G. Knepley         v += numComp;
228cb1e1211SMatthew G Knepley       }
229cb1e1211SMatthew G Knepley     }
23072f94c41SMatthew G. Knepley     ierr = DMPlexVecSetClosure(dm, section, localX, c, values, mode);CHKERRQ(ierr);
231cb1e1211SMatthew G Knepley   }
23272f94c41SMatthew G. Knepley   ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
2331f2da991SMatthew G. Knepley   ierr = PetscFree2(v0,J);CHKERRQ(ierr);
23472f94c41SMatthew G. Knepley   ierr = PetscFree(sp);CHKERRQ(ierr);
235cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
236cb1e1211SMatthew G Knepley }
237cb1e1211SMatthew G Knepley 
238cb1e1211SMatthew G Knepley #undef __FUNCT__
239cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexProjectFunction"
240cb1e1211SMatthew G Knepley /*@C
241cb1e1211SMatthew G Knepley   DMPlexProjectFunction - This projects the given function into the function space provided.
242cb1e1211SMatthew G Knepley 
243cb1e1211SMatthew G Knepley   Input Parameters:
244cb1e1211SMatthew G Knepley + dm      - The DM
24572f94c41SMatthew G. Knepley . fe      - The PetscFE associated with the field
24672f94c41SMatthew G. Knepley . funcs   - The coordinate functions to evaluate, one per field
247c110b1eeSGeoffrey Irving . ctxs    - Optional array of contexts to pass to each coordinate function.  ctxs itself may be null.
248cb1e1211SMatthew G Knepley - mode    - The insertion mode for values
249cb1e1211SMatthew G Knepley 
250cb1e1211SMatthew G Knepley   Output Parameter:
251cb1e1211SMatthew G Knepley . X - vector
252cb1e1211SMatthew G Knepley 
253cb1e1211SMatthew G Knepley   Level: developer
254cb1e1211SMatthew G Knepley 
255878cb397SSatish Balay .seealso: DMPlexComputeL2Diff()
256878cb397SSatish Balay @*/
257c110b1eeSGeoffrey Irving PetscErrorCode DMPlexProjectFunction(DM dm, PetscFE fe[], void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, InsertMode mode, Vec X)
258cb1e1211SMatthew G Knepley {
259cb1e1211SMatthew G Knepley   Vec            localX;
260cb1e1211SMatthew G Knepley   PetscErrorCode ierr;
261cb1e1211SMatthew G Knepley 
262cb1e1211SMatthew G Knepley   PetscFunctionBegin;
2639a800dd8SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
264cb1e1211SMatthew G Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
265c110b1eeSGeoffrey Irving   ierr = DMPlexProjectFunctionLocal(dm, fe, funcs, ctxs, mode, localX);CHKERRQ(ierr);
266cb1e1211SMatthew G Knepley   ierr = DMLocalToGlobalBegin(dm, localX, mode, X);CHKERRQ(ierr);
267cb1e1211SMatthew G Knepley   ierr = DMLocalToGlobalEnd(dm, localX, mode, X);CHKERRQ(ierr);
268cb1e1211SMatthew G Knepley   ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
269cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
270cb1e1211SMatthew G Knepley }
271cb1e1211SMatthew G Knepley 
272cb1e1211SMatthew G Knepley #undef __FUNCT__
273cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexComputeL2Diff"
274cb1e1211SMatthew G Knepley /*@C
275cb1e1211SMatthew G Knepley   DMPlexComputeL2Diff - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h.
276cb1e1211SMatthew G Knepley 
277cb1e1211SMatthew G Knepley   Input Parameters:
278cb1e1211SMatthew G Knepley + dm    - The DM
279c5bbbd5bSMatthew G. Knepley . fe    - The PetscFE object for each field
280cb1e1211SMatthew G Knepley . funcs - The functions to evaluate for each field component
28151259fa3SMatthew G. Knepley . ctxs  - Optional array of contexts to pass to each function, or NULL.
282cb1e1211SMatthew G Knepley - X     - The coefficient vector u_h
283cb1e1211SMatthew G Knepley 
284cb1e1211SMatthew G Knepley   Output Parameter:
285cb1e1211SMatthew G Knepley . diff - The diff ||u - u_h||_2
286cb1e1211SMatthew G Knepley 
287cb1e1211SMatthew G Knepley   Level: developer
288cb1e1211SMatthew G Knepley 
28923d86601SMatthew G. Knepley .seealso: DMPlexProjectFunction(), DMPlexComputeL2GradientDiff()
290878cb397SSatish Balay @*/
291c110b1eeSGeoffrey Irving PetscErrorCode DMPlexComputeL2Diff(DM dm, PetscFE fe[], void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, PetscReal *diff)
292cb1e1211SMatthew G Knepley {
293cb1e1211SMatthew G Knepley   const PetscInt  debug = 0;
294cb1e1211SMatthew G Knepley   PetscSection    section;
295c5bbbd5bSMatthew G. Knepley   PetscQuadrature quad;
296cb1e1211SMatthew G Knepley   Vec             localX;
29772f94c41SMatthew G. Knepley   PetscScalar    *funcVal;
298cb1e1211SMatthew G Knepley   PetscReal      *coords, *v0, *J, *invJ, detJ;
299cb1e1211SMatthew G Knepley   PetscReal       localDiff = 0.0;
300cb1e1211SMatthew G Knepley   PetscInt        dim, numFields, numComponents = 0, cStart, cEnd, c, field, fieldOffset, comp;
301cb1e1211SMatthew G Knepley   PetscErrorCode  ierr;
302cb1e1211SMatthew G Knepley 
303cb1e1211SMatthew G Knepley   PetscFunctionBegin;
304cb1e1211SMatthew G Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
305cb1e1211SMatthew G Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
306cb1e1211SMatthew G Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
307cb1e1211SMatthew G Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
308cb1e1211SMatthew G Knepley   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
309cb1e1211SMatthew G Knepley   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
310cb1e1211SMatthew G Knepley   for (field = 0; field < numFields; ++field) {
311c5bbbd5bSMatthew G. Knepley     PetscInt Nc;
312c5bbbd5bSMatthew G. Knepley 
313c5bbbd5bSMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe[field], &Nc);CHKERRQ(ierr);
314c5bbbd5bSMatthew G. Knepley     numComponents += Nc;
315cb1e1211SMatthew G Knepley   }
316c110b1eeSGeoffrey Irving   ierr = DMPlexProjectFunctionLocal(dm, fe, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr);
317dcca6d9dSJed Brown   ierr = PetscMalloc5(numComponents,&funcVal,dim,&coords,dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr);
318cb1e1211SMatthew G Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
319c5bbbd5bSMatthew G. Knepley   ierr = PetscFEGetQuadrature(fe[0], &quad);CHKERRQ(ierr);
320cb1e1211SMatthew G Knepley   for (c = cStart; c < cEnd; ++c) {
321a1e44745SMatthew G. Knepley     PetscScalar *x = NULL;
322cb1e1211SMatthew G Knepley     PetscReal    elemDiff = 0.0;
323cb1e1211SMatthew G Knepley 
324cb1e1211SMatthew G Knepley     ierr = DMPlexComputeCellGeometry(dm, c, v0, J, invJ, &detJ);CHKERRQ(ierr);
325cb1e1211SMatthew G Knepley     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
326cb1e1211SMatthew G Knepley     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
327cb1e1211SMatthew G Knepley 
328cb1e1211SMatthew G Knepley     for (field = 0, comp = 0, fieldOffset = 0; field < numFields; ++field) {
329c110b1eeSGeoffrey Irving       void * const     ctx = ctxs ? ctxs[field] : NULL;
33021454ff5SMatthew G. Knepley       const PetscReal *quadPoints, *quadWeights;
331c5bbbd5bSMatthew G. Knepley       PetscReal       *basis;
33221454ff5SMatthew G. Knepley       PetscInt         numQuadPoints, numBasisFuncs, numBasisComps, q, d, e, fc, f;
333cb1e1211SMatthew G Knepley 
33421454ff5SMatthew G. Knepley       ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr);
335c5bbbd5bSMatthew G. Knepley       ierr = PetscFEGetDimension(fe[field], &numBasisFuncs);CHKERRQ(ierr);
336c5bbbd5bSMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe[field], &numBasisComps);CHKERRQ(ierr);
337c5bbbd5bSMatthew G. Knepley       ierr = PetscFEGetDefaultTabulation(fe[field], &basis, NULL, NULL);CHKERRQ(ierr);
338cb1e1211SMatthew G Knepley       if (debug) {
339cb1e1211SMatthew G Knepley         char title[1024];
340cb1e1211SMatthew G Knepley         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
341cb1e1211SMatthew G Knepley         ierr = DMPrintCellVector(c, title, numBasisFuncs*numBasisComps, &x[fieldOffset]);CHKERRQ(ierr);
342cb1e1211SMatthew G Knepley       }
343cb1e1211SMatthew G Knepley       for (q = 0; q < numQuadPoints; ++q) {
344cb1e1211SMatthew G Knepley         for (d = 0; d < dim; d++) {
345cb1e1211SMatthew G Knepley           coords[d] = v0[d];
346cb1e1211SMatthew G Knepley           for (e = 0; e < dim; e++) {
347cb1e1211SMatthew G Knepley             coords[d] += J[d*dim+e]*(quadPoints[q*dim+e] + 1.0);
348cb1e1211SMatthew G Knepley           }
349cb1e1211SMatthew G Knepley         }
350c110b1eeSGeoffrey Irving         (*funcs[field])(coords, funcVal, ctx);
351cb1e1211SMatthew G Knepley         for (fc = 0; fc < numBasisComps; ++fc) {
352a1d24da5SMatthew G. Knepley           PetscScalar interpolant = 0.0;
353a1d24da5SMatthew G. Knepley 
354cb1e1211SMatthew G Knepley           for (f = 0; f < numBasisFuncs; ++f) {
355cb1e1211SMatthew G Knepley             const PetscInt fidx = f*numBasisComps+fc;
356a1d24da5SMatthew G. Knepley             interpolant += x[fieldOffset+fidx]*basis[q*numBasisFuncs*numBasisComps+fidx];
357cb1e1211SMatthew G Knepley           }
35872f94c41SMatthew 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);}
35972f94c41SMatthew G. Knepley           elemDiff += PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ;
360cb1e1211SMatthew G Knepley         }
361cb1e1211SMatthew G Knepley       }
362cb1e1211SMatthew G Knepley       comp        += numBasisComps;
363cb1e1211SMatthew G Knepley       fieldOffset += numBasisFuncs*numBasisComps;
364cb1e1211SMatthew G Knepley     }
365cb1e1211SMatthew G Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
366cb1e1211SMatthew G Knepley     if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);}
367cb1e1211SMatthew G Knepley     localDiff += elemDiff;
368cb1e1211SMatthew G Knepley   }
36972f94c41SMatthew G. Knepley   ierr  = PetscFree5(funcVal,coords,v0,J,invJ);CHKERRQ(ierr);
370cb1e1211SMatthew G Knepley   ierr  = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
37186a74ee0SMatthew G. Knepley   ierr  = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
372cb1e1211SMatthew G Knepley   *diff = PetscSqrtReal(*diff);
373cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
374cb1e1211SMatthew G Knepley }
375cb1e1211SMatthew G Knepley 
376cb1e1211SMatthew G Knepley #undef __FUNCT__
37740e14135SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeL2GradientDiff"
37840e14135SMatthew G. Knepley /*@C
37940e14135SMatthew 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.
38040e14135SMatthew G. Knepley 
38140e14135SMatthew G. Knepley   Input Parameters:
38240e14135SMatthew G. Knepley + dm    - The DM
38340e14135SMatthew G. Knepley . fe    - The PetscFE object for each field
38440e14135SMatthew G. Knepley . funcs - The gradient functions to evaluate for each field component
38551259fa3SMatthew G. Knepley . ctxs  - Optional array of contexts to pass to each function, or NULL.
38640e14135SMatthew G. Knepley . X     - The coefficient vector u_h
38740e14135SMatthew G. Knepley - n     - The vector to project along
38840e14135SMatthew G. Knepley 
38940e14135SMatthew G. Knepley   Output Parameter:
39040e14135SMatthew G. Knepley . diff - The diff ||(grad u - grad u_h) . n||_2
39140e14135SMatthew G. Knepley 
39240e14135SMatthew G. Knepley   Level: developer
39340e14135SMatthew G. Knepley 
39440e14135SMatthew G. Knepley .seealso: DMPlexProjectFunction(), DMPlexComputeL2Diff()
39540e14135SMatthew G. Knepley @*/
39651259fa3SMatthew 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)
397cb1e1211SMatthew G Knepley {
39840e14135SMatthew G. Knepley   const PetscInt  debug = 0;
399cb1e1211SMatthew G Knepley   PetscSection    section;
40040e14135SMatthew G. Knepley   PetscQuadrature quad;
40140e14135SMatthew G. Knepley   Vec             localX;
40240e14135SMatthew G. Knepley   PetscScalar    *funcVal, *interpolantVec;
40340e14135SMatthew G. Knepley   PetscReal      *coords, *realSpaceDer, *v0, *J, *invJ, detJ;
40440e14135SMatthew G. Knepley   PetscReal       localDiff = 0.0;
40540e14135SMatthew G. Knepley   PetscInt        dim, numFields, numComponents = 0, cStart, cEnd, c, field, fieldOffset, comp;
406cb1e1211SMatthew G Knepley   PetscErrorCode  ierr;
407cb1e1211SMatthew G Knepley 
408cb1e1211SMatthew G Knepley   PetscFunctionBegin;
40940e14135SMatthew G. Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
41040e14135SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
41140e14135SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
41240e14135SMatthew G. Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
41340e14135SMatthew G. Knepley   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
41440e14135SMatthew G. Knepley   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
415652b88e8SMatthew G. Knepley   for (field = 0; field < numFields; ++field) {
41640e14135SMatthew G. Knepley     PetscInt Nc;
417652b88e8SMatthew G. Knepley 
41840e14135SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe[field], &Nc);CHKERRQ(ierr);
41940e14135SMatthew G. Knepley     numComponents += Nc;
420652b88e8SMatthew G. Knepley   }
42140e14135SMatthew G. Knepley   /* ierr = DMPlexProjectFunctionLocal(dm, fe, funcs, INSERT_BC_VALUES, localX);CHKERRQ(ierr); */
42240e14135SMatthew G. Knepley   ierr = PetscMalloc7(numComponents,&funcVal,dim,&coords,dim,&realSpaceDer,dim,&v0,dim*dim,&J,dim*dim,&invJ,dim,&interpolantVec);CHKERRQ(ierr);
42340e14135SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
42440e14135SMatthew G. Knepley   ierr = PetscFEGetQuadrature(fe[0], &quad);CHKERRQ(ierr);
42540e14135SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
42640e14135SMatthew G. Knepley     PetscScalar *x = NULL;
42740e14135SMatthew G. Knepley     PetscReal    elemDiff = 0.0;
428652b88e8SMatthew G. Knepley 
42940e14135SMatthew G. Knepley     ierr = DMPlexComputeCellGeometry(dm, c, v0, J, invJ, &detJ);CHKERRQ(ierr);
43040e14135SMatthew G. Knepley     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
43140e14135SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
43240e14135SMatthew G. Knepley 
43340e14135SMatthew G. Knepley     for (field = 0, comp = 0, fieldOffset = 0; field < numFields; ++field) {
43451259fa3SMatthew G. Knepley       void * const     ctx = ctxs ? ctxs[field] : NULL;
43521454ff5SMatthew G. Knepley       const PetscReal *quadPoints, *quadWeights;
43640e14135SMatthew G. Knepley       PetscReal       *basisDer;
43721454ff5SMatthew G. Knepley       PetscInt         numQuadPoints, Nb, Ncomp, q, d, e, fc, f, g;
43840e14135SMatthew G. Knepley 
43921454ff5SMatthew G. Knepley       ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr);
44040e14135SMatthew G. Knepley       ierr = PetscFEGetDimension(fe[field], &Nb);CHKERRQ(ierr);
44140e14135SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe[field], &Ncomp);CHKERRQ(ierr);
44240e14135SMatthew G. Knepley       ierr = PetscFEGetDefaultTabulation(fe[field], NULL, &basisDer, NULL);CHKERRQ(ierr);
44340e14135SMatthew G. Knepley       if (debug) {
44440e14135SMatthew G. Knepley         char title[1024];
44540e14135SMatthew G. Knepley         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
44640e14135SMatthew G. Knepley         ierr = DMPrintCellVector(c, title, Nb*Ncomp, &x[fieldOffset]);CHKERRQ(ierr);
447652b88e8SMatthew G. Knepley       }
44840e14135SMatthew G. Knepley       for (q = 0; q < numQuadPoints; ++q) {
44940e14135SMatthew G. Knepley         for (d = 0; d < dim; d++) {
45040e14135SMatthew G. Knepley           coords[d] = v0[d];
45140e14135SMatthew G. Knepley           for (e = 0; e < dim; e++) {
45240e14135SMatthew G. Knepley             coords[d] += J[d*dim+e]*(quadPoints[q*dim+e] + 1.0);
453652b88e8SMatthew G. Knepley           }
45440e14135SMatthew G. Knepley         }
45551259fa3SMatthew G. Knepley         (*funcs[field])(coords, n, funcVal, ctx);
45640e14135SMatthew G. Knepley         for (fc = 0; fc < Ncomp; ++fc) {
45740e14135SMatthew G. Knepley           PetscScalar interpolant = 0.0;
45840e14135SMatthew G. Knepley 
45940e14135SMatthew G. Knepley           for (d = 0; d < dim; ++d) interpolantVec[d] = 0.0;
46040e14135SMatthew G. Knepley           for (f = 0; f < Nb; ++f) {
46140e14135SMatthew G. Knepley             const PetscInt fidx = f*Ncomp+fc;
46240e14135SMatthew G. Knepley 
46340e14135SMatthew G. Knepley             for (d = 0; d < dim; ++d) {
46440e14135SMatthew G. Knepley               realSpaceDer[d] = 0.0;
46540e14135SMatthew G. Knepley               for (g = 0; g < dim; ++g) {
46640e14135SMatthew G. Knepley                 realSpaceDer[d] += invJ[g*dim+d]*basisDer[(q*Nb*Ncomp+fidx)*dim+g];
46740e14135SMatthew G. Knepley               }
46840e14135SMatthew G. Knepley               interpolantVec[d] += x[fieldOffset+fidx]*realSpaceDer[d];
46940e14135SMatthew G. Knepley             }
47040e14135SMatthew G. Knepley           }
47140e14135SMatthew G. Knepley           for (d = 0; d < dim; ++d) interpolant += interpolantVec[d]*n[d];
47240e14135SMatthew 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);}
47340e14135SMatthew G. Knepley           elemDiff += PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ;
47440e14135SMatthew G. Knepley         }
47540e14135SMatthew G. Knepley       }
47640e14135SMatthew G. Knepley       comp        += Ncomp;
47740e14135SMatthew G. Knepley       fieldOffset += Nb*Ncomp;
47840e14135SMatthew G. Knepley     }
47940e14135SMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
48040e14135SMatthew G. Knepley     if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);}
48140e14135SMatthew G. Knepley     localDiff += elemDiff;
48240e14135SMatthew G. Knepley   }
48340e14135SMatthew G. Knepley   ierr  = PetscFree7(funcVal,coords,realSpaceDer,v0,J,invJ,interpolantVec);CHKERRQ(ierr);
48440e14135SMatthew G. Knepley   ierr  = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
48540e14135SMatthew G. Knepley   ierr  = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
48640e14135SMatthew G. Knepley   *diff = PetscSqrtReal(*diff);
487cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
488cb1e1211SMatthew G Knepley }
489cb1e1211SMatthew G Knepley 
490a0845e3aSMatthew G. Knepley #undef __FUNCT__
491a0845e3aSMatthew G. Knepley #define __FUNCT__ "DMPlexComputeResidualFEM"
492a0845e3aSMatthew G. Knepley /*@
493a0845e3aSMatthew G. Knepley   DMPlexComputeResidualFEM - Form the local residual F from the local input X using pointwise functions specified by the user
494a0845e3aSMatthew G. Knepley 
495a0845e3aSMatthew G. Knepley   Input Parameters:
496a0845e3aSMatthew G. Knepley + dm - The mesh
497a0845e3aSMatthew G. Knepley . X  - Local input vector
498a0845e3aSMatthew G. Knepley - user - The user context
499a0845e3aSMatthew G. Knepley 
500a0845e3aSMatthew G. Knepley   Output Parameter:
501a0845e3aSMatthew G. Knepley . F  - Local output vector
502a0845e3aSMatthew G. Knepley 
503a0845e3aSMatthew G. Knepley   Note:
5048026896cSMatthew G. Knepley   The first member of the user context must be an FEMContext.
505a0845e3aSMatthew G. Knepley 
506a0845e3aSMatthew G. Knepley   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
507a0845e3aSMatthew G. Knepley   like a GPU, or vectorize on a multicore machine.
508a0845e3aSMatthew G. Knepley 
509a0845e3aSMatthew G. Knepley   Level: developer
510a0845e3aSMatthew G. Knepley 
511a0845e3aSMatthew G. Knepley .seealso: DMPlexComputeJacobianActionFEM()
512a0845e3aSMatthew G. Knepley @*/
513a0845e3aSMatthew G. Knepley PetscErrorCode DMPlexComputeResidualFEM(DM dm, Vec X, Vec F, void *user)
514a0845e3aSMatthew G. Knepley {
515a0845e3aSMatthew G. Knepley   DM_Plex          *mesh  = (DM_Plex *) dm->data;
5169a559087SMatthew G. Knepley   PetscFEM         *fem   = (PetscFEM *) user;
517a0845e3aSMatthew G. Knepley   PetscFE          *fe    = fem->fe;
5189a559087SMatthew G. Knepley   PetscFE          *feAux = fem->feAux;
519f1ea0e2fSMatthew G. Knepley   PetscFE          *feBd  = fem->feBd;
520a0845e3aSMatthew G. Knepley   const char       *name  = "Residual";
5219a559087SMatthew G. Knepley   DM                dmAux;
5229a559087SMatthew G. Knepley   Vec               A;
523a0845e3aSMatthew G. Knepley   PetscQuadrature   q;
524a0845e3aSMatthew G. Knepley   PetscCellGeometry geom;
5259a559087SMatthew G. Knepley   PetscSection      section, sectionAux;
526a0845e3aSMatthew G. Knepley   PetscReal        *v0, *J, *invJ, *detJ;
52701599b20SMatthew G. Knepley   PetscScalar      *elemVec, *u, *a = NULL;
5289a559087SMatthew G. Knepley   PetscInt          dim, Nf, NfAux = 0, f, numCells, cStart, cEnd, c;
529a0845e3aSMatthew G. Knepley   PetscInt          cellDof = 0, numComponents = 0;
5309a559087SMatthew G. Knepley   PetscInt          cellDofAux = 0, numComponentsAux = 0;
531a0845e3aSMatthew G. Knepley   PetscErrorCode    ierr;
532a0845e3aSMatthew G. Knepley 
533a0845e3aSMatthew G. Knepley   PetscFunctionBegin;
534a0845e3aSMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_ResidualFEM,dm,0,0,0);CHKERRQ(ierr);
535a0845e3aSMatthew G. Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
536a0845e3aSMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
5379a559087SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
538a0845e3aSMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
539a0845e3aSMatthew G. Knepley   numCells = cEnd - cStart;
5409a559087SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
541a0845e3aSMatthew G. Knepley     PetscInt Nb, Nc;
542a0845e3aSMatthew G. Knepley 
543a0845e3aSMatthew G. Knepley     ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr);
544a0845e3aSMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe[f], &Nc);CHKERRQ(ierr);
545a0845e3aSMatthew G. Knepley     cellDof       += Nb*Nc;
546a0845e3aSMatthew G. Knepley     numComponents += Nc;
547a0845e3aSMatthew G. Knepley   }
5489a559087SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
5499a559087SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr);
5509a559087SMatthew G. Knepley   if (dmAux) {
5519a559087SMatthew G. Knepley     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
5529a559087SMatthew G. Knepley     ierr = PetscSectionGetNumFields(sectionAux, &NfAux);CHKERRQ(ierr);
5539a559087SMatthew G. Knepley   }
5549a559087SMatthew G. Knepley   for (f = 0; f < NfAux; ++f) {
5559a559087SMatthew G. Knepley     PetscInt Nb, Nc;
5569a559087SMatthew G. Knepley 
5579a559087SMatthew G. Knepley     ierr = PetscFEGetDimension(feAux[f], &Nb);CHKERRQ(ierr);
5589a559087SMatthew G. Knepley     ierr = PetscFEGetNumComponents(feAux[f], &Nc);CHKERRQ(ierr);
5599a559087SMatthew G. Knepley     cellDofAux       += Nb*Nc;
5609a559087SMatthew G. Knepley     numComponentsAux += Nc;
5619a559087SMatthew G. Knepley   }
562c110b1eeSGeoffrey Irving   ierr = DMPlexProjectFunctionLocal(dm, fe, fem->bcFuncs, fem->bcCtxs, INSERT_BC_VALUES, X);CHKERRQ(ierr);
563a0845e3aSMatthew G. Knepley   ierr = VecSet(F, 0.0);CHKERRQ(ierr);
564dcca6d9dSJed Brown   ierr = PetscMalloc6(numCells*cellDof,&u,numCells*dim,&v0,numCells*dim*dim,&J,numCells*dim*dim,&invJ,numCells,&detJ,numCells*cellDof,&elemVec);CHKERRQ(ierr);
565785e854fSJed Brown   if (dmAux) {ierr = PetscMalloc1(numCells*cellDofAux, &a);CHKERRQ(ierr);}
566a0845e3aSMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
567a0845e3aSMatthew G. Knepley     PetscScalar *x = NULL;
568a0845e3aSMatthew G. Knepley     PetscInt     i;
569a0845e3aSMatthew G. Knepley 
570a0845e3aSMatthew G. Knepley     ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr);
571a0845e3aSMatthew G. Knepley     if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c);
572a0845e3aSMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr);
573a0845e3aSMatthew G. Knepley     for (i = 0; i < cellDof; ++i) u[c*cellDof+i] = x[i];
574a0845e3aSMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr);
5759a559087SMatthew G. Knepley     if (dmAux) {
5769a559087SMatthew G. Knepley       ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
5779a559087SMatthew G. Knepley       for (i = 0; i < cellDofAux; ++i) a[c*cellDofAux+i] = x[i];
5789a559087SMatthew G. Knepley       ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
579a0845e3aSMatthew G. Knepley     }
5809a559087SMatthew G. Knepley   }
5819a559087SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
582c012ea0aSMatthew G. Knepley     void   (*f0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->f0Funcs[f];
583c012ea0aSMatthew G. Knepley     void   (*f1)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->f1Funcs[f];
58421454ff5SMatthew G. Knepley     PetscInt numQuadPoints, Nb;
585a0845e3aSMatthew G. Knepley     /* Conforming batches */
586f30c5766SMatthew G. Knepley     PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
587a0845e3aSMatthew G. Knepley     /* Remainder */
588a0845e3aSMatthew G. Knepley     PetscInt Nr, offset;
589a0845e3aSMatthew G. Knepley 
590a0845e3aSMatthew G. Knepley     ierr = PetscFEGetQuadrature(fe[f], &q);CHKERRQ(ierr);
591a0845e3aSMatthew G. Knepley     ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr);
592f30c5766SMatthew G. Knepley     ierr = PetscFEGetTileSizes(fe[f], NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
59321454ff5SMatthew G. Knepley     ierr = PetscQuadratureGetData(q, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr);
59421454ff5SMatthew G. Knepley     blockSize = Nb*numQuadPoints;
595a0845e3aSMatthew G. Knepley     batchSize = numBlocks * blockSize;
596f30c5766SMatthew G. Knepley     ierr =  PetscFESetTileSizes(fe[f], blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
597a0845e3aSMatthew G. Knepley     numChunks = numCells / (numBatches*batchSize);
598a0845e3aSMatthew G. Knepley     Ne        = numChunks*numBatches*batchSize;
599a0845e3aSMatthew G. Knepley     Nr        = numCells % (numBatches*batchSize);
600a0845e3aSMatthew G. Knepley     offset    = numCells - Nr;
601a0845e3aSMatthew G. Knepley     geom.v0   = v0;
602a0845e3aSMatthew G. Knepley     geom.J    = J;
603a0845e3aSMatthew G. Knepley     geom.invJ = invJ;
604a0845e3aSMatthew G. Knepley     geom.detJ = detJ;
6059a559087SMatthew G. Knepley     ierr = PetscFEIntegrateResidual(fe[f], Ne, Nf, fe, f, geom, u, NfAux, feAux, a, f0, f1, elemVec);CHKERRQ(ierr);
606a0845e3aSMatthew G. Knepley     geom.v0   = &v0[offset*dim];
607a0845e3aSMatthew G. Knepley     geom.J    = &J[offset*dim*dim];
608a0845e3aSMatthew G. Knepley     geom.invJ = &invJ[offset*dim*dim];
609a0845e3aSMatthew G. Knepley     geom.detJ = &detJ[offset];
6109a559087SMatthew 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);
611a0845e3aSMatthew G. Knepley   }
612a0845e3aSMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
613a0845e3aSMatthew G. Knepley     if (mesh->printFEM > 1) {ierr = DMPrintCellVector(c, name, cellDof, &elemVec[c*cellDof]);CHKERRQ(ierr);}
614a0845e3aSMatthew G. Knepley     ierr = DMPlexVecSetClosure(dm, section, F, c, &elemVec[c*cellDof], ADD_VALUES);CHKERRQ(ierr);
615a0845e3aSMatthew G. Knepley   }
616a0845e3aSMatthew G. Knepley   ierr = PetscFree6(u,v0,J,invJ,detJ,elemVec);CHKERRQ(ierr);
6179a559087SMatthew G. Knepley   if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);}
618f1ea0e2fSMatthew G. Knepley   if (feBd) {
619075da914SMatthew G. Knepley     DMLabel         label, depth;
620f1ea0e2fSMatthew G. Knepley     IS              pointIS;
621f1ea0e2fSMatthew G. Knepley     const PetscInt *points;
622075da914SMatthew G. Knepley     PetscInt        dep, numPoints, p, numFaces;
623f1ea0e2fSMatthew G. Knepley     PetscReal      *n;
624f1ea0e2fSMatthew G. Knepley 
625f1ea0e2fSMatthew G. Knepley     ierr = DMPlexGetLabel(dm, "boundary", &label);CHKERRQ(ierr);
626075da914SMatthew G. Knepley     ierr = DMPlexGetDepthLabel(dm, &depth);CHKERRQ(ierr);
627f1ea0e2fSMatthew G. Knepley     ierr = DMLabelGetStratumSize(label, 1, &numPoints);CHKERRQ(ierr);
628f1ea0e2fSMatthew G. Knepley     ierr = DMLabelGetStratumIS(label, 1, &pointIS);CHKERRQ(ierr);
629f1ea0e2fSMatthew G. Knepley     ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
630f1ea0e2fSMatthew G. Knepley     for (f = 0, cellDof = 0, numComponents = 0; f < Nf; ++f) {
631f1ea0e2fSMatthew G. Knepley       PetscInt Nb, Nc;
632f1ea0e2fSMatthew G. Knepley 
633f1ea0e2fSMatthew G. Knepley       ierr = PetscFEGetDimension(feBd[f], &Nb);CHKERRQ(ierr);
634f1ea0e2fSMatthew G. Knepley       ierr = PetscFEGetNumComponents(feBd[f], &Nc);CHKERRQ(ierr);
635f1ea0e2fSMatthew G. Knepley       cellDof       += Nb*Nc;
636f1ea0e2fSMatthew G. Knepley       numComponents += Nc;
637f1ea0e2fSMatthew G. Knepley     }
638075da914SMatthew G. Knepley     for (p = 0, numFaces = 0; p < numPoints; ++p) {
639075da914SMatthew G. Knepley       ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr);
640075da914SMatthew G. Knepley       if (dep == dim-1) ++numFaces;
641075da914SMatthew G. Knepley     }
642dcca6d9dSJed 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);
643075da914SMatthew G. Knepley     for (p = 0, f = 0; p < numPoints; ++p) {
644f1ea0e2fSMatthew G. Knepley       const PetscInt point = points[p];
645f1ea0e2fSMatthew G. Knepley       PetscScalar   *x     = NULL;
646f1ea0e2fSMatthew G. Knepley       PetscInt       i;
647f1ea0e2fSMatthew G. Knepley 
648075da914SMatthew G. Knepley       ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr);
649075da914SMatthew G. Knepley       if (dep != dim-1) continue;
650075da914SMatthew G. Knepley       ierr = DMPlexComputeCellGeometry(dm, point, &v0[f*dim], &J[f*dim*dim], &invJ[f*dim*dim], &detJ[f]);CHKERRQ(ierr);
651a8007bbfSMatthew G. Knepley       ierr = DMPlexComputeCellGeometryFVM(dm, point, NULL, NULL, &n[f*dim]);
652075da914SMatthew G. Knepley       if (detJ[f] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for face %d", detJ[f], point);
653f1ea0e2fSMatthew G. Knepley       ierr = DMPlexVecGetClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr);
654075da914SMatthew G. Knepley       for (i = 0; i < cellDof; ++i) u[f*cellDof+i] = x[i];
655f1ea0e2fSMatthew G. Knepley       ierr = DMPlexVecRestoreClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr);
656075da914SMatthew G. Knepley       ++f;
657f1ea0e2fSMatthew G. Knepley     }
658f1ea0e2fSMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
659f1ea0e2fSMatthew G. Knepley       void   (*f0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->f0BdFuncs[f];
660f1ea0e2fSMatthew G. Knepley       void   (*f1)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->f1BdFuncs[f];
66121454ff5SMatthew G. Knepley       PetscInt numQuadPoints, Nb;
662f1ea0e2fSMatthew G. Knepley       /* Conforming batches */
663f1ea0e2fSMatthew G. Knepley       PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
664f1ea0e2fSMatthew G. Knepley       /* Remainder */
665f1ea0e2fSMatthew G. Knepley       PetscInt Nr, offset;
666f1ea0e2fSMatthew G. Knepley 
667f1ea0e2fSMatthew G. Knepley       ierr = PetscFEGetQuadrature(feBd[f], &q);CHKERRQ(ierr);
668f1ea0e2fSMatthew G. Knepley       ierr = PetscFEGetDimension(feBd[f], &Nb);CHKERRQ(ierr);
669f1ea0e2fSMatthew G. Knepley       ierr = PetscFEGetTileSizes(feBd[f], NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
67021454ff5SMatthew G. Knepley       ierr = PetscQuadratureGetData(q, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr);
67121454ff5SMatthew G. Knepley       blockSize = Nb*numQuadPoints;
672f1ea0e2fSMatthew G. Knepley       batchSize = numBlocks * blockSize;
673f1ea0e2fSMatthew G. Knepley       ierr =  PetscFESetTileSizes(feBd[f], blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
674075da914SMatthew G. Knepley       numChunks = numFaces / (numBatches*batchSize);
675f1ea0e2fSMatthew G. Knepley       Ne        = numChunks*numBatches*batchSize;
676075da914SMatthew G. Knepley       Nr        = numFaces % (numBatches*batchSize);
677075da914SMatthew G. Knepley       offset    = numFaces - Nr;
678f1ea0e2fSMatthew G. Knepley       geom.v0   = v0;
679f1ea0e2fSMatthew G. Knepley       geom.n    = n;
680f1ea0e2fSMatthew G. Knepley       geom.J    = J;
681f1ea0e2fSMatthew G. Knepley       geom.invJ = invJ;
682f1ea0e2fSMatthew G. Knepley       geom.detJ = detJ;
683f1ea0e2fSMatthew G. Knepley       ierr = PetscFEIntegrateBdResidual(feBd[f], Ne, Nf, feBd, f, geom, u, 0, NULL, NULL, f0, f1, elemVec);CHKERRQ(ierr);
684f1ea0e2fSMatthew G. Knepley       geom.v0   = &v0[offset*dim];
685f1ea0e2fSMatthew G. Knepley       geom.n    = &n[offset*dim];
686f1ea0e2fSMatthew G. Knepley       geom.J    = &J[offset*dim*dim];
687f1ea0e2fSMatthew G. Knepley       geom.invJ = &invJ[offset*dim*dim];
688f1ea0e2fSMatthew G. Knepley       geom.detJ = &detJ[offset];
689f1ea0e2fSMatthew G. Knepley       ierr = PetscFEIntegrateBdResidual(feBd[f], Nr, Nf, feBd, f, geom, &u[offset*cellDof], 0, NULL, NULL, f0, f1, &elemVec[offset*cellDof]);CHKERRQ(ierr);
690f1ea0e2fSMatthew G. Knepley     }
691075da914SMatthew G. Knepley     for (p = 0, f = 0; p < numPoints; ++p) {
692f1ea0e2fSMatthew G. Knepley       const PetscInt point = points[p];
693f1ea0e2fSMatthew G. Knepley 
694075da914SMatthew G. Knepley       ierr = DMLabelGetValue(depth, point, &dep);CHKERRQ(ierr);
695075da914SMatthew G. Knepley       if (dep != dim-1) continue;
696075da914SMatthew G. Knepley       if (mesh->printFEM > 1) {ierr = DMPrintCellVector(point, "BdResidual", cellDof, &elemVec[f*cellDof]);CHKERRQ(ierr);}
697075da914SMatthew G. Knepley       ierr = DMPlexVecSetClosure(dm, NULL, F, point, &elemVec[f*cellDof], ADD_VALUES);CHKERRQ(ierr);
698075da914SMatthew G. Knepley       ++f;
699f1ea0e2fSMatthew G. Knepley     }
700f1ea0e2fSMatthew G. Knepley     ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
701f1ea0e2fSMatthew G. Knepley     ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
702f1ea0e2fSMatthew G. Knepley     ierr = PetscFree7(u,v0,n,J,invJ,detJ,elemVec);CHKERRQ(ierr);
703f1ea0e2fSMatthew G. Knepley   }
7046113b454SMatthew G. Knepley   if (mesh->printFEM) {ierr = DMPrintLocalVec(dm, name, mesh->printTol, F);CHKERRQ(ierr);}
705a0845e3aSMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_ResidualFEM,dm,0,0,0);CHKERRQ(ierr);
706a0845e3aSMatthew G. Knepley   PetscFunctionReturn(0);
707a0845e3aSMatthew G. Knepley }
708a0845e3aSMatthew G. Knepley 
709cb1e1211SMatthew G Knepley #undef __FUNCT__
710*af1eca97SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeIFunctionFEM"
711*af1eca97SMatthew G. Knepley /*@
712*af1eca97SMatthew G. Knepley   DMPlexComputeIFunctionFEM - Form the local implicit function F from the local input X, X_t using pointwise functions specified by the user
713*af1eca97SMatthew G. Knepley 
714*af1eca97SMatthew G. Knepley   Input Parameters:
715*af1eca97SMatthew G. Knepley + dm - The mesh
716*af1eca97SMatthew G. Knepley . X  - Local input vector
717*af1eca97SMatthew G. Knepley . X_t  - Time derivative of the local input vector
718*af1eca97SMatthew G. Knepley - user - The user context
719*af1eca97SMatthew G. Knepley 
720*af1eca97SMatthew G. Knepley   Output Parameter:
721*af1eca97SMatthew G. Knepley . F  - Local output vector
722*af1eca97SMatthew G. Knepley 
723*af1eca97SMatthew G. Knepley   Note:
724*af1eca97SMatthew G. Knepley   The first member of the user context must be an FEMContext.
725*af1eca97SMatthew G. Knepley 
726*af1eca97SMatthew G. Knepley   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
727*af1eca97SMatthew G. Knepley   like a GPU, or vectorize on a multicore machine.
728*af1eca97SMatthew G. Knepley 
729*af1eca97SMatthew G. Knepley   Level: developer
730*af1eca97SMatthew G. Knepley 
731*af1eca97SMatthew G. Knepley .seealso: DMPlexComputeResidualFEM()
732*af1eca97SMatthew G. Knepley @*/
733*af1eca97SMatthew G. Knepley PetscErrorCode DMPlexComputeIFunctionFEM(DM dm, Vec X, Vec X_t, Vec F, void *user)
734*af1eca97SMatthew G. Knepley {
735*af1eca97SMatthew G. Knepley   DM_Plex          *mesh  = (DM_Plex *) dm->data;
736*af1eca97SMatthew G. Knepley   PetscFEM         *fem   = (PetscFEM *) user;
737*af1eca97SMatthew G. Knepley   PetscFE          *fe    = fem->fe;
738*af1eca97SMatthew G. Knepley   PetscFE          *feAux = fem->feAux;
739*af1eca97SMatthew G. Knepley   PetscFE          *feBd  = fem->feBd;
740*af1eca97SMatthew G. Knepley   const char       *name  = "Residual";
741*af1eca97SMatthew G. Knepley   DM                dmAux;
742*af1eca97SMatthew G. Knepley   Vec               A;
743*af1eca97SMatthew G. Knepley   PetscQuadrature   q;
744*af1eca97SMatthew G. Knepley   PetscCellGeometry geom;
745*af1eca97SMatthew G. Knepley   PetscSection      section, sectionAux;
746*af1eca97SMatthew G. Knepley   PetscReal        *v0, *J, *invJ, *detJ;
747*af1eca97SMatthew G. Knepley   PetscScalar      *elemVec, *u, *u_t, *a = NULL;
748*af1eca97SMatthew G. Knepley   PetscInt          dim, Nf, NfAux = 0, f, numCells, cStart, cEnd, c;
749*af1eca97SMatthew G. Knepley   PetscInt          cellDof = 0, numComponents = 0;
750*af1eca97SMatthew G. Knepley   PetscInt          cellDofAux = 0, numComponentsAux = 0;
751*af1eca97SMatthew G. Knepley   PetscErrorCode    ierr;
752*af1eca97SMatthew G. Knepley 
753*af1eca97SMatthew G. Knepley   PetscFunctionBegin;
754*af1eca97SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_ResidualFEM,dm,0,0,0);CHKERRQ(ierr);
755*af1eca97SMatthew G. Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
756*af1eca97SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
757*af1eca97SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
758*af1eca97SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
759*af1eca97SMatthew G. Knepley   numCells = cEnd - cStart;
760*af1eca97SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
761*af1eca97SMatthew G. Knepley     PetscInt Nb, Nc;
762*af1eca97SMatthew G. Knepley 
763*af1eca97SMatthew G. Knepley     ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr);
764*af1eca97SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe[f], &Nc);CHKERRQ(ierr);
765*af1eca97SMatthew G. Knepley     cellDof       += Nb*Nc;
766*af1eca97SMatthew G. Knepley     numComponents += Nc;
767*af1eca97SMatthew G. Knepley   }
768*af1eca97SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
769*af1eca97SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr);
770*af1eca97SMatthew G. Knepley   if (dmAux) {
771*af1eca97SMatthew G. Knepley     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
772*af1eca97SMatthew G. Knepley     ierr = PetscSectionGetNumFields(sectionAux, &NfAux);CHKERRQ(ierr);
773*af1eca97SMatthew G. Knepley   }
774*af1eca97SMatthew G. Knepley   for (f = 0; f < NfAux; ++f) {
775*af1eca97SMatthew G. Knepley     PetscInt Nb, Nc;
776*af1eca97SMatthew G. Knepley 
777*af1eca97SMatthew G. Knepley     ierr = PetscFEGetDimension(feAux[f], &Nb);CHKERRQ(ierr);
778*af1eca97SMatthew G. Knepley     ierr = PetscFEGetNumComponents(feAux[f], &Nc);CHKERRQ(ierr);
779*af1eca97SMatthew G. Knepley     cellDofAux       += Nb*Nc;
780*af1eca97SMatthew G. Knepley     numComponentsAux += Nc;
781*af1eca97SMatthew G. Knepley   }
782*af1eca97SMatthew G. Knepley   ierr = DMPlexProjectFunctionLocal(dm, fe, fem->bcFuncs, fem->bcCtxs, INSERT_BC_VALUES, X);CHKERRQ(ierr);
783*af1eca97SMatthew G. Knepley   ierr = VecSet(F, 0.0);CHKERRQ(ierr);
784*af1eca97SMatthew 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);
785*af1eca97SMatthew G. Knepley   if (dmAux) {ierr = PetscMalloc1(numCells*cellDofAux, &a);CHKERRQ(ierr);}
786*af1eca97SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
787*af1eca97SMatthew G. Knepley     PetscScalar *x = NULL, *x_t = NULL;
788*af1eca97SMatthew G. Knepley     PetscInt     i;
789*af1eca97SMatthew G. Knepley 
790*af1eca97SMatthew G. Knepley     ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr);
791*af1eca97SMatthew G. Knepley     if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c);
792*af1eca97SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr);
793*af1eca97SMatthew G. Knepley     for (i = 0; i < cellDof; ++i) u[c*cellDof+i] = x[i];
794*af1eca97SMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr);
795*af1eca97SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, section, X_t, c, NULL, &x_t);CHKERRQ(ierr);
796*af1eca97SMatthew G. Knepley     for (i = 0; i < cellDof; ++i) u_t[c*cellDof+i] = x_t[i];
797*af1eca97SMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, section, X_t, c, NULL, &x_t);CHKERRQ(ierr);
798*af1eca97SMatthew G. Knepley     if (dmAux) {
799*af1eca97SMatthew G. Knepley       PetscScalar *x_a = NULL;
800*af1eca97SMatthew G. Knepley       ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x_a);CHKERRQ(ierr);
801*af1eca97SMatthew G. Knepley       for (i = 0; i < cellDofAux; ++i) a[c*cellDofAux+i] = x_a[i];
802*af1eca97SMatthew G. Knepley       ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x_a);CHKERRQ(ierr);
803*af1eca97SMatthew G. Knepley     }
804*af1eca97SMatthew G. Knepley   }
805*af1eca97SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
806*af1eca97SMatthew G. Knepley     void   (*f0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->f0IFuncs[f];
807*af1eca97SMatthew G. Knepley     void   (*f1)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->f1IFuncs[f];
808*af1eca97SMatthew G. Knepley     PetscInt numQuadPoints, Nb;
809*af1eca97SMatthew G. Knepley     /* Conforming batches */
810*af1eca97SMatthew G. Knepley     PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
811*af1eca97SMatthew G. Knepley     /* Remainder */
812*af1eca97SMatthew G. Knepley     PetscInt Nr, offset;
813*af1eca97SMatthew G. Knepley 
814*af1eca97SMatthew G. Knepley     ierr = PetscFEGetQuadrature(fe[f], &q);CHKERRQ(ierr);
815*af1eca97SMatthew G. Knepley     ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr);
816*af1eca97SMatthew G. Knepley     ierr = PetscFEGetTileSizes(fe[f], NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
817*af1eca97SMatthew G. Knepley     ierr = PetscQuadratureGetData(q, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr);
818*af1eca97SMatthew G. Knepley     blockSize = Nb*numQuadPoints;
819*af1eca97SMatthew G. Knepley     batchSize = numBlocks * blockSize;
820*af1eca97SMatthew G. Knepley     ierr =  PetscFESetTileSizes(fe[f], blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
821*af1eca97SMatthew G. Knepley     numChunks = numCells / (numBatches*batchSize);
822*af1eca97SMatthew G. Knepley     Ne        = numChunks*numBatches*batchSize;
823*af1eca97SMatthew G. Knepley     Nr        = numCells % (numBatches*batchSize);
824*af1eca97SMatthew G. Knepley     offset    = numCells - Nr;
825*af1eca97SMatthew G. Knepley     geom.v0   = v0;
826*af1eca97SMatthew G. Knepley     geom.J    = J;
827*af1eca97SMatthew G. Knepley     geom.invJ = invJ;
828*af1eca97SMatthew G. Knepley     geom.detJ = detJ;
829*af1eca97SMatthew G. Knepley     ierr = PetscFEIntegrateIFunction(fe[f], Ne, Nf, fe, f, geom, u, u_t, NfAux, feAux, a, f0, f1, elemVec);CHKERRQ(ierr);
830*af1eca97SMatthew G. Knepley     geom.v0   = &v0[offset*dim];
831*af1eca97SMatthew G. Knepley     geom.J    = &J[offset*dim*dim];
832*af1eca97SMatthew G. Knepley     geom.invJ = &invJ[offset*dim*dim];
833*af1eca97SMatthew G. Knepley     geom.detJ = &detJ[offset];
834*af1eca97SMatthew 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);
835*af1eca97SMatthew G. Knepley   }
836*af1eca97SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
837*af1eca97SMatthew G. Knepley     if (mesh->printFEM > 1) {ierr = DMPrintCellVector(c, name, cellDof, &elemVec[c*cellDof]);CHKERRQ(ierr);}
838*af1eca97SMatthew G. Knepley     ierr = DMPlexVecSetClosure(dm, section, F, c, &elemVec[c*cellDof], ADD_VALUES);CHKERRQ(ierr);
839*af1eca97SMatthew G. Knepley   }
840*af1eca97SMatthew G. Knepley   ierr = PetscFree7(u,u_t,v0,J,invJ,detJ,elemVec);CHKERRQ(ierr);
841*af1eca97SMatthew G. Knepley   if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);}
842*af1eca97SMatthew G. Knepley   if (feBd) {
843*af1eca97SMatthew G. Knepley     DMLabel         label, depth;
844*af1eca97SMatthew G. Knepley     IS              pointIS;
845*af1eca97SMatthew G. Knepley     const PetscInt *points;
846*af1eca97SMatthew G. Knepley     PetscInt        dep, numPoints, p, numFaces;
847*af1eca97SMatthew G. Knepley     PetscReal      *n;
848*af1eca97SMatthew G. Knepley 
849*af1eca97SMatthew G. Knepley     ierr = DMPlexGetLabel(dm, "boundary", &label);CHKERRQ(ierr);
850*af1eca97SMatthew G. Knepley     ierr = DMPlexGetDepthLabel(dm, &depth);CHKERRQ(ierr);
851*af1eca97SMatthew G. Knepley     ierr = DMLabelGetStratumSize(label, 1, &numPoints);CHKERRQ(ierr);
852*af1eca97SMatthew G. Knepley     ierr = DMLabelGetStratumIS(label, 1, &pointIS);CHKERRQ(ierr);
853*af1eca97SMatthew G. Knepley     ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
854*af1eca97SMatthew G. Knepley     for (f = 0, cellDof = 0, numComponents = 0; f < Nf; ++f) {
855*af1eca97SMatthew G. Knepley       PetscInt Nb, Nc;
856*af1eca97SMatthew G. Knepley 
857*af1eca97SMatthew G. Knepley       ierr = PetscFEGetDimension(feBd[f], &Nb);CHKERRQ(ierr);
858*af1eca97SMatthew G. Knepley       ierr = PetscFEGetNumComponents(feBd[f], &Nc);CHKERRQ(ierr);
859*af1eca97SMatthew G. Knepley       cellDof       += Nb*Nc;
860*af1eca97SMatthew G. Knepley       numComponents += Nc;
861*af1eca97SMatthew G. Knepley     }
862*af1eca97SMatthew G. Knepley     for (p = 0, numFaces = 0; p < numPoints; ++p) {
863*af1eca97SMatthew G. Knepley       ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr);
864*af1eca97SMatthew G. Knepley       if (dep == dim-1) ++numFaces;
865*af1eca97SMatthew G. Knepley     }
866*af1eca97SMatthew 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);
867*af1eca97SMatthew G. Knepley     ierr = PetscMalloc1(numFaces*cellDof,&u_t);CHKERRQ(ierr);
868*af1eca97SMatthew G. Knepley     for (p = 0, f = 0; p < numPoints; ++p) {
869*af1eca97SMatthew G. Knepley       const PetscInt point = points[p];
870*af1eca97SMatthew G. Knepley       PetscScalar   *x     = NULL;
871*af1eca97SMatthew G. Knepley       PetscInt       i;
872*af1eca97SMatthew G. Knepley 
873*af1eca97SMatthew G. Knepley       ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr);
874*af1eca97SMatthew G. Knepley       if (dep != dim-1) continue;
875*af1eca97SMatthew G. Knepley       ierr = DMPlexComputeCellGeometry(dm, point, &v0[f*dim], &J[f*dim*dim], &invJ[f*dim*dim], &detJ[f]);CHKERRQ(ierr);
876*af1eca97SMatthew G. Knepley       ierr = DMPlexComputeCellGeometryFVM(dm, point, NULL, NULL, &n[f*dim]);
877*af1eca97SMatthew G. Knepley       if (detJ[f] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for face %d", detJ[f], point);
878*af1eca97SMatthew G. Knepley       ierr = DMPlexVecGetClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr);
879*af1eca97SMatthew G. Knepley       for (i = 0; i < cellDof; ++i) u[f*cellDof+i] = x[i];
880*af1eca97SMatthew G. Knepley       ierr = DMPlexVecRestoreClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr);
881*af1eca97SMatthew G. Knepley       ierr = DMPlexVecGetClosure(dm, section, X_t, point, NULL, &x);CHKERRQ(ierr);
882*af1eca97SMatthew G. Knepley       for (i = 0; i < cellDof; ++i) u_t[f*cellDof+i] = x[i];
883*af1eca97SMatthew G. Knepley       ierr = DMPlexVecRestoreClosure(dm, section, X_t, point, NULL, &x);CHKERRQ(ierr);
884*af1eca97SMatthew G. Knepley       ++f;
885*af1eca97SMatthew G. Knepley     }
886*af1eca97SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
887*af1eca97SMatthew G. Knepley       void   (*f0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->f0BdIFuncs[f];
888*af1eca97SMatthew G. Knepley       void   (*f1)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->f1BdIFuncs[f];
889*af1eca97SMatthew G. Knepley       PetscInt numQuadPoints, Nb;
890*af1eca97SMatthew G. Knepley       /* Conforming batches */
891*af1eca97SMatthew G. Knepley       PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
892*af1eca97SMatthew G. Knepley       /* Remainder */
893*af1eca97SMatthew G. Knepley       PetscInt Nr, offset;
894*af1eca97SMatthew G. Knepley 
895*af1eca97SMatthew G. Knepley       ierr = PetscFEGetQuadrature(feBd[f], &q);CHKERRQ(ierr);
896*af1eca97SMatthew G. Knepley       ierr = PetscFEGetDimension(feBd[f], &Nb);CHKERRQ(ierr);
897*af1eca97SMatthew G. Knepley       ierr = PetscFEGetTileSizes(feBd[f], NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
898*af1eca97SMatthew G. Knepley       ierr = PetscQuadratureGetData(q, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr);
899*af1eca97SMatthew G. Knepley       blockSize = Nb*numQuadPoints;
900*af1eca97SMatthew G. Knepley       batchSize = numBlocks * blockSize;
901*af1eca97SMatthew G. Knepley       ierr =  PetscFESetTileSizes(feBd[f], blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
902*af1eca97SMatthew G. Knepley       numChunks = numFaces / (numBatches*batchSize);
903*af1eca97SMatthew G. Knepley       Ne        = numChunks*numBatches*batchSize;
904*af1eca97SMatthew G. Knepley       Nr        = numFaces % (numBatches*batchSize);
905*af1eca97SMatthew G. Knepley       offset    = numFaces - Nr;
906*af1eca97SMatthew G. Knepley       geom.v0   = v0;
907*af1eca97SMatthew G. Knepley       geom.n    = n;
908*af1eca97SMatthew G. Knepley       geom.J    = J;
909*af1eca97SMatthew G. Knepley       geom.invJ = invJ;
910*af1eca97SMatthew G. Knepley       geom.detJ = detJ;
911*af1eca97SMatthew G. Knepley       ierr = PetscFEIntegrateBdIFunction(feBd[f], Ne, Nf, feBd, f, geom, u, u_t, 0, NULL, NULL, f0, f1, elemVec);CHKERRQ(ierr);
912*af1eca97SMatthew G. Knepley       geom.v0   = &v0[offset*dim];
913*af1eca97SMatthew G. Knepley       geom.n    = &n[offset*dim];
914*af1eca97SMatthew G. Knepley       geom.J    = &J[offset*dim*dim];
915*af1eca97SMatthew G. Knepley       geom.invJ = &invJ[offset*dim*dim];
916*af1eca97SMatthew G. Knepley       geom.detJ = &detJ[offset];
917*af1eca97SMatthew 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);
918*af1eca97SMatthew G. Knepley     }
919*af1eca97SMatthew G. Knepley     for (p = 0, f = 0; p < numPoints; ++p) {
920*af1eca97SMatthew G. Knepley       const PetscInt point = points[p];
921*af1eca97SMatthew G. Knepley 
922*af1eca97SMatthew G. Knepley       ierr = DMLabelGetValue(depth, point, &dep);CHKERRQ(ierr);
923*af1eca97SMatthew G. Knepley       if (dep != dim-1) continue;
924*af1eca97SMatthew G. Knepley       if (mesh->printFEM > 1) {ierr = DMPrintCellVector(point, "BdResidual", cellDof, &elemVec[f*cellDof]);CHKERRQ(ierr);}
925*af1eca97SMatthew G. Knepley       ierr = DMPlexVecSetClosure(dm, NULL, F, point, &elemVec[f*cellDof], ADD_VALUES);CHKERRQ(ierr);
926*af1eca97SMatthew G. Knepley       ++f;
927*af1eca97SMatthew G. Knepley     }
928*af1eca97SMatthew G. Knepley     ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
929*af1eca97SMatthew G. Knepley     ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
930*af1eca97SMatthew G. Knepley     ierr = PetscFree7(u,v0,n,J,invJ,detJ,elemVec);CHKERRQ(ierr);
931*af1eca97SMatthew G. Knepley     ierr = PetscFree(u_t);CHKERRQ(ierr);
932*af1eca97SMatthew G. Knepley   }
933*af1eca97SMatthew G. Knepley   if (mesh->printFEM) {ierr = DMPrintLocalVec(dm, name, mesh->printTol, F);CHKERRQ(ierr);}
934*af1eca97SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_ResidualFEM,dm,0,0,0);CHKERRQ(ierr);
935*af1eca97SMatthew G. Knepley   PetscFunctionReturn(0);
936*af1eca97SMatthew G. Knepley }
937*af1eca97SMatthew G. Knepley 
938*af1eca97SMatthew G. Knepley #undef __FUNCT__
939cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexComputeJacobianActionFEM"
940cb1e1211SMatthew G Knepley /*@C
941cb1e1211SMatthew G Knepley   DMPlexComputeJacobianActionFEM - Form the local action of Jacobian J(u) on the local input X using pointwise functions specified by the user
942cb1e1211SMatthew G Knepley 
943cb1e1211SMatthew G Knepley   Input Parameters:
944cb1e1211SMatthew G Knepley + dm - The mesh
945cb1e1211SMatthew G Knepley . J  - The Jacobian shell matrix
946cb1e1211SMatthew G Knepley . X  - Local input vector
947cb1e1211SMatthew G Knepley - user - The user context
948cb1e1211SMatthew G Knepley 
949cb1e1211SMatthew G Knepley   Output Parameter:
950cb1e1211SMatthew G Knepley . F  - Local output vector
951cb1e1211SMatthew G Knepley 
952cb1e1211SMatthew G Knepley   Note:
9538026896cSMatthew G. Knepley   The first member of the user context must be an FEMContext.
954cb1e1211SMatthew G Knepley 
955cb1e1211SMatthew G Knepley   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
956cb1e1211SMatthew G Knepley   like a GPU, or vectorize on a multicore machine.
957cb1e1211SMatthew G Knepley 
9580059ad2aSSatish Balay   Level: developer
9590059ad2aSSatish Balay 
960cb1e1211SMatthew G Knepley .seealso: DMPlexComputeResidualFEM()
961878cb397SSatish Balay @*/
962cb1e1211SMatthew G Knepley PetscErrorCode DMPlexComputeJacobianActionFEM(DM dm, Mat Jac, Vec X, Vec F, void *user)
963cb1e1211SMatthew G Knepley {
964cb1e1211SMatthew G Knepley   DM_Plex          *mesh = (DM_Plex *) dm->data;
9659a559087SMatthew G. Knepley   PetscFEM         *fem  = (PetscFEM *) user;
9660483ade4SMatthew G. Knepley   PetscFE          *fe   = fem->fe;
9670483ade4SMatthew G. Knepley   PetscQuadrature   quad;
9680483ade4SMatthew G. Knepley   PetscCellGeometry geom;
969cb1e1211SMatthew G Knepley   PetscSection      section;
970cb1e1211SMatthew G Knepley   JacActionCtx     *jctx;
971cb1e1211SMatthew G Knepley   PetscReal        *v0, *J, *invJ, *detJ;
972cb1e1211SMatthew G Knepley   PetscScalar      *elemVec, *u, *a;
9730483ade4SMatthew G. Knepley   PetscInt          dim, numFields, field, numCells, cStart, cEnd, c;
974cb1e1211SMatthew G Knepley   PetscInt          cellDof = 0;
975cb1e1211SMatthew G Knepley   PetscErrorCode    ierr;
976cb1e1211SMatthew G Knepley 
977cb1e1211SMatthew G Knepley   PetscFunctionBegin;
9780483ade4SMatthew G. Knepley   /* ierr = PetscLogEventBegin(DMPLEX_JacobianActionFEM,dm,0,0,0);CHKERRQ(ierr); */
979cb1e1211SMatthew G Knepley   ierr = MatShellGetContext(Jac, &jctx);CHKERRQ(ierr);
980cb1e1211SMatthew G Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
981cb1e1211SMatthew G Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
982cb1e1211SMatthew G Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
983cb1e1211SMatthew G Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
984cb1e1211SMatthew G Knepley   numCells = cEnd - cStart;
985cb1e1211SMatthew G Knepley   for (field = 0; field < numFields; ++field) {
9860483ade4SMatthew G. Knepley     PetscInt Nb, Nc;
9870483ade4SMatthew G. Knepley 
9880483ade4SMatthew G. Knepley     ierr = PetscFEGetDimension(fe[field], &Nb);CHKERRQ(ierr);
9890483ade4SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe[field], &Nc);CHKERRQ(ierr);
9900483ade4SMatthew G. Knepley     cellDof += Nb*Nc;
991cb1e1211SMatthew G Knepley   }
992cb1e1211SMatthew G Knepley   ierr = VecSet(F, 0.0);CHKERRQ(ierr);
993dcca6d9dSJed 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);
994cb1e1211SMatthew G Knepley   for (c = cStart; c < cEnd; ++c) {
995a1e44745SMatthew G. Knepley     PetscScalar *x = NULL;
996cb1e1211SMatthew G Knepley     PetscInt     i;
997cb1e1211SMatthew G Knepley 
998cb1e1211SMatthew G Knepley     ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr);
999cb1e1211SMatthew G Knepley     if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c);
1000cb1e1211SMatthew G Knepley     ierr = DMPlexVecGetClosure(dm, NULL, jctx->u, c, NULL, &x);CHKERRQ(ierr);
1001cb1e1211SMatthew G Knepley     for (i = 0; i < cellDof; ++i) u[c*cellDof+i] = x[i];
1002cb1e1211SMatthew G Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, jctx->u, c, NULL, &x);CHKERRQ(ierr);
1003cb1e1211SMatthew G Knepley     ierr = DMPlexVecGetClosure(dm, NULL, X, c, NULL, &x);CHKERRQ(ierr);
1004cb1e1211SMatthew G Knepley     for (i = 0; i < cellDof; ++i) a[c*cellDof+i] = x[i];
1005cb1e1211SMatthew G Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, X, c, NULL, &x);CHKERRQ(ierr);
1006cb1e1211SMatthew G Knepley   }
1007cb1e1211SMatthew G Knepley   for (field = 0; field < numFields; ++field) {
100821454ff5SMatthew G. Knepley     PetscInt numQuadPoints, Nb;
1009cb1e1211SMatthew G Knepley     /* Conforming batches */
1010cb1e1211SMatthew G Knepley     PetscInt numBlocks  = 1;
10110483ade4SMatthew G. Knepley     PetscInt numBatches = 1;
10120483ade4SMatthew G. Knepley     PetscInt numChunks, Ne, blockSize, batchSize;
1013cb1e1211SMatthew G Knepley     /* Remainder */
10140483ade4SMatthew G. Knepley     PetscInt Nr, offset;
1015cb1e1211SMatthew G Knepley 
10160483ade4SMatthew G. Knepley     ierr = PetscFEGetQuadrature(fe[field], &quad);CHKERRQ(ierr);
10170483ade4SMatthew G. Knepley     ierr = PetscFEGetDimension(fe[field], &Nb);CHKERRQ(ierr);
101821454ff5SMatthew G. Knepley     ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr);
101921454ff5SMatthew G. Knepley     blockSize = Nb*numQuadPoints;
10200483ade4SMatthew G. Knepley     batchSize = numBlocks * blockSize;
10210483ade4SMatthew G. Knepley     numChunks = numCells / (numBatches*batchSize);
10220483ade4SMatthew G. Knepley     Ne        = numChunks*numBatches*batchSize;
10230483ade4SMatthew G. Knepley     Nr        = numCells % (numBatches*batchSize);
10240483ade4SMatthew G. Knepley     offset    = numCells - Nr;
10250483ade4SMatthew G. Knepley     geom.v0   = v0;
10260483ade4SMatthew G. Knepley     geom.J    = J;
10270483ade4SMatthew G. Knepley     geom.invJ = invJ;
10280483ade4SMatthew G. Knepley     geom.detJ = detJ;
10290483ade4SMatthew G. Knepley     ierr = PetscFEIntegrateJacobianAction(fe[field], Ne, numFields, fe, field, geom, u, a, fem->g0Funcs, fem->g1Funcs, fem->g2Funcs, fem->g3Funcs, elemVec);CHKERRQ(ierr);
10300483ade4SMatthew G. Knepley     geom.v0   = &v0[offset*dim];
10310483ade4SMatthew G. Knepley     geom.J    = &J[offset*dim*dim];
10320483ade4SMatthew G. Knepley     geom.invJ = &invJ[offset*dim*dim];
10330483ade4SMatthew G. Knepley     geom.detJ = &detJ[offset];
10340483ade4SMatthew G. Knepley     ierr = PetscFEIntegrateJacobianAction(fe[field], Nr, numFields, fe, field, geom, &u[offset*cellDof], &a[offset*cellDof],
1035cb1e1211SMatthew G Knepley                                           fem->g0Funcs, fem->g1Funcs, fem->g2Funcs, fem->g3Funcs, &elemVec[offset*cellDof]);CHKERRQ(ierr);
1036cb1e1211SMatthew G Knepley   }
1037cb1e1211SMatthew G Knepley   for (c = cStart; c < cEnd; ++c) {
1038cb1e1211SMatthew G Knepley     if (mesh->printFEM > 1) {ierr = DMPrintCellVector(c, "Jacobian Action", cellDof, &elemVec[c*cellDof]);CHKERRQ(ierr);}
1039cb1e1211SMatthew G Knepley     ierr = DMPlexVecSetClosure(dm, NULL, F, c, &elemVec[c*cellDof], ADD_VALUES);CHKERRQ(ierr);
1040cb1e1211SMatthew G Knepley   }
1041cb1e1211SMatthew G Knepley   ierr = PetscFree7(u,a,v0,J,invJ,detJ,elemVec);CHKERRQ(ierr);
1042cb1e1211SMatthew G Knepley   if (mesh->printFEM) {
1043cb1e1211SMatthew G Knepley     PetscMPIInt rank, numProcs;
1044cb1e1211SMatthew G Knepley     PetscInt    p;
1045cb1e1211SMatthew G Knepley 
1046cb1e1211SMatthew G Knepley     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr);
1047cb1e1211SMatthew G Knepley     ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm), &numProcs);CHKERRQ(ierr);
104886a74ee0SMatthew G. Knepley     ierr = PetscPrintf(PetscObjectComm((PetscObject)dm), "Jacobian Action:\n");CHKERRQ(ierr);
1049cb1e1211SMatthew G Knepley     for (p = 0; p < numProcs; ++p) {
1050cb1e1211SMatthew G Knepley       if (p == rank) {ierr = VecView(F, PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);}
1051cb1e1211SMatthew G Knepley       ierr = PetscBarrier((PetscObject) dm);CHKERRQ(ierr);
1052cb1e1211SMatthew G Knepley     }
1053cb1e1211SMatthew G Knepley   }
10540483ade4SMatthew G. Knepley   /* ierr = PetscLogEventEnd(DMPLEX_JacobianActionFEM,dm,0,0,0);CHKERRQ(ierr); */
1055cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
1056cb1e1211SMatthew G Knepley }
1057cb1e1211SMatthew G Knepley 
1058cb1e1211SMatthew G Knepley #undef __FUNCT__
1059cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexComputeJacobianFEM"
1060cb1e1211SMatthew G Knepley /*@
1061cb1e1211SMatthew G Knepley   DMPlexComputeJacobianFEM - Form the local portion of the Jacobian matrix J at the local solution X using pointwise functions specified by the user.
1062cb1e1211SMatthew G Knepley 
1063cb1e1211SMatthew G Knepley   Input Parameters:
1064cb1e1211SMatthew G Knepley + dm - The mesh
1065cb1e1211SMatthew G Knepley . X  - Local input vector
1066cb1e1211SMatthew G Knepley - user - The user context
1067cb1e1211SMatthew G Knepley 
1068cb1e1211SMatthew G Knepley   Output Parameter:
1069cb1e1211SMatthew G Knepley . Jac  - Jacobian matrix
1070cb1e1211SMatthew G Knepley 
1071cb1e1211SMatthew G Knepley   Note:
10728026896cSMatthew G. Knepley   The first member of the user context must be an FEMContext.
1073cb1e1211SMatthew G Knepley 
1074cb1e1211SMatthew G Knepley   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
1075cb1e1211SMatthew G Knepley   like a GPU, or vectorize on a multicore machine.
1076cb1e1211SMatthew G Knepley 
10770059ad2aSSatish Balay   Level: developer
10780059ad2aSSatish Balay 
1079cb1e1211SMatthew G Knepley .seealso: FormFunctionLocal()
1080878cb397SSatish Balay @*/
10810405ed22SMatthew G. Knepley PetscErrorCode DMPlexComputeJacobianFEM(DM dm, Vec X, Mat Jac, Mat JacP,void *user)
1082cb1e1211SMatthew G Knepley {
1083cb1e1211SMatthew G Knepley   DM_Plex          *mesh  = (DM_Plex *) dm->data;
10849a559087SMatthew G. Knepley   PetscFEM         *fem   = (PetscFEM *) user;
1085a319912fSMatthew G. Knepley   PetscFE          *fe    = fem->fe;
1086754551f4SMatthew G. Knepley   PetscFE          *feAux = fem->feAux;
1087a319912fSMatthew G. Knepley   const char       *name  = "Jacobian";
1088754551f4SMatthew G. Knepley   DM                dmAux;
1089754551f4SMatthew G. Knepley   Vec               A;
1090a319912fSMatthew G. Knepley   PetscQuadrature   quad;
1091a319912fSMatthew G. Knepley   PetscCellGeometry geom;
1092754551f4SMatthew G. Knepley   PetscSection      section, globalSection, sectionAux;
1093cb1e1211SMatthew G Knepley   PetscReal        *v0, *J, *invJ, *detJ;
1094754551f4SMatthew G. Knepley   PetscScalar      *elemMat, *u, *a;
1095754551f4SMatthew G. Knepley   PetscInt          dim, Nf, NfAux = 0, f, fieldI, fieldJ, numCells, cStart, cEnd, c;
1096cb1e1211SMatthew G Knepley   PetscInt          cellDof = 0, numComponents = 0;
1097754551f4SMatthew G. Knepley   PetscInt          cellDofAux = 0, numComponentsAux = 0;
1098cb1e1211SMatthew G Knepley   PetscBool         isShell;
1099cb1e1211SMatthew G Knepley   PetscErrorCode    ierr;
1100cb1e1211SMatthew G Knepley 
1101cb1e1211SMatthew G Knepley   PetscFunctionBegin;
1102a319912fSMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_JacobianFEM,dm,0,0,0);CHKERRQ(ierr);
1103cb1e1211SMatthew G Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
1104cb1e1211SMatthew G Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
1105a319912fSMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dm, &globalSection);CHKERRQ(ierr);
1106754551f4SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
1107cb1e1211SMatthew G Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
1108cb1e1211SMatthew G Knepley   numCells = cEnd - cStart;
1109754551f4SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
1110a319912fSMatthew G. Knepley     PetscInt Nb, Nc;
1111a319912fSMatthew G. Knepley 
1112a319912fSMatthew G. Knepley     ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr);
1113a319912fSMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe[f], &Nc);CHKERRQ(ierr);
1114a319912fSMatthew G. Knepley     cellDof       += Nb*Nc;
1115a319912fSMatthew G. Knepley     numComponents += Nc;
1116cb1e1211SMatthew G Knepley   }
1117754551f4SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
1118754551f4SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr);
1119754551f4SMatthew G. Knepley   if (dmAux) {
1120754551f4SMatthew G. Knepley     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
1121754551f4SMatthew G. Knepley     ierr = PetscSectionGetNumFields(sectionAux, &NfAux);CHKERRQ(ierr);
1122754551f4SMatthew G. Knepley   }
1123754551f4SMatthew G. Knepley   for (f = 0; f < NfAux; ++f) {
1124754551f4SMatthew G. Knepley     PetscInt Nb, Nc;
1125754551f4SMatthew G. Knepley 
1126754551f4SMatthew G. Knepley     ierr = PetscFEGetDimension(feAux[f], &Nb);CHKERRQ(ierr);
1127754551f4SMatthew G. Knepley     ierr = PetscFEGetNumComponents(feAux[f], &Nc);CHKERRQ(ierr);
1128754551f4SMatthew G. Knepley     cellDofAux       += Nb*Nc;
1129754551f4SMatthew G. Knepley     numComponentsAux += Nc;
1130754551f4SMatthew G. Knepley   }
1131c110b1eeSGeoffrey Irving   ierr = DMPlexProjectFunctionLocal(dm, fe, fem->bcFuncs, fem->bcCtxs, INSERT_BC_VALUES, X);CHKERRQ(ierr);
1132cb1e1211SMatthew G Knepley   ierr = MatZeroEntries(JacP);CHKERRQ(ierr);
1133dcca6d9dSJed 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);
1134785e854fSJed Brown   if (dmAux) {ierr = PetscMalloc1(numCells*cellDofAux, &a);CHKERRQ(ierr);}
1135cb1e1211SMatthew G Knepley   for (c = cStart; c < cEnd; ++c) {
1136a1e44745SMatthew G. Knepley     PetscScalar *x = NULL;
1137cb1e1211SMatthew G Knepley     PetscInt     i;
1138cb1e1211SMatthew G Knepley 
1139cb1e1211SMatthew G Knepley     ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr);
1140cb1e1211SMatthew G Knepley     if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c);
1141a319912fSMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr);
1142cb1e1211SMatthew G Knepley     for (i = 0; i < cellDof; ++i) u[c*cellDof+i] = x[i];
1143a319912fSMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr);
1144754551f4SMatthew G. Knepley     if (dmAux) {
1145754551f4SMatthew G. Knepley       ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
1146754551f4SMatthew G. Knepley       for (i = 0; i < cellDofAux; ++i) a[c*cellDofAux+i] = x[i];
1147754551f4SMatthew G. Knepley       ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
1148754551f4SMatthew G. Knepley     }
1149cb1e1211SMatthew G Knepley   }
1150cb1e1211SMatthew G Knepley   ierr = PetscMemzero(elemMat, numCells*cellDof*cellDof * sizeof(PetscScalar));CHKERRQ(ierr);
1151754551f4SMatthew G. Knepley   for (fieldI = 0; fieldI < Nf; ++fieldI) {
115221454ff5SMatthew G. Knepley     PetscInt numQuadPoints, Nb;
1153cb1e1211SMatthew G Knepley     /* Conforming batches */
1154754551f4SMatthew G. Knepley     PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
1155cb1e1211SMatthew G Knepley     /* Remainder */
1156a319912fSMatthew G. Knepley     PetscInt Nr, offset;
1157cb1e1211SMatthew G Knepley 
1158754551f4SMatthew G. Knepley     ierr = PetscFEGetQuadrature(fe[fieldI], &quad);CHKERRQ(ierr);
1159754551f4SMatthew G. Knepley     ierr = PetscFEGetDimension(fe[fieldI], &Nb);CHKERRQ(ierr);
1160754551f4SMatthew G. Knepley     ierr = PetscFEGetTileSizes(fe[fieldI], NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
116121454ff5SMatthew G. Knepley     ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr);
116221454ff5SMatthew G. Knepley     blockSize = Nb*numQuadPoints;
1163a319912fSMatthew G. Knepley     batchSize = numBlocks * blockSize;
1164754551f4SMatthew G. Knepley     ierr = PetscFESetTileSizes(fe[fieldI], blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
1165a319912fSMatthew G. Knepley     numChunks = numCells / (numBatches*batchSize);
1166a319912fSMatthew G. Knepley     Ne        = numChunks*numBatches*batchSize;
1167a319912fSMatthew G. Knepley     Nr        = numCells % (numBatches*batchSize);
1168a319912fSMatthew G. Knepley     offset    = numCells - Nr;
1169754551f4SMatthew G. Knepley     for (fieldJ = 0; fieldJ < Nf; ++fieldJ) {
1170754551f4SMatthew G. Knepley       void   (*g0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g0Funcs[fieldI*Nf+fieldJ];
1171754551f4SMatthew G. Knepley       void   (*g1)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g1Funcs[fieldI*Nf+fieldJ];
1172754551f4SMatthew G. Knepley       void   (*g2)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g2Funcs[fieldI*Nf+fieldJ];
1173754551f4SMatthew G. Knepley       void   (*g3)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g3Funcs[fieldI*Nf+fieldJ];
1174754551f4SMatthew G. Knepley 
1175a319912fSMatthew G. Knepley       geom.v0   = v0;
1176a319912fSMatthew G. Knepley       geom.J    = J;
1177a319912fSMatthew G. Knepley       geom.invJ = invJ;
1178a319912fSMatthew G. Knepley       geom.detJ = detJ;
1179754551f4SMatthew G. Knepley       ierr = PetscFEIntegrateJacobian(fe[fieldI], Ne, Nf, fe, fieldI, fieldJ, geom, u, NfAux, feAux, a, g0, g1, g2, g3, elemMat);CHKERRQ(ierr);
1180a319912fSMatthew G. Knepley       geom.v0   = &v0[offset*dim];
1181a319912fSMatthew G. Knepley       geom.J    = &J[offset*dim*dim];
1182a319912fSMatthew G. Knepley       geom.invJ = &invJ[offset*dim*dim];
1183a319912fSMatthew G. Knepley       geom.detJ = &detJ[offset];
1184754551f4SMatthew 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);
1185cb1e1211SMatthew G Knepley     }
1186cb1e1211SMatthew G Knepley   }
1187cb1e1211SMatthew G Knepley   for (c = cStart; c < cEnd; ++c) {
1188a319912fSMatthew G. Knepley     if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(c, name, cellDof, cellDof, &elemMat[c*cellDof*cellDof]);CHKERRQ(ierr);}
1189a319912fSMatthew G. Knepley     ierr = DMPlexMatSetClosure(dm, section, globalSection, JacP, c, &elemMat[c*cellDof*cellDof], ADD_VALUES);CHKERRQ(ierr);
1190cb1e1211SMatthew G Knepley   }
1191cb1e1211SMatthew G Knepley   ierr = PetscFree6(u,v0,J,invJ,detJ,elemMat);CHKERRQ(ierr);
1192754551f4SMatthew G. Knepley   if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);}
1193cb1e1211SMatthew G Knepley   ierr = MatAssemblyBegin(JacP, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1194cb1e1211SMatthew G Knepley   ierr = MatAssemblyEnd(JacP, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1195cb1e1211SMatthew G Knepley   if (mesh->printFEM) {
1196a319912fSMatthew G. Knepley     ierr = PetscPrintf(PETSC_COMM_WORLD, "%s:\n", name);CHKERRQ(ierr);
1197cb1e1211SMatthew G Knepley     ierr = MatChop(JacP, 1.0e-10);CHKERRQ(ierr);
1198cb1e1211SMatthew G Knepley     ierr = MatView(JacP, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1199cb1e1211SMatthew G Knepley   }
1200a319912fSMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_JacobianFEM,dm,0,0,0);CHKERRQ(ierr);
1201cb1e1211SMatthew G Knepley   ierr = PetscObjectTypeCompare((PetscObject) Jac, MATSHELL, &isShell);CHKERRQ(ierr);
1202cb1e1211SMatthew G Knepley   if (isShell) {
1203cb1e1211SMatthew G Knepley     JacActionCtx *jctx;
1204cb1e1211SMatthew G Knepley 
1205cb1e1211SMatthew G Knepley     ierr = MatShellGetContext(Jac, &jctx);CHKERRQ(ierr);
1206cb1e1211SMatthew G Knepley     ierr = VecCopy(X, jctx->u);CHKERRQ(ierr);
1207cb1e1211SMatthew G Knepley   }
1208cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
1209cb1e1211SMatthew G Knepley }
1210bceba477SMatthew G. Knepley 
1211d69c5d34SMatthew G. Knepley #if 0
1212d69c5d34SMatthew G. Knepley 
121332356553SMatthew 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[])
121432356553SMatthew G. Knepley {
121532356553SMatthew G. Knepley   g3[0] = 1.0;
121632356553SMatthew G. Knepley }
121732356553SMatthew G. Knepley 
121832356553SMatthew 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[])
121932356553SMatthew G. Knepley {
122032356553SMatthew G. Knepley   g3[0] = g3[3] = 1.0;
122132356553SMatthew G. Knepley }
122232356553SMatthew G. Knepley 
122332356553SMatthew 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[])
122432356553SMatthew G. Knepley {
122532356553SMatthew G. Knepley   g3[0] = g3[4] = g3[8] = 1.0;
122632356553SMatthew G. Knepley }
122732356553SMatthew G. Knepley 
1228bceba477SMatthew G. Knepley #undef __FUNCT__
1229d69c5d34SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeInterpolatorFEMBroken"
1230d69c5d34SMatthew G. Knepley PetscErrorCode DMPlexComputeInterpolatorFEMBroken(DM dmc, DM dmf, Mat I, void *user)
1231bceba477SMatthew G. Knepley {
1232bceba477SMatthew G. Knepley   DM_Plex          *mesh  = (DM_Plex *) dmc->data;
1233bceba477SMatthew G. Knepley   PetscFEM         *fem   = (PetscFEM *) user;
1234bceba477SMatthew G. Knepley   PetscFE          *fe    = fem->fe;
1235bceba477SMatthew G. Knepley   const char       *name  = "Interpolator";
123632356553SMatthew G. Knepley   PetscFE          *feRef;
1237bceba477SMatthew G. Knepley   PetscQuadrature   quad, quadOld;
1238bceba477SMatthew G. Knepley   PetscCellGeometry geom;
123932356553SMatthew G. Knepley   PetscSection      fsection, fglobalSection;
124032356553SMatthew G. Knepley   PetscSection      csection, cglobalSection;
1241bceba477SMatthew G. Knepley   PetscReal        *v0, *J, *invJ, *detJ;
1242bceba477SMatthew G. Knepley   PetscScalar      *elemMat;
1243bceba477SMatthew G. Knepley   PetscInt          dim, Nf, f, fieldI, fieldJ, numCells, cStart, cEnd, c;
124432356553SMatthew G. Knepley   PetscInt          rCellDof = 0, cCellDof = 0, numComponents = 0;
1245bceba477SMatthew G. Knepley   PetscErrorCode    ierr;
1246bceba477SMatthew G. Knepley 
1247bceba477SMatthew G. Knepley   PetscFunctionBegin;
1248bceba477SMatthew G. Knepley #if 0
1249bceba477SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1250bceba477SMatthew G. Knepley #endif
1251bceba477SMatthew G. Knepley   ierr = DMPlexGetDimension(dmf, &dim);CHKERRQ(ierr);
125232356553SMatthew G. Knepley   ierr = DMGetDefaultSection(dmf, &fsection);CHKERRQ(ierr);
125332356553SMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dmf, &fglobalSection);CHKERRQ(ierr);
125432356553SMatthew G. Knepley   ierr = DMGetDefaultSection(dmc, &csection);CHKERRQ(ierr);
125532356553SMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dmc, &cglobalSection);CHKERRQ(ierr);
125632356553SMatthew G. Knepley   ierr = PetscSectionGetNumFields(fsection, &Nf);CHKERRQ(ierr);
1257bceba477SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dmc, 0, &cStart, &cEnd);CHKERRQ(ierr);
1258bceba477SMatthew G. Knepley   numCells = cEnd - cStart;
125932356553SMatthew G. Knepley   ierr = PetscMalloc1(Nf,&feRef);CHKERRQ(ierr);
126032356553SMatthew G. Knepley   for (fieldI = 0; fieldI < Nf; ++fieldI) {
126132356553SMatthew G. Knepley     ierr = PetscFERefine(fe[fieldI], &feRef[fieldI]);CHKERRQ(ierr);
126232356553SMatthew G. Knepley   }
1263bceba477SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
126432356553SMatthew G. Knepley     PetscInt rNb, cNb, Nc;
1265bceba477SMatthew G. Knepley 
126632356553SMatthew G. Knepley     ierr = PetscFEGetDimension(feRef[f], &rNb);CHKERRQ(ierr);
126732356553SMatthew G. Knepley     ierr = PetscFEGetDimension(fe[f], &cNb);CHKERRQ(ierr);
1268bceba477SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe[f], &Nc);CHKERRQ(ierr);
1269bceba477SMatthew G. Knepley     numComponents += Nc;
127032356553SMatthew G. Knepley     rCellDof += rNb*Nc;
127132356553SMatthew G. Knepley     cCellDof += cNb*Nc;
1272bceba477SMatthew G. Knepley   }
1273bceba477SMatthew G. Knepley   ierr = MatZeroEntries(I);CHKERRQ(ierr);
127432356553SMatthew G. Knepley   ierr = PetscMalloc5(numCells*dim,&v0,numCells*dim*dim,&J,numCells*dim*dim,&invJ,numCells,&detJ,numCells*rCellDof*cCellDof,&elemMat);CHKERRQ(ierr);
1275bceba477SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
1276bceba477SMatthew G. Knepley     ierr = DMPlexComputeCellGeometry(dmc, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr);
1277bceba477SMatthew G. Knepley     if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c);
1278bceba477SMatthew G. Knepley   }
127932356553SMatthew G. Knepley   ierr = PetscMemzero(elemMat, numCells*rCellDof*cCellDof * sizeof(PetscScalar));CHKERRQ(ierr);
1280bceba477SMatthew G. Knepley   for (fieldI = 0; fieldI < Nf; ++fieldI) {
1281bceba477SMatthew G. Knepley     PetscInt numQuadPoints, Nb;
1282bceba477SMatthew G. Knepley     /* Conforming batches */
1283bceba477SMatthew G. Knepley     PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
1284bceba477SMatthew G. Knepley     /* Remainder */
1285bceba477SMatthew G. Knepley     PetscInt Nr, offset;
1286bceba477SMatthew G. Knepley 
1287bceba477SMatthew G. Knepley     /* Make new fine FE which refines the ref cell and the quadrature rule */
128832356553SMatthew G. Knepley     ierr = PetscFEGetQuadrature(feRef[fieldI], &quad);CHKERRQ(ierr);
128932356553SMatthew G. Knepley     ierr = PetscFEGetDimension(feRef[fieldI], &Nb);CHKERRQ(ierr);
129032356553SMatthew G. Knepley     ierr = PetscFEGetTileSizes(feRef[fieldI], NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
1291bceba477SMatthew G. Knepley     ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr);
1292bceba477SMatthew G. Knepley     blockSize = Nb*numQuadPoints;
1293bceba477SMatthew G. Knepley     batchSize = numBlocks * blockSize;
129432356553SMatthew G. Knepley     ierr = PetscFESetTileSizes(feRef[fieldI], blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
1295bceba477SMatthew G. Knepley     numChunks = numCells / (numBatches*batchSize);
1296bceba477SMatthew G. Knepley     Ne        = numChunks*numBatches*batchSize;
1297bceba477SMatthew G. Knepley     Nr        = numCells % (numBatches*batchSize);
1298bceba477SMatthew G. Knepley     offset    = numCells - Nr;
1299d69c5d34SMatthew G. Knepley 
1300bceba477SMatthew G. Knepley     for (fieldJ = 0; fieldJ < Nf; ++fieldJ) {
130132356553SMatthew G. Knepley       /* void   (*g0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g0Funcs[fieldI*Nf+fieldJ]; */
130232356553SMatthew G. Knepley       void   (*g0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = g0_identity_2d_static;
1303bceba477SMatthew G. Knepley 
1304bceba477SMatthew G. Knepley       /* Replace quadrature in coarse FE with refined quadrature */
1305bceba477SMatthew G. Knepley       ierr = PetscFEGetQuadrature(fe[fieldJ], &quadOld);CHKERRQ(ierr);
1306bceba477SMatthew G. Knepley       ierr = PetscObjectReference((PetscObject) quadOld);CHKERRQ(ierr);
1307bceba477SMatthew G. Knepley       ierr = PetscFESetQuadrature(fe[fieldJ], quad);CHKERRQ(ierr);
1308bceba477SMatthew G. Knepley       geom.v0   = v0;
1309bceba477SMatthew G. Knepley       geom.J    = J;
1310bceba477SMatthew G. Knepley       geom.invJ = invJ;
1311bceba477SMatthew G. Knepley       geom.detJ = detJ;
131232356553SMatthew G. Knepley       ierr = PetscFEIntegrateInterpolator_Basic(feRef[fieldI], Ne, Nf, feRef, fieldI, fe, fieldJ, geom, g0, elemMat);CHKERRQ(ierr);
1313bceba477SMatthew G. Knepley       geom.v0   = &v0[offset*dim];
1314bceba477SMatthew G. Knepley       geom.J    = &J[offset*dim*dim];
1315bceba477SMatthew G. Knepley       geom.invJ = &invJ[offset*dim*dim];
1316bceba477SMatthew G. Knepley       geom.detJ = &detJ[offset];
131732356553SMatthew G. Knepley       ierr = PetscFEIntegrateInterpolator_Basic(feRef[fieldI], Nr, Nf, feRef, fieldI, fe, fieldJ, geom, g0, &elemMat[offset*rCellDof*cCellDof]);CHKERRQ(ierr);
1318bceba477SMatthew G. Knepley       ierr = PetscFESetQuadrature(fe[fieldJ], quadOld);CHKERRQ(ierr);
1319bceba477SMatthew G. Knepley     }
1320bceba477SMatthew G. Knepley   }
1321bceba477SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
132232356553SMatthew G. Knepley     if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(c, name, rCellDof, cCellDof, &elemMat[c*rCellDof*cCellDof]);CHKERRQ(ierr);}
132332356553SMatthew G. Knepley     ierr = DMPlexMatSetClosureRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, I, c, &elemMat[c*rCellDof*cCellDof], ADD_VALUES);CHKERRQ(ierr);
1324bceba477SMatthew G. Knepley   }
1325bceba477SMatthew G. Knepley   ierr = PetscFree5(v0,J,invJ,detJ,elemMat);CHKERRQ(ierr);
132632356553SMatthew G. Knepley   ierr = PetscFree(feRef);CHKERRQ(ierr);
1327bceba477SMatthew G. Knepley   ierr = MatAssemblyBegin(I, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1328bceba477SMatthew G. Knepley   ierr = MatAssemblyEnd(I, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1329bceba477SMatthew G. Knepley   if (mesh->printFEM) {
1330bceba477SMatthew G. Knepley     ierr = PetscPrintf(PETSC_COMM_WORLD, "%s:\n", name);CHKERRQ(ierr);
1331bceba477SMatthew G. Knepley     ierr = MatChop(I, 1.0e-10);CHKERRQ(ierr);
1332bceba477SMatthew G. Knepley     ierr = MatView(I, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1333bceba477SMatthew G. Knepley   }
1334bceba477SMatthew G. Knepley #if 0
1335bceba477SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1336bceba477SMatthew G. Knepley #endif
1337bceba477SMatthew G. Knepley   PetscFunctionReturn(0);
1338bceba477SMatthew G. Knepley }
1339d69c5d34SMatthew G. Knepley #endif
1340d69c5d34SMatthew G. Knepley 
1341d69c5d34SMatthew G. Knepley #undef __FUNCT__
1342d69c5d34SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeInterpolatorFEM"
1343d69c5d34SMatthew G. Knepley /*@
1344d69c5d34SMatthew G. Knepley   DMPlexComputeInterpolatorFEM - Form the local portion of the interpolation matrix I from the coarse DM to the uniformly refined DM.
1345d69c5d34SMatthew G. Knepley 
1346d69c5d34SMatthew G. Knepley   Input Parameters:
1347d69c5d34SMatthew G. Knepley + dmf  - The fine mesh
1348d69c5d34SMatthew G. Knepley . dmc  - The coarse mesh
1349d69c5d34SMatthew G. Knepley - user - The user context
1350d69c5d34SMatthew G. Knepley 
1351d69c5d34SMatthew G. Knepley   Output Parameter:
1352934789fcSMatthew G. Knepley . In  - The interpolation matrix
1353d69c5d34SMatthew G. Knepley 
1354d69c5d34SMatthew G. Knepley   Note:
1355d69c5d34SMatthew G. Knepley   The first member of the user context must be an FEMContext.
1356d69c5d34SMatthew G. Knepley 
1357d69c5d34SMatthew G. Knepley   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
1358d69c5d34SMatthew G. Knepley   like a GPU, or vectorize on a multicore machine.
1359d69c5d34SMatthew G. Knepley 
1360d69c5d34SMatthew G. Knepley   Level: developer
1361d69c5d34SMatthew G. Knepley 
1362d69c5d34SMatthew G. Knepley .seealso: DMPlexComputeJacobianFEM()
1363d69c5d34SMatthew G. Knepley @*/
1364934789fcSMatthew G. Knepley PetscErrorCode DMPlexComputeInterpolatorFEM(DM dmc, DM dmf, Mat In, void *user)
1365d69c5d34SMatthew G. Knepley {
1366d69c5d34SMatthew G. Knepley   DM_Plex          *mesh  = (DM_Plex *) dmc->data;
1367d69c5d34SMatthew G. Knepley   PetscFEM         *fem   = (PetscFEM *) user;
1368d69c5d34SMatthew G. Knepley   PetscFE          *fe    = fem->fe;
1369d69c5d34SMatthew G. Knepley   const char       *name  = "Interpolator";
1370d69c5d34SMatthew G. Knepley   PetscFE          *feRef;
1371d69c5d34SMatthew G. Knepley   PetscSection      fsection, fglobalSection;
1372d69c5d34SMatthew G. Knepley   PetscSection      csection, cglobalSection;
1373d69c5d34SMatthew G. Knepley   PetscScalar      *elemMat;
1374942a7a06SMatthew G. Knepley   PetscInt          dim, Nf, f, fieldI, fieldJ, offsetI, offsetJ, cStart, cEnd, c;
1375d69c5d34SMatthew G. Knepley   PetscInt          rCellDof = 0, cCellDof = 0;
1376d69c5d34SMatthew G. Knepley   PetscErrorCode    ierr;
1377d69c5d34SMatthew G. Knepley 
1378d69c5d34SMatthew G. Knepley   PetscFunctionBegin;
1379d69c5d34SMatthew G. Knepley #if 0
1380d69c5d34SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1381d69c5d34SMatthew G. Knepley #endif
1382d69c5d34SMatthew G. Knepley   ierr = DMPlexGetDimension(dmf, &dim);CHKERRQ(ierr);
1383d69c5d34SMatthew G. Knepley   ierr = DMGetDefaultSection(dmf, &fsection);CHKERRQ(ierr);
1384d69c5d34SMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dmf, &fglobalSection);CHKERRQ(ierr);
1385d69c5d34SMatthew G. Knepley   ierr = DMGetDefaultSection(dmc, &csection);CHKERRQ(ierr);
1386d69c5d34SMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dmc, &cglobalSection);CHKERRQ(ierr);
1387d69c5d34SMatthew G. Knepley   ierr = PetscSectionGetNumFields(fsection, &Nf);CHKERRQ(ierr);
1388d69c5d34SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dmc, 0, &cStart, &cEnd);CHKERRQ(ierr);
1389d69c5d34SMatthew G. Knepley   ierr = PetscMalloc1(Nf,&feRef);CHKERRQ(ierr);
1390d69c5d34SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
1391d69c5d34SMatthew G. Knepley     PetscInt rNb, cNb, Nc;
1392d69c5d34SMatthew G. Knepley 
1393d69c5d34SMatthew G. Knepley     ierr = PetscFERefine(fe[f], &feRef[f]);CHKERRQ(ierr);
1394d69c5d34SMatthew G. Knepley     ierr = PetscFEGetDimension(feRef[f], &rNb);CHKERRQ(ierr);
1395d69c5d34SMatthew G. Knepley     ierr = PetscFEGetDimension(fe[f], &cNb);CHKERRQ(ierr);
1396d69c5d34SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe[f], &Nc);CHKERRQ(ierr);
1397d69c5d34SMatthew G. Knepley     rCellDof += rNb*Nc;
1398d69c5d34SMatthew G. Knepley     cCellDof += cNb*Nc;
1399d69c5d34SMatthew G. Knepley   }
1400934789fcSMatthew G. Knepley   ierr = MatZeroEntries(In);CHKERRQ(ierr);
1401d69c5d34SMatthew G. Knepley   ierr = PetscMalloc1(rCellDof*cCellDof,&elemMat);CHKERRQ(ierr);
1402d69c5d34SMatthew G. Knepley   ierr = PetscMemzero(elemMat, rCellDof*cCellDof * sizeof(PetscScalar));CHKERRQ(ierr);
1403d69c5d34SMatthew G. Knepley   for (fieldI = 0, offsetI = 0; fieldI < Nf; ++fieldI) {
1404d69c5d34SMatthew G. Knepley     PetscDualSpace   Qref;
1405d69c5d34SMatthew G. Knepley     PetscQuadrature  f;
1406d69c5d34SMatthew G. Knepley     const PetscReal *qpoints, *qweights;
1407d69c5d34SMatthew G. Knepley     PetscReal       *points;
1408d69c5d34SMatthew G. Knepley     PetscInt         npoints = 0, Nc, Np, fpdim, i, k, p, d;
1409d69c5d34SMatthew G. Knepley 
1410d69c5d34SMatthew G. Knepley     /* Compose points from all dual basis functionals */
1411d69c5d34SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe[fieldI], &Nc);CHKERRQ(ierr);
1412d69c5d34SMatthew G. Knepley     ierr = PetscFEGetDualSpace(feRef[fieldI], &Qref);CHKERRQ(ierr);
1413d69c5d34SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(Qref, &fpdim);CHKERRQ(ierr);
1414d69c5d34SMatthew G. Knepley     for (i = 0; i < fpdim; ++i) {
1415d69c5d34SMatthew G. Knepley       ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1416d69c5d34SMatthew G. Knepley       ierr = PetscQuadratureGetData(f, NULL, &Np, NULL, NULL);CHKERRQ(ierr);
1417d69c5d34SMatthew G. Knepley       npoints += Np;
1418d69c5d34SMatthew G. Knepley     }
1419d69c5d34SMatthew G. Knepley     ierr = PetscMalloc1(npoints*dim,&points);CHKERRQ(ierr);
1420d69c5d34SMatthew G. Knepley     for (i = 0, k = 0; i < fpdim; ++i) {
1421d69c5d34SMatthew G. Knepley       ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1422d69c5d34SMatthew G. Knepley       ierr = PetscQuadratureGetData(f, NULL, &Np, &qpoints, NULL);CHKERRQ(ierr);
1423d69c5d34SMatthew G. Knepley       for (p = 0; p < Np; ++p, ++k) for (d = 0; d < dim; ++d) points[k*dim+d] = qpoints[p*dim+d];
1424d69c5d34SMatthew G. Knepley     }
1425d69c5d34SMatthew G. Knepley 
1426d69c5d34SMatthew G. Knepley     for (fieldJ = 0, offsetJ = 0; fieldJ < Nf; ++fieldJ) {
1427d69c5d34SMatthew G. Knepley       PetscReal *B;
142836a6d9c0SMatthew G. Knepley       PetscInt   NcJ, cpdim, j;
1429d69c5d34SMatthew G. Knepley 
1430d69c5d34SMatthew G. Knepley       /* Evaluate basis at points */
143136a6d9c0SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe[fieldJ], &NcJ);CHKERRQ(ierr);
143236a6d9c0SMatthew 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);
1433d69c5d34SMatthew G. Knepley       ierr = PetscFEGetDimension(fe[fieldJ], &cpdim);CHKERRQ(ierr);
1434ffe73a53SMatthew G. Knepley       /* For now, fields only interpolate themselves */
1435ffe73a53SMatthew G. Knepley       if (fieldI == fieldJ) {
1436d69c5d34SMatthew G. Knepley         ierr = PetscFEGetTabulation(fe[fieldJ], npoints, points, &B, NULL, NULL);CHKERRQ(ierr);
1437d69c5d34SMatthew G. Knepley         for (i = 0, k = 0; i < fpdim; ++i) {
1438d69c5d34SMatthew G. Knepley           ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1439d69c5d34SMatthew G. Knepley           ierr = PetscQuadratureGetData(f, NULL, &Np, NULL, &qweights);CHKERRQ(ierr);
1440d69c5d34SMatthew G. Knepley           for (p = 0; p < Np; ++p, ++k) {
144136a6d9c0SMatthew G. Knepley             for (j = 0; j < cpdim; ++j) {
144236a6d9c0SMatthew 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];
144336a6d9c0SMatthew G. Knepley             }
1444d69c5d34SMatthew G. Knepley           }
1445d69c5d34SMatthew G. Knepley         }
1446d69c5d34SMatthew G. Knepley         ierr = PetscFERestoreTabulation(fe[fieldJ], npoints, points, &B, NULL, NULL);CHKERRQ(ierr);CHKERRQ(ierr);
1447ffe73a53SMatthew G. Knepley       }
144836a6d9c0SMatthew G. Knepley       offsetJ += cpdim*NcJ;
1449d69c5d34SMatthew G. Knepley     }
1450d69c5d34SMatthew G. Knepley     offsetI += fpdim*Nc;
1451549a8adaSMatthew G. Knepley     ierr = PetscFree(points);CHKERRQ(ierr);
1452d69c5d34SMatthew G. Knepley   }
1453d69c5d34SMatthew G. Knepley   if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(0, name, rCellDof, cCellDof, elemMat);CHKERRQ(ierr);}
1454d69c5d34SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
1455934789fcSMatthew G. Knepley     ierr = DMPlexMatSetClosureRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, In, c, elemMat, INSERT_VALUES);CHKERRQ(ierr);
1456d69c5d34SMatthew G. Knepley   }
1457549a8adaSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {ierr = PetscFEDestroy(&feRef[f]);CHKERRQ(ierr);}
1458d69c5d34SMatthew G. Knepley   ierr = PetscFree(feRef);CHKERRQ(ierr);
1459549a8adaSMatthew G. Knepley   ierr = PetscFree(elemMat);CHKERRQ(ierr);
1460934789fcSMatthew G. Knepley   ierr = MatAssemblyBegin(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1461934789fcSMatthew G. Knepley   ierr = MatAssemblyEnd(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1462d69c5d34SMatthew G. Knepley   if (mesh->printFEM) {
1463d69c5d34SMatthew G. Knepley     ierr = PetscPrintf(PETSC_COMM_WORLD, "%s:\n", name);CHKERRQ(ierr);
1464934789fcSMatthew G. Knepley     ierr = MatChop(In, 1.0e-10);CHKERRQ(ierr);
1465934789fcSMatthew G. Knepley     ierr = MatView(In, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1466d69c5d34SMatthew G. Knepley   }
1467d69c5d34SMatthew G. Knepley #if 0
1468d69c5d34SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1469d69c5d34SMatthew G. Knepley #endif
1470d69c5d34SMatthew G. Knepley   PetscFunctionReturn(0);
1471d69c5d34SMatthew G. Knepley }
14726c73c22cSMatthew G. Knepley 
14736c73c22cSMatthew G. Knepley #undef __FUNCT__
14746c73c22cSMatthew G. Knepley #define __FUNCT__ "DMPlexAddBoundary"
14756c73c22cSMatthew G. Knepley /* The ids can be overridden by the command line option -bc_<boundary name> */
14766c73c22cSMatthew G. Knepley PetscErrorCode DMPlexAddBoundary(DM dm, PetscBool isEssential, const char name[], PetscInt field, void (*bcFunc)(), PetscInt numids, const PetscInt *ids, void *ctx)
14776c73c22cSMatthew G. Knepley {
14786c73c22cSMatthew G. Knepley   DM_Plex       *mesh = (DM_Plex *) dm->data;
14796c73c22cSMatthew G. Knepley   DMBoundary     b;
14806c73c22cSMatthew G. Knepley   PetscErrorCode ierr;
14816c73c22cSMatthew G. Knepley 
14826c73c22cSMatthew G. Knepley   PetscFunctionBegin;
14836c73c22cSMatthew G. Knepley   ierr = PetscNew(&b);CHKERRQ(ierr);
14846c73c22cSMatthew G. Knepley   ierr = PetscStrallocpy(name, (char **) &b->name);CHKERRQ(ierr);
14856c73c22cSMatthew G. Knepley   ierr = PetscMalloc1(numids, &b->ids);CHKERRQ(ierr);
14866c73c22cSMatthew G. Knepley   ierr = PetscMemcpy(b->ids, ids, numids*sizeof(PetscInt));CHKERRQ(ierr);
14876c73c22cSMatthew G. Knepley   b->essential   = isEssential;
14886c73c22cSMatthew G. Knepley   b->field       = field;
14896c73c22cSMatthew G. Knepley   b->func        = bcFunc;
14906c73c22cSMatthew G. Knepley   b->numids      = numids;
14916c73c22cSMatthew G. Knepley   b->ctx         = ctx;
14926c73c22cSMatthew G. Knepley   b->next        = mesh->boundary;
14936c73c22cSMatthew G. Knepley   mesh->boundary = b;
14946c73c22cSMatthew G. Knepley   PetscFunctionReturn(0);
14956c73c22cSMatthew G. Knepley }
14966c73c22cSMatthew G. Knepley 
14976c73c22cSMatthew G. Knepley #undef __FUNCT__
14986c73c22cSMatthew G. Knepley #define __FUNCT__ "DMPlexGetNumBoundary"
14996c73c22cSMatthew G. Knepley PetscErrorCode DMPlexGetNumBoundary(DM dm, PetscInt *numBd)
15006c73c22cSMatthew G. Knepley {
15016c73c22cSMatthew G. Knepley   DM_Plex   *mesh = (DM_Plex *) dm->data;
15026c73c22cSMatthew G. Knepley   DMBoundary b    = mesh->boundary;
15036c73c22cSMatthew G. Knepley 
15046c73c22cSMatthew G. Knepley   PetscFunctionBegin;
15056c73c22cSMatthew G. Knepley   *numBd = 0;
15066c73c22cSMatthew G. Knepley   while (b) {++(*numBd); b = b->next;}
15076c73c22cSMatthew G. Knepley   PetscFunctionReturn(0);
15086c73c22cSMatthew G. Knepley }
15096c73c22cSMatthew G. Knepley 
15106c73c22cSMatthew G. Knepley #undef __FUNCT__
15116c73c22cSMatthew G. Knepley #define __FUNCT__ "DMPlexGetBoundary"
15126c73c22cSMatthew G. Knepley PetscErrorCode DMPlexGetBoundary(DM dm, PetscInt bd, PetscBool *isEssential, const char **name, PetscInt *field, void (**func)(), PetscInt *numids, const PetscInt **ids, void **ctx)
15136c73c22cSMatthew G. Knepley {
15146c73c22cSMatthew G. Knepley   DM_Plex   *mesh = (DM_Plex *) dm->data;
15156c73c22cSMatthew G. Knepley   DMBoundary b    = mesh->boundary;
15166c73c22cSMatthew G. Knepley   PetscInt   n    = 0;
15176c73c22cSMatthew G. Knepley 
15186c73c22cSMatthew G. Knepley   PetscFunctionBegin;
15196c73c22cSMatthew G. Knepley   while (b) {
15206c73c22cSMatthew G. Knepley     if (n == bd) break;
15216c73c22cSMatthew G. Knepley     b = b->next;
15226c73c22cSMatthew G. Knepley     ++n;
15236c73c22cSMatthew G. Knepley   }
15246c73c22cSMatthew G. Knepley   if (n != bd) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Boundary %d is not in [0, %d)", bd, n);
15256c73c22cSMatthew G. Knepley   if (isEssential) {
15266c73c22cSMatthew G. Knepley     PetscValidPointer(isEssential, 3);
15276c73c22cSMatthew G. Knepley     *isEssential = b->essential;
15286c73c22cSMatthew G. Knepley   }
15296c73c22cSMatthew G. Knepley   if (name) {
15306c73c22cSMatthew G. Knepley     PetscValidPointer(name, 4);
15316c73c22cSMatthew G. Knepley     *name = b->name;
15326c73c22cSMatthew G. Knepley   }
15336c73c22cSMatthew G. Knepley   if (field) {
15346c73c22cSMatthew G. Knepley     PetscValidPointer(field, 5);
15356c73c22cSMatthew G. Knepley     *field = b->field;
15366c73c22cSMatthew G. Knepley   }
15376c73c22cSMatthew G. Knepley   if (func) {
15386c73c22cSMatthew G. Knepley     PetscValidPointer(func, 6);
15396c73c22cSMatthew G. Knepley     *func = b->func;
15406c73c22cSMatthew G. Knepley   }
15416c73c22cSMatthew G. Knepley   if (numids) {
15426c73c22cSMatthew G. Knepley     PetscValidPointer(numids, 7);
15436c73c22cSMatthew G. Knepley     *numids = b->numids;
15446c73c22cSMatthew G. Knepley   }
15456c73c22cSMatthew G. Knepley   if (ids) {
15466c73c22cSMatthew G. Knepley     PetscValidPointer(ids, 8);
15476c73c22cSMatthew G. Knepley     *ids = b->ids;
15486c73c22cSMatthew G. Knepley   }
15496c73c22cSMatthew G. Knepley   if (ctx) {
15506c73c22cSMatthew G. Knepley     PetscValidPointer(ctx, 9);
15516c73c22cSMatthew G. Knepley     *ctx = b->ctx;
15526c73c22cSMatthew G. Knepley   }
15536c73c22cSMatthew G. Knepley   PetscFunctionReturn(0);
15546c73c22cSMatthew G. Knepley }
1555