xref: /petsc/src/dm/impls/plex/plexfem.c (revision c110b1eed3fc6d435822e10a4aa5917ffc8efa18)
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"
185*c110b1eeSGeoffrey 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;
19172f94c41SMatthew G. Knepley   PetscInt        numFields, numComp, dim, spDim, totDim = 0, numValues, cStart, cEnd, c, f, d, v;
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) {
218*c110b1eeSGeoffrey 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) {
222*c110b1eeSGeoffrey Irving         ierr = PetscDualSpaceApply(sp[f], d, geom, numComp, funcs[f], ctx, &values[v]);CHKERRQ(ierr);
22372f94c41SMatthew G. Knepley         v += numComp;
224cb1e1211SMatthew G Knepley       }
225cb1e1211SMatthew G Knepley     }
22672f94c41SMatthew G. Knepley     ierr = DMPlexVecSetClosure(dm, section, localX, c, values, mode);CHKERRQ(ierr);
227cb1e1211SMatthew G Knepley   }
22872f94c41SMatthew G. Knepley   ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
2291f2da991SMatthew G. Knepley   ierr = PetscFree2(v0,J);CHKERRQ(ierr);
23072f94c41SMatthew G. Knepley   ierr = PetscFree(sp);CHKERRQ(ierr);
231cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
232cb1e1211SMatthew G Knepley }
233cb1e1211SMatthew G Knepley 
234cb1e1211SMatthew G Knepley #undef __FUNCT__
235cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexProjectFunction"
236cb1e1211SMatthew G Knepley /*@C
237cb1e1211SMatthew G Knepley   DMPlexProjectFunction - This projects the given function into the function space provided.
238cb1e1211SMatthew G Knepley 
239cb1e1211SMatthew G Knepley   Input Parameters:
240cb1e1211SMatthew G Knepley + dm      - The DM
24172f94c41SMatthew G. Knepley . fe      - The PetscFE associated with the field
24272f94c41SMatthew G. Knepley . funcs   - The coordinate functions to evaluate, one per field
243*c110b1eeSGeoffrey Irving . ctxs    - Optional array of contexts to pass to each coordinate function.  ctxs itself may be null.
244cb1e1211SMatthew G Knepley - mode    - The insertion mode for values
245cb1e1211SMatthew G Knepley 
246cb1e1211SMatthew G Knepley   Output Parameter:
247cb1e1211SMatthew G Knepley . X - vector
248cb1e1211SMatthew G Knepley 
249cb1e1211SMatthew G Knepley   Level: developer
250cb1e1211SMatthew G Knepley 
251878cb397SSatish Balay .seealso: DMPlexComputeL2Diff()
252878cb397SSatish Balay @*/
253*c110b1eeSGeoffrey Irving PetscErrorCode DMPlexProjectFunction(DM dm, PetscFE fe[], void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, InsertMode mode, Vec X)
254cb1e1211SMatthew G Knepley {
255cb1e1211SMatthew G Knepley   Vec            localX;
256cb1e1211SMatthew G Knepley   PetscErrorCode ierr;
257cb1e1211SMatthew G Knepley 
258cb1e1211SMatthew G Knepley   PetscFunctionBegin;
2599a800dd8SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
260cb1e1211SMatthew G Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
261*c110b1eeSGeoffrey Irving   ierr = DMPlexProjectFunctionLocal(dm, fe, funcs, ctxs, mode, localX);CHKERRQ(ierr);
262cb1e1211SMatthew G Knepley   ierr = DMLocalToGlobalBegin(dm, localX, mode, X);CHKERRQ(ierr);
263cb1e1211SMatthew G Knepley   ierr = DMLocalToGlobalEnd(dm, localX, mode, X);CHKERRQ(ierr);
264cb1e1211SMatthew G Knepley   ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
265cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
266cb1e1211SMatthew G Knepley }
267cb1e1211SMatthew G Knepley 
268cb1e1211SMatthew G Knepley #undef __FUNCT__
269cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexComputeL2Diff"
270cb1e1211SMatthew G Knepley /*@C
271cb1e1211SMatthew G Knepley   DMPlexComputeL2Diff - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h.
272cb1e1211SMatthew G Knepley 
273cb1e1211SMatthew G Knepley   Input Parameters:
274cb1e1211SMatthew G Knepley + dm    - The DM
275c5bbbd5bSMatthew G. Knepley . fe    - The PetscFE object for each field
276cb1e1211SMatthew G Knepley . funcs - The functions to evaluate for each field component
277*c110b1eeSGeoffrey Irving . ctxs  - Optional array of contexts to pass to each function.  ctxs itself may be null.
278cb1e1211SMatthew G Knepley - X     - The coefficient vector u_h
279cb1e1211SMatthew G Knepley 
280cb1e1211SMatthew G Knepley   Output Parameter:
281cb1e1211SMatthew G Knepley . diff - The diff ||u - u_h||_2
282cb1e1211SMatthew G Knepley 
283cb1e1211SMatthew G Knepley   Level: developer
284cb1e1211SMatthew G Knepley 
285cb1e1211SMatthew G Knepley .seealso: DMPlexProjectFunction()
286878cb397SSatish Balay @*/
287*c110b1eeSGeoffrey Irving PetscErrorCode DMPlexComputeL2Diff(DM dm, PetscFE fe[], void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, PetscReal *diff)
288cb1e1211SMatthew G Knepley {
289cb1e1211SMatthew G Knepley   const PetscInt  debug = 0;
290cb1e1211SMatthew G Knepley   PetscSection    section;
291c5bbbd5bSMatthew G. Knepley   PetscQuadrature quad;
292cb1e1211SMatthew G Knepley   Vec             localX;
29372f94c41SMatthew G. Knepley   PetscScalar    *funcVal;
294cb1e1211SMatthew G Knepley   PetscReal      *coords, *v0, *J, *invJ, detJ;
295cb1e1211SMatthew G Knepley   PetscReal       localDiff = 0.0;
296cb1e1211SMatthew G Knepley   PetscInt        dim, numFields, numComponents = 0, cStart, cEnd, c, field, fieldOffset, comp;
297cb1e1211SMatthew G Knepley   PetscErrorCode  ierr;
298cb1e1211SMatthew G Knepley 
299cb1e1211SMatthew G Knepley   PetscFunctionBegin;
300cb1e1211SMatthew G Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
301cb1e1211SMatthew G Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
302cb1e1211SMatthew G Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
303cb1e1211SMatthew G Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
304cb1e1211SMatthew G Knepley   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
305cb1e1211SMatthew G Knepley   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
306cb1e1211SMatthew G Knepley   for (field = 0; field < numFields; ++field) {
307c5bbbd5bSMatthew G. Knepley     PetscInt Nc;
308c5bbbd5bSMatthew G. Knepley 
309c5bbbd5bSMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe[field], &Nc);CHKERRQ(ierr);
310c5bbbd5bSMatthew G. Knepley     numComponents += Nc;
311cb1e1211SMatthew G Knepley   }
312*c110b1eeSGeoffrey Irving   ierr = DMPlexProjectFunctionLocal(dm, fe, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr);
313dcca6d9dSJed Brown   ierr = PetscMalloc5(numComponents,&funcVal,dim,&coords,dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr);
314cb1e1211SMatthew G Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
315c5bbbd5bSMatthew G. Knepley   ierr = PetscFEGetQuadrature(fe[0], &quad);CHKERRQ(ierr);
316cb1e1211SMatthew G Knepley   for (c = cStart; c < cEnd; ++c) {
317a1e44745SMatthew G. Knepley     PetscScalar *x = NULL;
318cb1e1211SMatthew G Knepley     PetscReal    elemDiff = 0.0;
319cb1e1211SMatthew G Knepley 
320cb1e1211SMatthew G Knepley     ierr = DMPlexComputeCellGeometry(dm, c, v0, J, invJ, &detJ);CHKERRQ(ierr);
321cb1e1211SMatthew G Knepley     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
322cb1e1211SMatthew G Knepley     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
323cb1e1211SMatthew G Knepley 
324cb1e1211SMatthew G Knepley     for (field = 0, comp = 0, fieldOffset = 0; field < numFields; ++field) {
325*c110b1eeSGeoffrey Irving       void * const     ctx           = ctxs ? ctxs[field] : NULL;
326f9fd7fdbSMatthew G. Knepley       const PetscInt   numQuadPoints = quad.numPoints;
327f9fd7fdbSMatthew G. Knepley       const PetscReal *quadPoints    = quad.points;
328f9fd7fdbSMatthew G. Knepley       const PetscReal *quadWeights   = quad.weights;
329c5bbbd5bSMatthew G. Knepley       PetscReal       *basis;
330c5bbbd5bSMatthew G. Knepley       PetscInt         numBasisFuncs, numBasisComps, q, d, e, fc, f;
331cb1e1211SMatthew G Knepley 
332c5bbbd5bSMatthew G. Knepley       ierr = PetscFEGetDimension(fe[field], &numBasisFuncs);CHKERRQ(ierr);
333c5bbbd5bSMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe[field], &numBasisComps);CHKERRQ(ierr);
334c5bbbd5bSMatthew G. Knepley       ierr = PetscFEGetDefaultTabulation(fe[field], &basis, NULL, NULL);CHKERRQ(ierr);
335cb1e1211SMatthew G Knepley       if (debug) {
336cb1e1211SMatthew G Knepley         char title[1024];
337cb1e1211SMatthew G Knepley         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
338cb1e1211SMatthew G Knepley         ierr = DMPrintCellVector(c, title, numBasisFuncs*numBasisComps, &x[fieldOffset]);CHKERRQ(ierr);
339cb1e1211SMatthew G Knepley       }
340cb1e1211SMatthew G Knepley       for (q = 0; q < numQuadPoints; ++q) {
341cb1e1211SMatthew G Knepley         for (d = 0; d < dim; d++) {
342cb1e1211SMatthew G Knepley           coords[d] = v0[d];
343cb1e1211SMatthew G Knepley           for (e = 0; e < dim; e++) {
344cb1e1211SMatthew G Knepley             coords[d] += J[d*dim+e]*(quadPoints[q*dim+e] + 1.0);
345cb1e1211SMatthew G Knepley           }
346cb1e1211SMatthew G Knepley         }
347*c110b1eeSGeoffrey Irving         (*funcs[field])(coords, funcVal, ctx);
348cb1e1211SMatthew G Knepley         for (fc = 0; fc < numBasisComps; ++fc) {
349a1d24da5SMatthew G. Knepley           PetscScalar interpolant = 0.0;
350a1d24da5SMatthew G. Knepley 
351cb1e1211SMatthew G Knepley           for (f = 0; f < numBasisFuncs; ++f) {
352cb1e1211SMatthew G Knepley             const PetscInt fidx = f*numBasisComps+fc;
353a1d24da5SMatthew G. Knepley             interpolant += x[fieldOffset+fidx]*basis[q*numBasisFuncs*numBasisComps+fidx];
354cb1e1211SMatthew G Knepley           }
35572f94c41SMatthew 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);}
35672f94c41SMatthew G. Knepley           elemDiff += PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ;
357cb1e1211SMatthew G Knepley         }
358cb1e1211SMatthew G Knepley       }
359cb1e1211SMatthew G Knepley       comp        += numBasisComps;
360cb1e1211SMatthew G Knepley       fieldOffset += numBasisFuncs*numBasisComps;
361cb1e1211SMatthew G Knepley     }
362cb1e1211SMatthew G Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
363cb1e1211SMatthew G Knepley     if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);}
364cb1e1211SMatthew G Knepley     localDiff += elemDiff;
365cb1e1211SMatthew G Knepley   }
36672f94c41SMatthew G. Knepley   ierr  = PetscFree5(funcVal,coords,v0,J,invJ);CHKERRQ(ierr);
367cb1e1211SMatthew G Knepley   ierr  = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
36886a74ee0SMatthew G. Knepley   ierr  = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
369cb1e1211SMatthew G Knepley   *diff = PetscSqrtReal(*diff);
370cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
371cb1e1211SMatthew G Knepley }
372cb1e1211SMatthew G Knepley 
373a0845e3aSMatthew G. Knepley #if 0
374a0845e3aSMatthew G. Knepley 
375cb1e1211SMatthew G Knepley #undef __FUNCT__
376cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexComputeResidualFEM"
377cb1e1211SMatthew G Knepley PetscErrorCode DMPlexComputeResidualFEM(DM dm, Vec X, Vec F, void *user)
378cb1e1211SMatthew G Knepley {
379cb1e1211SMatthew G Knepley   DM_Plex         *mesh   = (DM_Plex *) dm->data;
3809a559087SMatthew G. Knepley   PetscFEM        *fem    = (PetscFEM *) user;
381cb1e1211SMatthew G Knepley   PetscQuadrature *quad   = fem->quad;
382652b88e8SMatthew G. Knepley   PetscQuadrature *quadBd = fem->quadBd;
383cb1e1211SMatthew G Knepley   PetscSection     section;
38402a80efeSMatthew G. Knepley   PetscReal       *v0, *n, *J, *invJ, *detJ;
385cb1e1211SMatthew G Knepley   PetscScalar     *elemVec, *u;
386cb1e1211SMatthew G Knepley   PetscInt         dim, numFields, field, numBatchesTmp = 1, numCells, cStart, cEnd, c;
387652b88e8SMatthew G. Knepley   PetscInt         cellDof, numComponents;
388652b88e8SMatthew G. Knepley   PetscBool        has;
389cb1e1211SMatthew G Knepley   PetscErrorCode   ierr;
390cb1e1211SMatthew G Knepley 
391cb1e1211SMatthew G Knepley   PetscFunctionBegin;
392652b88e8SMatthew G. Knepley   if (has && quadBd) {
393652b88e8SMatthew G. Knepley     DMLabel         label;
394652b88e8SMatthew G. Knepley     IS              pointIS;
395652b88e8SMatthew G. Knepley     const PetscInt *points;
396652b88e8SMatthew G. Knepley     PetscInt        numPoints, p;
397652b88e8SMatthew G. Knepley 
398652b88e8SMatthew G. Knepley     ierr = DMPlexGetLabel(dm, "boundary", &label);CHKERRQ(ierr);
399652b88e8SMatthew G. Knepley     ierr = DMLabelGetStratumSize(label, 1, &numPoints);CHKERRQ(ierr);
400652b88e8SMatthew G. Knepley     ierr = DMLabelGetStratumIS(label, 1, &pointIS);CHKERRQ(ierr);
401652b88e8SMatthew G. Knepley     ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
402652b88e8SMatthew G. Knepley     for (field = 0, cellDof = 0, numComponents = 0; field < numFields; ++field) {
403652b88e8SMatthew G. Knepley       cellDof       += quadBd[field].numBasisFuncs*quadBd[field].numComponents;
404652b88e8SMatthew G. Knepley       numComponents += quadBd[field].numComponents;
405652b88e8SMatthew G. Knepley     }
406dcca6d9dSJed Brown     ierr = PetscMalloc7(numPoints*cellDof,&u,numPoints*dim,&v0,numPoints*dim,&n,numPoints*dim*dim,&J,numPoints*dim*dim,&invJ,numPoints,&detJ,numPoints*cellDof,&elemVec);CHKERRQ(ierr);
407652b88e8SMatthew G. Knepley     for (p = 0; p < numPoints; ++p) {
408652b88e8SMatthew G. Knepley       const PetscInt point = points[p];
409652b88e8SMatthew G. Knepley       PetscScalar   *x;
410652b88e8SMatthew G. Knepley       PetscInt       i;
411652b88e8SMatthew G. Knepley 
41202a80efeSMatthew G. Knepley       /* TODO: Add normal determination here */
413652b88e8SMatthew G. Knepley       ierr = DMPlexComputeCellGeometry(dm, point, &v0[p*dim], &J[p*dim*dim], &invJ[p*dim*dim], &detJ[p]);CHKERRQ(ierr);
4141d930511SMatthew G. Knepley       if (detJ[p] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for face %d", detJ[p], point);
415652b88e8SMatthew G. Knepley       ierr = DMPlexVecGetClosure(dm, NULL, X, point, NULL, &x);CHKERRQ(ierr);
416652b88e8SMatthew G. Knepley 
417652b88e8SMatthew G. Knepley       for (i = 0; i < cellDof; ++i) u[p*cellDof+i] = x[i];
418652b88e8SMatthew G. Knepley       ierr = DMPlexVecRestoreClosure(dm, NULL, X, point, NULL, &x);CHKERRQ(ierr);
419652b88e8SMatthew G. Knepley     }
420652b88e8SMatthew G. Knepley     for (field = 0; field < numFields; ++field) {
421652b88e8SMatthew G. Knepley       const PetscInt numQuadPoints = quadBd[field].numQuadPoints;
422652b88e8SMatthew G. Knepley       const PetscInt numBasisFuncs = quadBd[field].numBasisFuncs;
423a9dc2124SMatthew G. Knepley       void           (*f0)(const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->f0BdFuncs[field];
424a9dc2124SMatthew G. Knepley       void           (*f1)(const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->f1BdFuncs[field];
425652b88e8SMatthew G. Knepley       /* Conforming batches */
426652b88e8SMatthew G. Knepley       PetscInt blockSize  = numBasisFuncs*numQuadPoints;
427652b88e8SMatthew G. Knepley       PetscInt numBlocks  = 1;
428652b88e8SMatthew G. Knepley       PetscInt batchSize  = numBlocks * blockSize;
429652b88e8SMatthew G. Knepley       PetscInt numBatches = numBatchesTmp;
430652b88e8SMatthew G. Knepley       PetscInt numChunks  = numPoints / (numBatches*batchSize);
431652b88e8SMatthew G. Knepley       /* Remainder */
432652b88e8SMatthew G. Knepley       PetscInt numRemainder = numPoints % (numBatches * batchSize);
433652b88e8SMatthew G. Knepley       PetscInt offset       = numPoints - numRemainder;
434652b88e8SMatthew G. Knepley 
43502a80efeSMatthew G. Knepley       ierr = (*mesh->integrateBdResidualFEM)(numChunks*numBatches*batchSize, numFields, field, quadBd, u, v0, n, J, invJ, detJ, f0, f1, elemVec);CHKERRQ(ierr);
43602a80efeSMatthew G. Knepley       ierr = (*mesh->integrateBdResidualFEM)(numRemainder, numFields, field, quadBd, &u[offset*cellDof], &v0[offset*dim], &n[offset*dim], &J[offset*dim*dim], &invJ[offset*dim*dim], &detJ[offset],
437652b88e8SMatthew G. Knepley                                              f0, f1, &elemVec[offset*cellDof]);CHKERRQ(ierr);
438652b88e8SMatthew G. Knepley     }
439652b88e8SMatthew G. Knepley     for (p = 0; p < numPoints; ++p) {
440652b88e8SMatthew G. Knepley       const PetscInt point = points[p];
441652b88e8SMatthew G. Knepley 
442652b88e8SMatthew G. Knepley       if (mesh->printFEM > 1) {ierr = DMPrintCellVector(point, "Residual", cellDof, &elemVec[p*cellDof]);CHKERRQ(ierr);}
443652b88e8SMatthew G. Knepley       ierr = DMPlexVecSetClosure(dm, NULL, F, point, &elemVec[p*cellDof], ADD_VALUES);CHKERRQ(ierr);
444652b88e8SMatthew G. Knepley     }
445652b88e8SMatthew G. Knepley     ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
446652b88e8SMatthew G. Knepley     ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
44702a80efeSMatthew G. Knepley     ierr = PetscFree7(u,v0,n,J,invJ,detJ,elemVec);CHKERRQ(ierr);
448652b88e8SMatthew G. Knepley   }
449cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
450cb1e1211SMatthew G Knepley }
451cb1e1211SMatthew G Knepley 
452a0845e3aSMatthew G. Knepley #else
453a0845e3aSMatthew G. Knepley 
454a0845e3aSMatthew G. Knepley #undef __FUNCT__
455a0845e3aSMatthew G. Knepley #define __FUNCT__ "DMPlexComputeResidualFEM"
456a0845e3aSMatthew G. Knepley /*@
457a0845e3aSMatthew G. Knepley   DMPlexComputeResidualFEM - Form the local residual F from the local input X using pointwise functions specified by the user
458a0845e3aSMatthew G. Knepley 
459a0845e3aSMatthew G. Knepley   Input Parameters:
460a0845e3aSMatthew G. Knepley + dm - The mesh
461a0845e3aSMatthew G. Knepley . X  - Local input vector
462a0845e3aSMatthew G. Knepley - user - The user context
463a0845e3aSMatthew G. Knepley 
464a0845e3aSMatthew G. Knepley   Output Parameter:
465a0845e3aSMatthew G. Knepley . F  - Local output vector
466a0845e3aSMatthew G. Knepley 
467a0845e3aSMatthew G. Knepley   Note:
468a0845e3aSMatthew G. Knepley   The second member of the user context must be an FEMContext.
469a0845e3aSMatthew G. Knepley 
470a0845e3aSMatthew G. Knepley   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
471a0845e3aSMatthew G. Knepley   like a GPU, or vectorize on a multicore machine.
472a0845e3aSMatthew G. Knepley 
473a0845e3aSMatthew G. Knepley   Level: developer
474a0845e3aSMatthew G. Knepley 
475a0845e3aSMatthew G. Knepley .seealso: DMPlexComputeJacobianActionFEM()
476a0845e3aSMatthew G. Knepley @*/
477a0845e3aSMatthew G. Knepley PetscErrorCode DMPlexComputeResidualFEM(DM dm, Vec X, Vec F, void *user)
478a0845e3aSMatthew G. Knepley {
479a0845e3aSMatthew G. Knepley   DM_Plex          *mesh  = (DM_Plex *) dm->data;
4809a559087SMatthew G. Knepley   PetscFEM         *fem   = (PetscFEM *) user;
481a0845e3aSMatthew G. Knepley   PetscFE          *fe    = fem->fe;
4829a559087SMatthew G. Knepley   PetscFE          *feAux = fem->feAux;
483f1ea0e2fSMatthew G. Knepley   PetscFE          *feBd  = fem->feBd;
484a0845e3aSMatthew G. Knepley   const char       *name  = "Residual";
4859a559087SMatthew G. Knepley   DM                dmAux;
4869a559087SMatthew G. Knepley   Vec               A;
487a0845e3aSMatthew G. Knepley   PetscQuadrature   q;
488a0845e3aSMatthew G. Knepley   PetscCellGeometry geom;
4899a559087SMatthew G. Knepley   PetscSection      section, sectionAux;
490a0845e3aSMatthew G. Knepley   PetscReal        *v0, *J, *invJ, *detJ;
4919a559087SMatthew G. Knepley   PetscScalar      *elemVec, *u, *a;
4929a559087SMatthew G. Knepley   PetscInt          dim, Nf, NfAux = 0, f, numCells, cStart, cEnd, c;
493a0845e3aSMatthew G. Knepley   PetscInt          cellDof = 0, numComponents = 0;
4949a559087SMatthew G. Knepley   PetscInt          cellDofAux = 0, numComponentsAux = 0;
495a0845e3aSMatthew G. Knepley   PetscErrorCode    ierr;
496a0845e3aSMatthew G. Knepley 
497a0845e3aSMatthew G. Knepley   PetscFunctionBegin;
498a0845e3aSMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_ResidualFEM,dm,0,0,0);CHKERRQ(ierr);
499a0845e3aSMatthew G. Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
500a0845e3aSMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
5019a559087SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
502a0845e3aSMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
503a0845e3aSMatthew G. Knepley   numCells = cEnd - cStart;
5049a559087SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
505a0845e3aSMatthew G. Knepley     PetscInt Nb, Nc;
506a0845e3aSMatthew G. Knepley 
507a0845e3aSMatthew G. Knepley     ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr);
508a0845e3aSMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe[f], &Nc);CHKERRQ(ierr);
509a0845e3aSMatthew G. Knepley     cellDof       += Nb*Nc;
510a0845e3aSMatthew G. Knepley     numComponents += Nc;
511a0845e3aSMatthew G. Knepley   }
5129a559087SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
5139a559087SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr);
5149a559087SMatthew G. Knepley   if (dmAux) {
5159a559087SMatthew G. Knepley     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
5169a559087SMatthew G. Knepley     ierr = PetscSectionGetNumFields(sectionAux, &NfAux);CHKERRQ(ierr);
5179a559087SMatthew G. Knepley   }
5189a559087SMatthew G. Knepley   for (f = 0; f < NfAux; ++f) {
5199a559087SMatthew G. Knepley     PetscInt Nb, Nc;
5209a559087SMatthew G. Knepley 
5219a559087SMatthew G. Knepley     ierr = PetscFEGetDimension(feAux[f], &Nb);CHKERRQ(ierr);
5229a559087SMatthew G. Knepley     ierr = PetscFEGetNumComponents(feAux[f], &Nc);CHKERRQ(ierr);
5239a559087SMatthew G. Knepley     cellDofAux       += Nb*Nc;
5249a559087SMatthew G. Knepley     numComponentsAux += Nc;
5259a559087SMatthew G. Knepley   }
526*c110b1eeSGeoffrey Irving   ierr = DMPlexProjectFunctionLocal(dm, fe, fem->bcFuncs, fem->bcCtxs, INSERT_BC_VALUES, X);CHKERRQ(ierr);
527a0845e3aSMatthew G. Knepley   ierr = VecSet(F, 0.0);CHKERRQ(ierr);
528dcca6d9dSJed Brown   ierr = PetscMalloc6(numCells*cellDof,&u,numCells*dim,&v0,numCells*dim*dim,&J,numCells*dim*dim,&invJ,numCells,&detJ,numCells*cellDof,&elemVec);CHKERRQ(ierr);
529785e854fSJed Brown   if (dmAux) {ierr = PetscMalloc1(numCells*cellDofAux, &a);CHKERRQ(ierr);}
530a0845e3aSMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
531a0845e3aSMatthew G. Knepley     PetscScalar *x = NULL;
532a0845e3aSMatthew G. Knepley     PetscInt     i;
533a0845e3aSMatthew G. Knepley 
534a0845e3aSMatthew G. Knepley     ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr);
535a0845e3aSMatthew G. Knepley     if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c);
536a0845e3aSMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr);
537a0845e3aSMatthew G. Knepley     for (i = 0; i < cellDof; ++i) u[c*cellDof+i] = x[i];
538a0845e3aSMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr);
5399a559087SMatthew G. Knepley     if (dmAux) {
5409a559087SMatthew G. Knepley       ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
5419a559087SMatthew G. Knepley       for (i = 0; i < cellDofAux; ++i) a[c*cellDofAux+i] = x[i];
5429a559087SMatthew G. Knepley       ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
543a0845e3aSMatthew G. Knepley     }
5449a559087SMatthew G. Knepley   }
5459a559087SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
546c012ea0aSMatthew G. Knepley     void   (*f0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->f0Funcs[f];
547c012ea0aSMatthew G. Knepley     void   (*f1)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->f1Funcs[f];
548a0845e3aSMatthew G. Knepley     PetscInt Nb;
549a0845e3aSMatthew G. Knepley     /* Conforming batches */
550f30c5766SMatthew G. Knepley     PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
551a0845e3aSMatthew G. Knepley     /* Remainder */
552a0845e3aSMatthew G. Knepley     PetscInt Nr, offset;
553a0845e3aSMatthew G. Knepley 
554a0845e3aSMatthew G. Knepley     ierr = PetscFEGetQuadrature(fe[f], &q);CHKERRQ(ierr);
555a0845e3aSMatthew G. Knepley     ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr);
556f30c5766SMatthew G. Knepley     ierr = PetscFEGetTileSizes(fe[f], NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
557f9fd7fdbSMatthew G. Knepley     blockSize = Nb*q.numPoints;
558a0845e3aSMatthew G. Knepley     batchSize = numBlocks * blockSize;
559f30c5766SMatthew G. Knepley     ierr =  PetscFESetTileSizes(fe[f], blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
560a0845e3aSMatthew G. Knepley     numChunks = numCells / (numBatches*batchSize);
561a0845e3aSMatthew G. Knepley     Ne        = numChunks*numBatches*batchSize;
562a0845e3aSMatthew G. Knepley     Nr        = numCells % (numBatches*batchSize);
563a0845e3aSMatthew G. Knepley     offset    = numCells - Nr;
564a0845e3aSMatthew G. Knepley     geom.v0   = v0;
565a0845e3aSMatthew G. Knepley     geom.J    = J;
566a0845e3aSMatthew G. Knepley     geom.invJ = invJ;
567a0845e3aSMatthew G. Knepley     geom.detJ = detJ;
5689a559087SMatthew G. Knepley     ierr = PetscFEIntegrateResidual(fe[f], Ne, Nf, fe, f, geom, u, NfAux, feAux, a, f0, f1, elemVec);CHKERRQ(ierr);
569a0845e3aSMatthew G. Knepley     geom.v0   = &v0[offset*dim];
570a0845e3aSMatthew G. Knepley     geom.J    = &J[offset*dim*dim];
571a0845e3aSMatthew G. Knepley     geom.invJ = &invJ[offset*dim*dim];
572a0845e3aSMatthew G. Knepley     geom.detJ = &detJ[offset];
5739a559087SMatthew 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);
574a0845e3aSMatthew G. Knepley   }
575a0845e3aSMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
576a0845e3aSMatthew G. Knepley     if (mesh->printFEM > 1) {ierr = DMPrintCellVector(c, name, cellDof, &elemVec[c*cellDof]);CHKERRQ(ierr);}
577a0845e3aSMatthew G. Knepley     ierr = DMPlexVecSetClosure(dm, section, F, c, &elemVec[c*cellDof], ADD_VALUES);CHKERRQ(ierr);
578a0845e3aSMatthew G. Knepley   }
579a0845e3aSMatthew G. Knepley   ierr = PetscFree6(u,v0,J,invJ,detJ,elemVec);CHKERRQ(ierr);
5809a559087SMatthew G. Knepley   if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);}
581f1ea0e2fSMatthew G. Knepley   if (feBd) {
582075da914SMatthew G. Knepley     DMLabel         label, depth;
583f1ea0e2fSMatthew G. Knepley     IS              pointIS;
584f1ea0e2fSMatthew G. Knepley     const PetscInt *points;
585075da914SMatthew G. Knepley     PetscInt        dep, numPoints, p, numFaces;
586f1ea0e2fSMatthew G. Knepley     PetscReal      *n;
587f1ea0e2fSMatthew G. Knepley 
588f1ea0e2fSMatthew G. Knepley     ierr = DMPlexGetLabel(dm, "boundary", &label);CHKERRQ(ierr);
589075da914SMatthew G. Knepley     ierr = DMPlexGetDepthLabel(dm, &depth);CHKERRQ(ierr);
590f1ea0e2fSMatthew G. Knepley     ierr = DMLabelGetStratumSize(label, 1, &numPoints);CHKERRQ(ierr);
591f1ea0e2fSMatthew G. Knepley     ierr = DMLabelGetStratumIS(label, 1, &pointIS);CHKERRQ(ierr);
592f1ea0e2fSMatthew G. Knepley     ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
593f1ea0e2fSMatthew G. Knepley     for (f = 0, cellDof = 0, numComponents = 0; f < Nf; ++f) {
594f1ea0e2fSMatthew G. Knepley       PetscInt Nb, Nc;
595f1ea0e2fSMatthew G. Knepley 
596f1ea0e2fSMatthew G. Knepley       ierr = PetscFEGetDimension(feBd[f], &Nb);CHKERRQ(ierr);
597f1ea0e2fSMatthew G. Knepley       ierr = PetscFEGetNumComponents(feBd[f], &Nc);CHKERRQ(ierr);
598f1ea0e2fSMatthew G. Knepley       cellDof       += Nb*Nc;
599f1ea0e2fSMatthew G. Knepley       numComponents += Nc;
600f1ea0e2fSMatthew G. Knepley     }
601075da914SMatthew G. Knepley     for (p = 0, numFaces = 0; p < numPoints; ++p) {
602075da914SMatthew G. Knepley       ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr);
603075da914SMatthew G. Knepley       if (dep == dim-1) ++numFaces;
604075da914SMatthew G. Knepley     }
605dcca6d9dSJed 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);
606075da914SMatthew G. Knepley     for (p = 0, f = 0; p < numPoints; ++p) {
607f1ea0e2fSMatthew G. Knepley       const PetscInt point = points[p];
608f1ea0e2fSMatthew G. Knepley       PetscScalar   *x     = NULL;
609f1ea0e2fSMatthew G. Knepley       PetscInt       i;
610f1ea0e2fSMatthew G. Knepley 
611075da914SMatthew G. Knepley       ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr);
612075da914SMatthew G. Knepley       if (dep != dim-1) continue;
613075da914SMatthew G. Knepley       ierr = DMPlexComputeCellGeometry(dm, point, &v0[f*dim], &J[f*dim*dim], &invJ[f*dim*dim], &detJ[f]);CHKERRQ(ierr);
614a8007bbfSMatthew G. Knepley       ierr = DMPlexComputeCellGeometryFVM(dm, point, NULL, NULL, &n[f*dim]);
615075da914SMatthew G. Knepley       if (detJ[f] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for face %d", detJ[f], point);
616f1ea0e2fSMatthew G. Knepley       ierr = DMPlexVecGetClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr);
617075da914SMatthew G. Knepley       for (i = 0; i < cellDof; ++i) u[f*cellDof+i] = x[i];
618f1ea0e2fSMatthew G. Knepley       ierr = DMPlexVecRestoreClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr);
619075da914SMatthew G. Knepley       ++f;
620f1ea0e2fSMatthew G. Knepley     }
621f1ea0e2fSMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
622f1ea0e2fSMatthew G. Knepley       void   (*f0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->f0BdFuncs[f];
623f1ea0e2fSMatthew G. Knepley       void   (*f1)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->f1BdFuncs[f];
624f1ea0e2fSMatthew G. Knepley       PetscInt Nb;
625f1ea0e2fSMatthew G. Knepley       /* Conforming batches */
626f1ea0e2fSMatthew G. Knepley       PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
627f1ea0e2fSMatthew G. Knepley       /* Remainder */
628f1ea0e2fSMatthew G. Knepley       PetscInt Nr, offset;
629f1ea0e2fSMatthew G. Knepley 
630f1ea0e2fSMatthew G. Knepley       ierr = PetscFEGetQuadrature(feBd[f], &q);CHKERRQ(ierr);
631f1ea0e2fSMatthew G. Knepley       ierr = PetscFEGetDimension(feBd[f], &Nb);CHKERRQ(ierr);
632f1ea0e2fSMatthew G. Knepley       ierr = PetscFEGetTileSizes(feBd[f], NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
633f9fd7fdbSMatthew G. Knepley       blockSize = Nb*q.numPoints;
634f1ea0e2fSMatthew G. Knepley       batchSize = numBlocks * blockSize;
635f1ea0e2fSMatthew G. Knepley       ierr =  PetscFESetTileSizes(feBd[f], blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
636075da914SMatthew G. Knepley       numChunks = numFaces / (numBatches*batchSize);
637f1ea0e2fSMatthew G. Knepley       Ne        = numChunks*numBatches*batchSize;
638075da914SMatthew G. Knepley       Nr        = numFaces % (numBatches*batchSize);
639075da914SMatthew G. Knepley       offset    = numFaces - Nr;
640f1ea0e2fSMatthew G. Knepley       geom.v0   = v0;
641f1ea0e2fSMatthew G. Knepley       geom.n    = n;
642f1ea0e2fSMatthew G. Knepley       geom.J    = J;
643f1ea0e2fSMatthew G. Knepley       geom.invJ = invJ;
644f1ea0e2fSMatthew G. Knepley       geom.detJ = detJ;
645f1ea0e2fSMatthew G. Knepley       ierr = PetscFEIntegrateBdResidual(feBd[f], Ne, Nf, feBd, f, geom, u, 0, NULL, NULL, f0, f1, elemVec);CHKERRQ(ierr);
646f1ea0e2fSMatthew G. Knepley       geom.v0   = &v0[offset*dim];
647f1ea0e2fSMatthew G. Knepley       geom.n    = &n[offset*dim];
648f1ea0e2fSMatthew G. Knepley       geom.J    = &J[offset*dim*dim];
649f1ea0e2fSMatthew G. Knepley       geom.invJ = &invJ[offset*dim*dim];
650f1ea0e2fSMatthew G. Knepley       geom.detJ = &detJ[offset];
651f1ea0e2fSMatthew G. Knepley       ierr = PetscFEIntegrateBdResidual(feBd[f], Nr, Nf, feBd, f, geom, &u[offset*cellDof], 0, NULL, NULL, f0, f1, &elemVec[offset*cellDof]);CHKERRQ(ierr);
652f1ea0e2fSMatthew G. Knepley     }
653075da914SMatthew G. Knepley     for (p = 0, f = 0; p < numPoints; ++p) {
654f1ea0e2fSMatthew G. Knepley       const PetscInt point = points[p];
655f1ea0e2fSMatthew G. Knepley 
656075da914SMatthew G. Knepley       ierr = DMLabelGetValue(depth, point, &dep);CHKERRQ(ierr);
657075da914SMatthew G. Knepley       if (dep != dim-1) continue;
658075da914SMatthew G. Knepley       if (mesh->printFEM > 1) {ierr = DMPrintCellVector(point, "BdResidual", cellDof, &elemVec[f*cellDof]);CHKERRQ(ierr);}
659075da914SMatthew G. Knepley       ierr = DMPlexVecSetClosure(dm, NULL, F, point, &elemVec[f*cellDof], ADD_VALUES);CHKERRQ(ierr);
660075da914SMatthew G. Knepley       ++f;
661f1ea0e2fSMatthew G. Knepley     }
662f1ea0e2fSMatthew G. Knepley     ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
663f1ea0e2fSMatthew G. Knepley     ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
664f1ea0e2fSMatthew G. Knepley     ierr = PetscFree7(u,v0,n,J,invJ,detJ,elemVec);CHKERRQ(ierr);
665f1ea0e2fSMatthew G. Knepley   }
6666113b454SMatthew G. Knepley   if (mesh->printFEM) {ierr = DMPrintLocalVec(dm, name, mesh->printTol, F);CHKERRQ(ierr);}
667a0845e3aSMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_ResidualFEM,dm,0,0,0);CHKERRQ(ierr);
668a0845e3aSMatthew G. Knepley   PetscFunctionReturn(0);
669a0845e3aSMatthew G. Knepley }
670a0845e3aSMatthew G. Knepley 
671a0845e3aSMatthew G. Knepley #endif
672a0845e3aSMatthew G. Knepley 
673cb1e1211SMatthew G Knepley #undef __FUNCT__
674cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexComputeJacobianActionFEM"
675cb1e1211SMatthew G Knepley /*@C
676cb1e1211SMatthew G Knepley   DMPlexComputeJacobianActionFEM - Form the local action of Jacobian J(u) on the local input X using pointwise functions specified by the user
677cb1e1211SMatthew G Knepley 
678cb1e1211SMatthew G Knepley   Input Parameters:
679cb1e1211SMatthew G Knepley + dm - The mesh
680cb1e1211SMatthew G Knepley . J  - The Jacobian shell matrix
681cb1e1211SMatthew G Knepley . X  - Local input vector
682cb1e1211SMatthew G Knepley - user - The user context
683cb1e1211SMatthew G Knepley 
684cb1e1211SMatthew G Knepley   Output Parameter:
685cb1e1211SMatthew G Knepley . F  - Local output vector
686cb1e1211SMatthew G Knepley 
687cb1e1211SMatthew G Knepley   Note:
688cb1e1211SMatthew G Knepley   The second member of the user context must be an FEMContext.
689cb1e1211SMatthew G Knepley 
690cb1e1211SMatthew G Knepley   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
691cb1e1211SMatthew G Knepley   like a GPU, or vectorize on a multicore machine.
692cb1e1211SMatthew G Knepley 
6930059ad2aSSatish Balay   Level: developer
6940059ad2aSSatish Balay 
695cb1e1211SMatthew G Knepley .seealso: DMPlexComputeResidualFEM()
696878cb397SSatish Balay @*/
697cb1e1211SMatthew G Knepley PetscErrorCode DMPlexComputeJacobianActionFEM(DM dm, Mat Jac, Vec X, Vec F, void *user)
698cb1e1211SMatthew G Knepley {
699cb1e1211SMatthew G Knepley   DM_Plex          *mesh = (DM_Plex *) dm->data;
7009a559087SMatthew G. Knepley   PetscFEM         *fem  = (PetscFEM *) user;
7010483ade4SMatthew G. Knepley   PetscFE          *fe   = fem->fe;
7020483ade4SMatthew G. Knepley   PetscQuadrature   quad;
7030483ade4SMatthew G. Knepley   PetscCellGeometry geom;
704cb1e1211SMatthew G Knepley   PetscSection      section;
705cb1e1211SMatthew G Knepley   JacActionCtx     *jctx;
706cb1e1211SMatthew G Knepley   PetscReal        *v0, *J, *invJ, *detJ;
707cb1e1211SMatthew G Knepley   PetscScalar      *elemVec, *u, *a;
7080483ade4SMatthew G. Knepley   PetscInt          dim, numFields, field, numCells, cStart, cEnd, c;
709cb1e1211SMatthew G Knepley   PetscInt          cellDof = 0;
710cb1e1211SMatthew G Knepley   PetscErrorCode    ierr;
711cb1e1211SMatthew G Knepley 
712cb1e1211SMatthew G Knepley   PetscFunctionBegin;
7130483ade4SMatthew G. Knepley   /* ierr = PetscLogEventBegin(DMPLEX_JacobianActionFEM,dm,0,0,0);CHKERRQ(ierr); */
714cb1e1211SMatthew G Knepley   ierr = MatShellGetContext(Jac, &jctx);CHKERRQ(ierr);
715cb1e1211SMatthew G Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
716cb1e1211SMatthew G Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
717cb1e1211SMatthew G Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
718cb1e1211SMatthew G Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
719cb1e1211SMatthew G Knepley   numCells = cEnd - cStart;
720cb1e1211SMatthew G Knepley   for (field = 0; field < numFields; ++field) {
7210483ade4SMatthew G. Knepley     PetscInt Nb, Nc;
7220483ade4SMatthew G. Knepley 
7230483ade4SMatthew G. Knepley     ierr = PetscFEGetDimension(fe[field], &Nb);CHKERRQ(ierr);
7240483ade4SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe[field], &Nc);CHKERRQ(ierr);
7250483ade4SMatthew G. Knepley     cellDof += Nb*Nc;
726cb1e1211SMatthew G Knepley   }
727cb1e1211SMatthew G Knepley   ierr = VecSet(F, 0.0);CHKERRQ(ierr);
728dcca6d9dSJed 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);
729cb1e1211SMatthew G Knepley   for (c = cStart; c < cEnd; ++c) {
730a1e44745SMatthew G. Knepley     PetscScalar *x = NULL;
731cb1e1211SMatthew G Knepley     PetscInt     i;
732cb1e1211SMatthew G Knepley 
733cb1e1211SMatthew G Knepley     ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr);
734cb1e1211SMatthew G Knepley     if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c);
735cb1e1211SMatthew G Knepley     ierr = DMPlexVecGetClosure(dm, NULL, jctx->u, c, NULL, &x);CHKERRQ(ierr);
736cb1e1211SMatthew G Knepley     for (i = 0; i < cellDof; ++i) u[c*cellDof+i] = x[i];
737cb1e1211SMatthew G Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, jctx->u, c, NULL, &x);CHKERRQ(ierr);
738cb1e1211SMatthew G Knepley     ierr = DMPlexVecGetClosure(dm, NULL, X, c, NULL, &x);CHKERRQ(ierr);
739cb1e1211SMatthew G Knepley     for (i = 0; i < cellDof; ++i) a[c*cellDof+i] = x[i];
740cb1e1211SMatthew G Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, X, c, NULL, &x);CHKERRQ(ierr);
741cb1e1211SMatthew G Knepley   }
742cb1e1211SMatthew G Knepley   for (field = 0; field < numFields; ++field) {
7430483ade4SMatthew G. Knepley     PetscInt Nb;
744cb1e1211SMatthew G Knepley     /* Conforming batches */
745cb1e1211SMatthew G Knepley     PetscInt numBlocks  = 1;
7460483ade4SMatthew G. Knepley     PetscInt numBatches = 1;
7470483ade4SMatthew G. Knepley     PetscInt numChunks, Ne, blockSize, batchSize;
748cb1e1211SMatthew G Knepley     /* Remainder */
7490483ade4SMatthew G. Knepley     PetscInt Nr, offset;
750cb1e1211SMatthew G Knepley 
7510483ade4SMatthew G. Knepley     ierr = PetscFEGetQuadrature(fe[field], &quad);CHKERRQ(ierr);
7520483ade4SMatthew G. Knepley     ierr = PetscFEGetDimension(fe[field], &Nb);CHKERRQ(ierr);
753f9fd7fdbSMatthew G. Knepley     blockSize = Nb*quad.numPoints;
7540483ade4SMatthew G. Knepley     batchSize = numBlocks * blockSize;
7550483ade4SMatthew G. Knepley     numChunks = numCells / (numBatches*batchSize);
7560483ade4SMatthew G. Knepley     Ne        = numChunks*numBatches*batchSize;
7570483ade4SMatthew G. Knepley     Nr        = numCells % (numBatches*batchSize);
7580483ade4SMatthew G. Knepley     offset    = numCells - Nr;
7590483ade4SMatthew G. Knepley     geom.v0   = v0;
7600483ade4SMatthew G. Knepley     geom.J    = J;
7610483ade4SMatthew G. Knepley     geom.invJ = invJ;
7620483ade4SMatthew G. Knepley     geom.detJ = detJ;
7630483ade4SMatthew G. Knepley     ierr = PetscFEIntegrateJacobianAction(fe[field], Ne, numFields, fe, field, geom, u, a, fem->g0Funcs, fem->g1Funcs, fem->g2Funcs, fem->g3Funcs, elemVec);CHKERRQ(ierr);
7640483ade4SMatthew G. Knepley     geom.v0   = &v0[offset*dim];
7650483ade4SMatthew G. Knepley     geom.J    = &J[offset*dim*dim];
7660483ade4SMatthew G. Knepley     geom.invJ = &invJ[offset*dim*dim];
7670483ade4SMatthew G. Knepley     geom.detJ = &detJ[offset];
7680483ade4SMatthew G. Knepley     ierr = PetscFEIntegrateJacobianAction(fe[field], Nr, numFields, fe, field, geom, &u[offset*cellDof], &a[offset*cellDof],
769cb1e1211SMatthew G Knepley                                           fem->g0Funcs, fem->g1Funcs, fem->g2Funcs, fem->g3Funcs, &elemVec[offset*cellDof]);CHKERRQ(ierr);
770cb1e1211SMatthew G Knepley   }
771cb1e1211SMatthew G Knepley   for (c = cStart; c < cEnd; ++c) {
772cb1e1211SMatthew G Knepley     if (mesh->printFEM > 1) {ierr = DMPrintCellVector(c, "Jacobian Action", cellDof, &elemVec[c*cellDof]);CHKERRQ(ierr);}
773cb1e1211SMatthew G Knepley     ierr = DMPlexVecSetClosure(dm, NULL, F, c, &elemVec[c*cellDof], ADD_VALUES);CHKERRQ(ierr);
774cb1e1211SMatthew G Knepley   }
775cb1e1211SMatthew G Knepley   ierr = PetscFree7(u,a,v0,J,invJ,detJ,elemVec);CHKERRQ(ierr);
776cb1e1211SMatthew G Knepley   if (mesh->printFEM) {
777cb1e1211SMatthew G Knepley     PetscMPIInt rank, numProcs;
778cb1e1211SMatthew G Knepley     PetscInt    p;
779cb1e1211SMatthew G Knepley 
780cb1e1211SMatthew G Knepley     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr);
781cb1e1211SMatthew G Knepley     ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm), &numProcs);CHKERRQ(ierr);
78286a74ee0SMatthew G. Knepley     ierr = PetscPrintf(PetscObjectComm((PetscObject)dm), "Jacobian Action:\n");CHKERRQ(ierr);
783cb1e1211SMatthew G Knepley     for (p = 0; p < numProcs; ++p) {
784cb1e1211SMatthew G Knepley       if (p == rank) {ierr = VecView(F, PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);}
785cb1e1211SMatthew G Knepley       ierr = PetscBarrier((PetscObject) dm);CHKERRQ(ierr);
786cb1e1211SMatthew G Knepley     }
787cb1e1211SMatthew G Knepley   }
7880483ade4SMatthew G. Knepley   /* ierr = PetscLogEventEnd(DMPLEX_JacobianActionFEM,dm,0,0,0);CHKERRQ(ierr); */
789cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
790cb1e1211SMatthew G Knepley }
791cb1e1211SMatthew G Knepley 
792cb1e1211SMatthew G Knepley #undef __FUNCT__
793cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexComputeJacobianFEM"
794cb1e1211SMatthew G Knepley /*@
795cb1e1211SMatthew G Knepley   DMPlexComputeJacobianFEM - Form the local portion of the Jacobian matrix J at the local solution X using pointwise functions specified by the user.
796cb1e1211SMatthew G Knepley 
797cb1e1211SMatthew G Knepley   Input Parameters:
798cb1e1211SMatthew G Knepley + dm - The mesh
799cb1e1211SMatthew G Knepley . X  - Local input vector
800cb1e1211SMatthew G Knepley - user - The user context
801cb1e1211SMatthew G Knepley 
802cb1e1211SMatthew G Knepley   Output Parameter:
803cb1e1211SMatthew G Knepley . Jac  - Jacobian matrix
804cb1e1211SMatthew G Knepley 
805cb1e1211SMatthew G Knepley   Note:
806cb1e1211SMatthew G Knepley   The second member of the user context must be an FEMContext.
807cb1e1211SMatthew G Knepley 
808cb1e1211SMatthew G Knepley   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
809cb1e1211SMatthew G Knepley   like a GPU, or vectorize on a multicore machine.
810cb1e1211SMatthew G Knepley 
8110059ad2aSSatish Balay   Level: developer
8120059ad2aSSatish Balay 
813cb1e1211SMatthew G Knepley .seealso: FormFunctionLocal()
814878cb397SSatish Balay @*/
815cb1e1211SMatthew G Knepley PetscErrorCode DMPlexComputeJacobianFEM(DM dm, Vec X, Mat Jac, Mat JacP, MatStructure *str,void *user)
816cb1e1211SMatthew G Knepley {
817cb1e1211SMatthew G Knepley   DM_Plex          *mesh  = (DM_Plex *) dm->data;
8189a559087SMatthew G. Knepley   PetscFEM         *fem   = (PetscFEM *) user;
819a319912fSMatthew G. Knepley   PetscFE          *fe    = fem->fe;
820754551f4SMatthew G. Knepley   PetscFE          *feAux = fem->feAux;
821a319912fSMatthew G. Knepley   const char       *name  = "Jacobian";
822754551f4SMatthew G. Knepley   DM                dmAux;
823754551f4SMatthew G. Knepley   Vec               A;
824a319912fSMatthew G. Knepley   PetscQuadrature   quad;
825a319912fSMatthew G. Knepley   PetscCellGeometry geom;
826754551f4SMatthew G. Knepley   PetscSection      section, globalSection, sectionAux;
827cb1e1211SMatthew G Knepley   PetscReal        *v0, *J, *invJ, *detJ;
828754551f4SMatthew G. Knepley   PetscScalar      *elemMat, *u, *a;
829754551f4SMatthew G. Knepley   PetscInt          dim, Nf, NfAux = 0, f, fieldI, fieldJ, numCells, cStart, cEnd, c;
830cb1e1211SMatthew G Knepley   PetscInt          cellDof = 0, numComponents = 0;
831754551f4SMatthew G. Knepley   PetscInt          cellDofAux = 0, numComponentsAux = 0;
832cb1e1211SMatthew G Knepley   PetscBool         isShell;
833cb1e1211SMatthew G Knepley   PetscErrorCode    ierr;
834cb1e1211SMatthew G Knepley 
835cb1e1211SMatthew G Knepley   PetscFunctionBegin;
836a319912fSMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_JacobianFEM,dm,0,0,0);CHKERRQ(ierr);
837cb1e1211SMatthew G Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
838cb1e1211SMatthew G Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
839a319912fSMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dm, &globalSection);CHKERRQ(ierr);
840754551f4SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
841cb1e1211SMatthew G Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
842cb1e1211SMatthew G Knepley   numCells = cEnd - cStart;
843754551f4SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
844a319912fSMatthew G. Knepley     PetscInt Nb, Nc;
845a319912fSMatthew G. Knepley 
846a319912fSMatthew G. Knepley     ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr);
847a319912fSMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe[f], &Nc);CHKERRQ(ierr);
848a319912fSMatthew G. Knepley     cellDof       += Nb*Nc;
849a319912fSMatthew G. Knepley     numComponents += Nc;
850cb1e1211SMatthew G Knepley   }
851754551f4SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
852754551f4SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr);
853754551f4SMatthew G. Knepley   if (dmAux) {
854754551f4SMatthew G. Knepley     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
855754551f4SMatthew G. Knepley     ierr = PetscSectionGetNumFields(sectionAux, &NfAux);CHKERRQ(ierr);
856754551f4SMatthew G. Knepley   }
857754551f4SMatthew G. Knepley   for (f = 0; f < NfAux; ++f) {
858754551f4SMatthew G. Knepley     PetscInt Nb, Nc;
859754551f4SMatthew G. Knepley 
860754551f4SMatthew G. Knepley     ierr = PetscFEGetDimension(feAux[f], &Nb);CHKERRQ(ierr);
861754551f4SMatthew G. Knepley     ierr = PetscFEGetNumComponents(feAux[f], &Nc);CHKERRQ(ierr);
862754551f4SMatthew G. Knepley     cellDofAux       += Nb*Nc;
863754551f4SMatthew G. Knepley     numComponentsAux += Nc;
864754551f4SMatthew G. Knepley   }
865*c110b1eeSGeoffrey Irving   ierr = DMPlexProjectFunctionLocal(dm, fe, fem->bcFuncs, fem->bcCtxs, INSERT_BC_VALUES, X);CHKERRQ(ierr);
866cb1e1211SMatthew G Knepley   ierr = MatZeroEntries(JacP);CHKERRQ(ierr);
867dcca6d9dSJed 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);
868785e854fSJed Brown   if (dmAux) {ierr = PetscMalloc1(numCells*cellDofAux, &a);CHKERRQ(ierr);}
869cb1e1211SMatthew G Knepley   for (c = cStart; c < cEnd; ++c) {
870a1e44745SMatthew G. Knepley     PetscScalar *x = NULL;
871cb1e1211SMatthew G Knepley     PetscInt     i;
872cb1e1211SMatthew G Knepley 
873cb1e1211SMatthew G Knepley     ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr);
874cb1e1211SMatthew G Knepley     if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c);
875a319912fSMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr);
876cb1e1211SMatthew G Knepley     for (i = 0; i < cellDof; ++i) u[c*cellDof+i] = x[i];
877a319912fSMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr);
878754551f4SMatthew G. Knepley     if (dmAux) {
879754551f4SMatthew G. Knepley       ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
880754551f4SMatthew G. Knepley       for (i = 0; i < cellDofAux; ++i) a[c*cellDofAux+i] = x[i];
881754551f4SMatthew G. Knepley       ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
882754551f4SMatthew G. Knepley     }
883cb1e1211SMatthew G Knepley   }
884cb1e1211SMatthew G Knepley   ierr = PetscMemzero(elemMat, numCells*cellDof*cellDof * sizeof(PetscScalar));CHKERRQ(ierr);
885754551f4SMatthew G. Knepley   for (fieldI = 0; fieldI < Nf; ++fieldI) {
886a319912fSMatthew G. Knepley     PetscInt Nb;
887cb1e1211SMatthew G Knepley     /* Conforming batches */
888754551f4SMatthew G. Knepley     PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
889cb1e1211SMatthew G Knepley     /* Remainder */
890a319912fSMatthew G. Knepley     PetscInt Nr, offset;
891cb1e1211SMatthew G Knepley 
892754551f4SMatthew G. Knepley     ierr = PetscFEGetQuadrature(fe[fieldI], &quad);CHKERRQ(ierr);
893754551f4SMatthew G. Knepley     ierr = PetscFEGetDimension(fe[fieldI], &Nb);CHKERRQ(ierr);
894754551f4SMatthew G. Knepley     ierr = PetscFEGetTileSizes(fe[fieldI], NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
895f9fd7fdbSMatthew G. Knepley     blockSize = Nb*quad.numPoints;
896a319912fSMatthew G. Knepley     batchSize = numBlocks * blockSize;
897754551f4SMatthew G. Knepley     ierr = PetscFESetTileSizes(fe[fieldI], blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
898a319912fSMatthew G. Knepley     numChunks = numCells / (numBatches*batchSize);
899a319912fSMatthew G. Knepley     Ne        = numChunks*numBatches*batchSize;
900a319912fSMatthew G. Knepley     Nr        = numCells % (numBatches*batchSize);
901a319912fSMatthew G. Knepley     offset    = numCells - Nr;
902754551f4SMatthew G. Knepley     for (fieldJ = 0; fieldJ < Nf; ++fieldJ) {
903754551f4SMatthew G. Knepley       void   (*g0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g0Funcs[fieldI*Nf+fieldJ];
904754551f4SMatthew G. Knepley       void   (*g1)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g1Funcs[fieldI*Nf+fieldJ];
905754551f4SMatthew G. Knepley       void   (*g2)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g2Funcs[fieldI*Nf+fieldJ];
906754551f4SMatthew G. Knepley       void   (*g3)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g3Funcs[fieldI*Nf+fieldJ];
907754551f4SMatthew G. Knepley 
908a319912fSMatthew G. Knepley       geom.v0   = v0;
909a319912fSMatthew G. Knepley       geom.J    = J;
910a319912fSMatthew G. Knepley       geom.invJ = invJ;
911a319912fSMatthew G. Knepley       geom.detJ = detJ;
912754551f4SMatthew G. Knepley       ierr = PetscFEIntegrateJacobian(fe[fieldI], Ne, Nf, fe, fieldI, fieldJ, geom, u, NfAux, feAux, a, g0, g1, g2, g3, elemMat);CHKERRQ(ierr);
913a319912fSMatthew G. Knepley       geom.v0   = &v0[offset*dim];
914a319912fSMatthew G. Knepley       geom.J    = &J[offset*dim*dim];
915a319912fSMatthew G. Knepley       geom.invJ = &invJ[offset*dim*dim];
916a319912fSMatthew G. Knepley       geom.detJ = &detJ[offset];
917754551f4SMatthew 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);
918cb1e1211SMatthew G Knepley     }
919cb1e1211SMatthew G Knepley   }
920cb1e1211SMatthew G Knepley   for (c = cStart; c < cEnd; ++c) {
921a319912fSMatthew G. Knepley     if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(c, name, cellDof, cellDof, &elemMat[c*cellDof*cellDof]);CHKERRQ(ierr);}
922a319912fSMatthew G. Knepley     ierr = DMPlexMatSetClosure(dm, section, globalSection, JacP, c, &elemMat[c*cellDof*cellDof], ADD_VALUES);CHKERRQ(ierr);
923cb1e1211SMatthew G Knepley   }
924cb1e1211SMatthew G Knepley   ierr = PetscFree6(u,v0,J,invJ,detJ,elemMat);CHKERRQ(ierr);
925754551f4SMatthew G. Knepley   if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);}
926cb1e1211SMatthew G Knepley   ierr = MatAssemblyBegin(JacP, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
927cb1e1211SMatthew G Knepley   ierr = MatAssemblyEnd(JacP, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
928cb1e1211SMatthew G Knepley   if (mesh->printFEM) {
929a319912fSMatthew G. Knepley     ierr = PetscPrintf(PETSC_COMM_WORLD, "%s:\n", name);CHKERRQ(ierr);
930cb1e1211SMatthew G Knepley     ierr = MatChop(JacP, 1.0e-10);CHKERRQ(ierr);
931cb1e1211SMatthew G Knepley     ierr = MatView(JacP, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
932cb1e1211SMatthew G Knepley   }
933a319912fSMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_JacobianFEM,dm,0,0,0);CHKERRQ(ierr);
934cb1e1211SMatthew G Knepley   ierr = PetscObjectTypeCompare((PetscObject) Jac, MATSHELL, &isShell);CHKERRQ(ierr);
935cb1e1211SMatthew G Knepley   if (isShell) {
936cb1e1211SMatthew G Knepley     JacActionCtx *jctx;
937cb1e1211SMatthew G Knepley 
938cb1e1211SMatthew G Knepley     ierr = MatShellGetContext(Jac, &jctx);CHKERRQ(ierr);
939cb1e1211SMatthew G Knepley     ierr = VecCopy(X, jctx->u);CHKERRQ(ierr);
940cb1e1211SMatthew G Knepley   }
941cb1e1211SMatthew G Knepley   *str = SAME_NONZERO_PATTERN;
942cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
943cb1e1211SMatthew G Knepley }
944