xref: /petsc/src/dm/impls/plex/plexfem.c (revision 72f94c41d68956ae88eca928477c9616bca58bb7)
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);
126cb1e1211SMatthew G Knepley   ierr = DMPlexGetCoordinateSection(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*72f94c41SMatthew G. Knepley PetscErrorCode DMPlexProjectFunctionLocal(DM dm, PetscFE fe[], void (**funcs)(const PetscReal [], PetscScalar *), InsertMode mode, Vec localX)
186cb1e1211SMatthew G Knepley {
187*72f94c41SMatthew G. Knepley   PetscDualSpace *sp;
188*72f94c41SMatthew G. Knepley   PetscSection    section;
189*72f94c41SMatthew G. Knepley   PetscScalar    *values;
190*72f94c41SMatthew G. Knepley   PetscReal      *v0, *J, detJ;
191*72f94c41SMatthew 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);
196*72f94c41SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
197*72f94c41SMatthew G. Knepley   ierr = PetscMalloc(numFields * sizeof(PetscDualSpace), &sp);CHKERRQ(ierr);
198*72f94c41SMatthew G. Knepley   for (f = 0; f < numFields; ++f) {
199*72f94c41SMatthew G. Knepley     ierr = PetscFEGetDualSpace(fe[f], &sp[f]);CHKERRQ(ierr);
200*72f94c41SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe[f], &numComp);CHKERRQ(ierr);
201*72f94c41SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
202*72f94c41SMatthew G. Knepley     totDim += spDim*numComp;
203cb1e1211SMatthew G Knepley   }
204*72f94c41SMatthew G. Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
205*72f94c41SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
206*72f94c41SMatthew G. Knepley   ierr = DMPlexVecGetClosure(dm, section, localX, cStart, &numValues, NULL);CHKERRQ(ierr);
207*72f94c41SMatthew 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);
208*72f94c41SMatthew G. Knepley   ierr = DMGetWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
209cb1e1211SMatthew G Knepley   ierr = PetscMalloc2(dim,PetscReal,&v0,dim*dim,PetscReal,&J);CHKERRQ(ierr);
210*72f94c41SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
211*72f94c41SMatthew G. Knepley     PetscCellGeometry geom;
212cb1e1211SMatthew G Knepley 
213cb1e1211SMatthew G Knepley     ierr = DMPlexComputeCellGeometry(dm, c, v0, J, NULL, &detJ);CHKERRQ(ierr);
214*72f94c41SMatthew G. Knepley     geom.v0   = v0;
215*72f94c41SMatthew G. Knepley     geom.J    = J;
216*72f94c41SMatthew G. Knepley     geom.detJ = &detJ;
217*72f94c41SMatthew G. Knepley     for (f = 0, v = 0; f < numFields; ++f) {
218*72f94c41SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe[f], &numComp);CHKERRQ(ierr);
219*72f94c41SMatthew G. Knepley       ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
220*72f94c41SMatthew G. Knepley       for (d = 0; d < spDim; ++d) {
221*72f94c41SMatthew G. Knepley         ierr = PetscDualSpaceApply(sp[f], d, geom, numComp, funcs[f], &values[v]);CHKERRQ(ierr);
222*72f94c41SMatthew G. Knepley         v += numComp;
223cb1e1211SMatthew G Knepley       }
224cb1e1211SMatthew G Knepley     }
225*72f94c41SMatthew G. Knepley     ierr = DMPlexVecSetClosure(dm, section, localX, c, values, mode);CHKERRQ(ierr);
226cb1e1211SMatthew G Knepley   }
227*72f94c41SMatthew G. Knepley   ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
228*72f94c41SMatthew G. Knepley   ierr = PetscFree(sp);CHKERRQ(ierr);
229cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
230cb1e1211SMatthew G Knepley }
231cb1e1211SMatthew G Knepley 
232cb1e1211SMatthew G Knepley #undef __FUNCT__
233cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexProjectFunction"
234cb1e1211SMatthew G Knepley /*@C
235cb1e1211SMatthew G Knepley   DMPlexProjectFunction - This projects the given function into the function space provided.
236cb1e1211SMatthew G Knepley 
237cb1e1211SMatthew G Knepley   Input Parameters:
238cb1e1211SMatthew G Knepley + dm      - The DM
239*72f94c41SMatthew G. Knepley . fe      - The PetscFE associated with the field
240*72f94c41SMatthew G. Knepley . funcs   - The coordinate functions to evaluate, one per field
241cb1e1211SMatthew G Knepley - mode    - The insertion mode for values
242cb1e1211SMatthew G Knepley 
243cb1e1211SMatthew G Knepley   Output Parameter:
244cb1e1211SMatthew G Knepley . X - vector
245cb1e1211SMatthew G Knepley 
246cb1e1211SMatthew G Knepley   Level: developer
247cb1e1211SMatthew G Knepley 
248878cb397SSatish Balay .seealso: DMPlexComputeL2Diff()
249878cb397SSatish Balay @*/
250*72f94c41SMatthew G. Knepley PetscErrorCode DMPlexProjectFunction(DM dm, PetscFE fe[], void (**funcs)(const PetscReal [], PetscScalar *), InsertMode mode, Vec X)
251cb1e1211SMatthew G Knepley {
252cb1e1211SMatthew G Knepley   Vec            localX;
253cb1e1211SMatthew G Knepley   PetscErrorCode ierr;
254cb1e1211SMatthew G Knepley 
255cb1e1211SMatthew G Knepley   PetscFunctionBegin;
256cb1e1211SMatthew G Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
257*72f94c41SMatthew G. Knepley   ierr = DMPlexProjectFunctionLocal(dm, fe, funcs, mode, localX);CHKERRQ(ierr);
258cb1e1211SMatthew G Knepley   ierr = DMLocalToGlobalBegin(dm, localX, mode, X);CHKERRQ(ierr);
259cb1e1211SMatthew G Knepley   ierr = DMLocalToGlobalEnd(dm, localX, mode, X);CHKERRQ(ierr);
260cb1e1211SMatthew G Knepley   ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
261cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
262cb1e1211SMatthew G Knepley }
263cb1e1211SMatthew G Knepley 
264cb1e1211SMatthew G Knepley #undef __FUNCT__
265cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexComputeL2Diff"
266cb1e1211SMatthew G Knepley /*@C
267cb1e1211SMatthew G Knepley   DMPlexComputeL2Diff - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h.
268cb1e1211SMatthew G Knepley 
269cb1e1211SMatthew G Knepley   Input Parameters:
270cb1e1211SMatthew G Knepley + dm    - The DM
271c5bbbd5bSMatthew G. Knepley . fe    - The PetscFE object for each field
272cb1e1211SMatthew G Knepley . funcs - The functions to evaluate for each field component
273cb1e1211SMatthew G Knepley - X     - The coefficient vector u_h
274cb1e1211SMatthew G Knepley 
275cb1e1211SMatthew G Knepley   Output Parameter:
276cb1e1211SMatthew G Knepley . diff - The diff ||u - u_h||_2
277cb1e1211SMatthew G Knepley 
278cb1e1211SMatthew G Knepley   Level: developer
279cb1e1211SMatthew G Knepley 
280cb1e1211SMatthew G Knepley .seealso: DMPlexProjectFunction()
281878cb397SSatish Balay @*/
282c5bbbd5bSMatthew G. Knepley PetscErrorCode DMPlexComputeL2Diff(DM dm, PetscFE fe[], void (**funcs)(const PetscReal [], PetscScalar *), Vec X, PetscReal *diff)
283cb1e1211SMatthew G Knepley {
284cb1e1211SMatthew G Knepley   const PetscInt  debug = 0;
285cb1e1211SMatthew G Knepley   PetscSection    section;
286c5bbbd5bSMatthew G. Knepley   PetscQuadrature quad;
287cb1e1211SMatthew G Knepley   Vec             localX;
288*72f94c41SMatthew G. Knepley   PetscScalar    *funcVal;
289cb1e1211SMatthew G Knepley   PetscReal      *coords, *v0, *J, *invJ, detJ;
290cb1e1211SMatthew G Knepley   PetscReal       localDiff = 0.0;
291cb1e1211SMatthew G Knepley   PetscInt        dim, numFields, numComponents = 0, cStart, cEnd, c, field, fieldOffset, comp;
292cb1e1211SMatthew G Knepley   PetscErrorCode  ierr;
293cb1e1211SMatthew G Knepley 
294cb1e1211SMatthew G Knepley   PetscFunctionBegin;
295cb1e1211SMatthew G Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
296cb1e1211SMatthew G Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
297cb1e1211SMatthew G Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
298cb1e1211SMatthew G Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
299cb1e1211SMatthew G Knepley   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
300cb1e1211SMatthew G Knepley   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
301cb1e1211SMatthew G Knepley   for (field = 0; field < numFields; ++field) {
302c5bbbd5bSMatthew G. Knepley     PetscInt Nc;
303c5bbbd5bSMatthew G. Knepley 
304c5bbbd5bSMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe[field], &Nc);CHKERRQ(ierr);
305c5bbbd5bSMatthew G. Knepley     numComponents += Nc;
306cb1e1211SMatthew G Knepley   }
307*72f94c41SMatthew G. Knepley   ierr = DMPlexProjectFunctionLocal(dm, fe, funcs, INSERT_BC_VALUES, localX);CHKERRQ(ierr);
308*72f94c41SMatthew G. Knepley   ierr = PetscMalloc5(numComponents,PetscScalar,&funcVal,dim,PetscReal,&coords,dim,PetscReal,&v0,dim*dim,PetscReal,&J,dim*dim,PetscReal,&invJ);CHKERRQ(ierr);
309cb1e1211SMatthew G Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
310c5bbbd5bSMatthew G. Knepley   ierr = PetscFEGetQuadrature(fe[0], &quad);CHKERRQ(ierr);
311cb1e1211SMatthew G Knepley   for (c = cStart; c < cEnd; ++c) {
312a1e44745SMatthew G. Knepley     PetscScalar *x = NULL;
313cb1e1211SMatthew G Knepley     PetscReal    elemDiff = 0.0;
314cb1e1211SMatthew G Knepley 
315cb1e1211SMatthew G Knepley     ierr = DMPlexComputeCellGeometry(dm, c, v0, J, invJ, &detJ);CHKERRQ(ierr);
316cb1e1211SMatthew G Knepley     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
317cb1e1211SMatthew G Knepley     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
318cb1e1211SMatthew G Knepley 
319cb1e1211SMatthew G Knepley     for (field = 0, comp = 0, fieldOffset = 0; field < numFields; ++field) {
320c5bbbd5bSMatthew G. Knepley       const PetscInt   numQuadPoints = quad.numQuadPoints;
321c5bbbd5bSMatthew G. Knepley       const PetscReal *quadPoints    = quad.quadPoints;
322c5bbbd5bSMatthew G. Knepley       const PetscReal *quadWeights   = quad.quadWeights;
323c5bbbd5bSMatthew G. Knepley       PetscReal       *basis;
324c5bbbd5bSMatthew G. Knepley       PetscInt         numBasisFuncs, numBasisComps, q, d, e, fc, f;
325cb1e1211SMatthew G Knepley 
326c5bbbd5bSMatthew G. Knepley       ierr = PetscFEGetDimension(fe[field], &numBasisFuncs);CHKERRQ(ierr);
327c5bbbd5bSMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe[field], &numBasisComps);CHKERRQ(ierr);
328c5bbbd5bSMatthew G. Knepley       ierr = PetscFEGetDefaultTabulation(fe[field], &basis, NULL, NULL);CHKERRQ(ierr);
329cb1e1211SMatthew G Knepley       if (debug) {
330cb1e1211SMatthew G Knepley         char title[1024];
331cb1e1211SMatthew G Knepley         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
332cb1e1211SMatthew G Knepley         ierr = DMPrintCellVector(c, title, numBasisFuncs*numBasisComps, &x[fieldOffset]);CHKERRQ(ierr);
333cb1e1211SMatthew G Knepley       }
334cb1e1211SMatthew G Knepley       for (q = 0; q < numQuadPoints; ++q) {
335cb1e1211SMatthew G Knepley         for (d = 0; d < dim; d++) {
336cb1e1211SMatthew G Knepley           coords[d] = v0[d];
337cb1e1211SMatthew G Knepley           for (e = 0; e < dim; e++) {
338cb1e1211SMatthew G Knepley             coords[d] += J[d*dim+e]*(quadPoints[q*dim+e] + 1.0);
339cb1e1211SMatthew G Knepley           }
340cb1e1211SMatthew G Knepley         }
341*72f94c41SMatthew G. Knepley         (*funcs[field])(coords, funcVal);
342cb1e1211SMatthew G Knepley         for (fc = 0; fc < numBasisComps; ++fc) {
343a1d24da5SMatthew G. Knepley           PetscScalar interpolant = 0.0;
344a1d24da5SMatthew G. Knepley 
345cb1e1211SMatthew G Knepley           for (f = 0; f < numBasisFuncs; ++f) {
346cb1e1211SMatthew G Knepley             const PetscInt fidx = f*numBasisComps+fc;
347a1d24da5SMatthew G. Knepley             interpolant += x[fieldOffset+fidx]*basis[q*numBasisFuncs*numBasisComps+fidx];
348cb1e1211SMatthew G Knepley           }
349*72f94c41SMatthew 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);}
350*72f94c41SMatthew G. Knepley           elemDiff += PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ;
351cb1e1211SMatthew G Knepley         }
352cb1e1211SMatthew G Knepley       }
353cb1e1211SMatthew G Knepley       comp        += numBasisComps;
354cb1e1211SMatthew G Knepley       fieldOffset += numBasisFuncs*numBasisComps;
355cb1e1211SMatthew G Knepley     }
356cb1e1211SMatthew G Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
357cb1e1211SMatthew G Knepley     if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);}
358cb1e1211SMatthew G Knepley     localDiff += elemDiff;
359cb1e1211SMatthew G Knepley   }
360*72f94c41SMatthew G. Knepley   ierr  = PetscFree5(funcVal,coords,v0,J,invJ);CHKERRQ(ierr);
361cb1e1211SMatthew G Knepley   ierr  = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
36286a74ee0SMatthew G. Knepley   ierr  = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
363cb1e1211SMatthew G Knepley   *diff = PetscSqrtReal(*diff);
364cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
365cb1e1211SMatthew G Knepley }
366cb1e1211SMatthew G Knepley 
367a0845e3aSMatthew G. Knepley #if 0
368a0845e3aSMatthew G. Knepley 
369cb1e1211SMatthew G Knepley #undef __FUNCT__
370cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexComputeResidualFEM"
371cb1e1211SMatthew G Knepley PetscErrorCode DMPlexComputeResidualFEM(DM dm, Vec X, Vec F, void *user)
372cb1e1211SMatthew G Knepley {
373cb1e1211SMatthew G Knepley   DM_Plex         *mesh   = (DM_Plex *) dm->data;
3749a559087SMatthew G. Knepley   PetscFEM        *fem    = (PetscFEM *) user;
375cb1e1211SMatthew G Knepley   PetscQuadrature *quad   = fem->quad;
376652b88e8SMatthew G. Knepley   PetscQuadrature *quadBd = fem->quadBd;
377cb1e1211SMatthew G Knepley   PetscSection     section;
37802a80efeSMatthew G. Knepley   PetscReal       *v0, *n, *J, *invJ, *detJ;
379cb1e1211SMatthew G Knepley   PetscScalar     *elemVec, *u;
380cb1e1211SMatthew G Knepley   PetscInt         dim, numFields, field, numBatchesTmp = 1, numCells, cStart, cEnd, c;
381652b88e8SMatthew G. Knepley   PetscInt         cellDof, numComponents;
382652b88e8SMatthew G. Knepley   PetscBool        has;
383cb1e1211SMatthew G Knepley   PetscErrorCode   ierr;
384cb1e1211SMatthew G Knepley 
385cb1e1211SMatthew G Knepley   PetscFunctionBegin;
386652b88e8SMatthew G. Knepley   if (has && quadBd) {
387652b88e8SMatthew G. Knepley     DMLabel         label;
388652b88e8SMatthew G. Knepley     IS              pointIS;
389652b88e8SMatthew G. Knepley     const PetscInt *points;
390652b88e8SMatthew G. Knepley     PetscInt        numPoints, p;
391652b88e8SMatthew G. Knepley 
392652b88e8SMatthew G. Knepley     ierr = DMPlexGetLabel(dm, "boundary", &label);CHKERRQ(ierr);
393652b88e8SMatthew G. Knepley     ierr = DMLabelGetStratumSize(label, 1, &numPoints);CHKERRQ(ierr);
394652b88e8SMatthew G. Knepley     ierr = DMLabelGetStratumIS(label, 1, &pointIS);CHKERRQ(ierr);
395652b88e8SMatthew G. Knepley     ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
396652b88e8SMatthew G. Knepley     for (field = 0, cellDof = 0, numComponents = 0; field < numFields; ++field) {
397652b88e8SMatthew G. Knepley       cellDof       += quadBd[field].numBasisFuncs*quadBd[field].numComponents;
398652b88e8SMatthew G. Knepley       numComponents += quadBd[field].numComponents;
399652b88e8SMatthew G. Knepley     }
40002a80efeSMatthew G. Knepley     ierr = PetscMalloc7(numPoints*cellDof,PetscScalar,&u,numPoints*dim,PetscReal,&v0,numPoints*dim,PetscReal,&n,numPoints*dim*dim,PetscReal,&J,numPoints*dim*dim,PetscReal,&invJ,numPoints,PetscReal,&detJ,numPoints*cellDof,PetscScalar,&elemVec);CHKERRQ(ierr);
401652b88e8SMatthew G. Knepley     for (p = 0; p < numPoints; ++p) {
402652b88e8SMatthew G. Knepley       const PetscInt point = points[p];
403652b88e8SMatthew G. Knepley       PetscScalar   *x;
404652b88e8SMatthew G. Knepley       PetscInt       i;
405652b88e8SMatthew G. Knepley 
40602a80efeSMatthew G. Knepley       /* TODO: Add normal determination here */
407652b88e8SMatthew G. Knepley       ierr = DMPlexComputeCellGeometry(dm, point, &v0[p*dim], &J[p*dim*dim], &invJ[p*dim*dim], &detJ[p]);CHKERRQ(ierr);
4081d930511SMatthew G. Knepley       if (detJ[p] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for face %d", detJ[p], point);
409652b88e8SMatthew G. Knepley       ierr = DMPlexVecGetClosure(dm, NULL, X, point, NULL, &x);CHKERRQ(ierr);
410652b88e8SMatthew G. Knepley 
411652b88e8SMatthew G. Knepley       for (i = 0; i < cellDof; ++i) u[p*cellDof+i] = x[i];
412652b88e8SMatthew G. Knepley       ierr = DMPlexVecRestoreClosure(dm, NULL, X, point, NULL, &x);CHKERRQ(ierr);
413652b88e8SMatthew G. Knepley     }
414652b88e8SMatthew G. Knepley     for (field = 0; field < numFields; ++field) {
415652b88e8SMatthew G. Knepley       const PetscInt numQuadPoints = quadBd[field].numQuadPoints;
416652b88e8SMatthew G. Knepley       const PetscInt numBasisFuncs = quadBd[field].numBasisFuncs;
417a9dc2124SMatthew G. Knepley       void           (*f0)(const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->f0BdFuncs[field];
418a9dc2124SMatthew G. Knepley       void           (*f1)(const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->f1BdFuncs[field];
419652b88e8SMatthew G. Knepley       /* Conforming batches */
420652b88e8SMatthew G. Knepley       PetscInt blockSize  = numBasisFuncs*numQuadPoints;
421652b88e8SMatthew G. Knepley       PetscInt numBlocks  = 1;
422652b88e8SMatthew G. Knepley       PetscInt batchSize  = numBlocks * blockSize;
423652b88e8SMatthew G. Knepley       PetscInt numBatches = numBatchesTmp;
424652b88e8SMatthew G. Knepley       PetscInt numChunks  = numPoints / (numBatches*batchSize);
425652b88e8SMatthew G. Knepley       /* Remainder */
426652b88e8SMatthew G. Knepley       PetscInt numRemainder = numPoints % (numBatches * batchSize);
427652b88e8SMatthew G. Knepley       PetscInt offset       = numPoints - numRemainder;
428652b88e8SMatthew G. Knepley 
42902a80efeSMatthew G. Knepley       ierr = (*mesh->integrateBdResidualFEM)(numChunks*numBatches*batchSize, numFields, field, quadBd, u, v0, n, J, invJ, detJ, f0, f1, elemVec);CHKERRQ(ierr);
43002a80efeSMatthew 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],
431652b88e8SMatthew G. Knepley                                              f0, f1, &elemVec[offset*cellDof]);CHKERRQ(ierr);
432652b88e8SMatthew G. Knepley     }
433652b88e8SMatthew G. Knepley     for (p = 0; p < numPoints; ++p) {
434652b88e8SMatthew G. Knepley       const PetscInt point = points[p];
435652b88e8SMatthew G. Knepley 
436652b88e8SMatthew G. Knepley       if (mesh->printFEM > 1) {ierr = DMPrintCellVector(point, "Residual", cellDof, &elemVec[p*cellDof]);CHKERRQ(ierr);}
437652b88e8SMatthew G. Knepley       ierr = DMPlexVecSetClosure(dm, NULL, F, point, &elemVec[p*cellDof], ADD_VALUES);CHKERRQ(ierr);
438652b88e8SMatthew G. Knepley     }
439652b88e8SMatthew G. Knepley     ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
440652b88e8SMatthew G. Knepley     ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
44102a80efeSMatthew G. Knepley     ierr = PetscFree7(u,v0,n,J,invJ,detJ,elemVec);CHKERRQ(ierr);
442652b88e8SMatthew G. Knepley   }
443cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
444cb1e1211SMatthew G Knepley }
445cb1e1211SMatthew G Knepley 
446a0845e3aSMatthew G. Knepley #else
447a0845e3aSMatthew G. Knepley 
448a0845e3aSMatthew G. Knepley #undef __FUNCT__
449a0845e3aSMatthew G. Knepley #define __FUNCT__ "DMPlexComputeResidualFEM"
450a0845e3aSMatthew G. Knepley /*@
451a0845e3aSMatthew G. Knepley   DMPlexComputeResidualFEM - Form the local residual F from the local input X using pointwise functions specified by the user
452a0845e3aSMatthew G. Knepley 
453a0845e3aSMatthew G. Knepley   Input Parameters:
454a0845e3aSMatthew G. Knepley + dm - The mesh
455a0845e3aSMatthew G. Knepley . X  - Local input vector
456a0845e3aSMatthew G. Knepley - user - The user context
457a0845e3aSMatthew G. Knepley 
458a0845e3aSMatthew G. Knepley   Output Parameter:
459a0845e3aSMatthew G. Knepley . F  - Local output vector
460a0845e3aSMatthew G. Knepley 
461a0845e3aSMatthew G. Knepley   Note:
462a0845e3aSMatthew G. Knepley   The second member of the user context must be an FEMContext.
463a0845e3aSMatthew G. Knepley 
464a0845e3aSMatthew G. Knepley   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
465a0845e3aSMatthew G. Knepley   like a GPU, or vectorize on a multicore machine.
466a0845e3aSMatthew G. Knepley 
467a0845e3aSMatthew G. Knepley   Level: developer
468a0845e3aSMatthew G. Knepley 
469a0845e3aSMatthew G. Knepley .seealso: DMPlexComputeJacobianActionFEM()
470a0845e3aSMatthew G. Knepley @*/
471a0845e3aSMatthew G. Knepley PetscErrorCode DMPlexComputeResidualFEM(DM dm, Vec X, Vec F, void *user)
472a0845e3aSMatthew G. Knepley {
473a0845e3aSMatthew G. Knepley   DM_Plex          *mesh  = (DM_Plex *) dm->data;
4749a559087SMatthew G. Knepley   PetscFEM         *fem   = (PetscFEM *) user;
475a0845e3aSMatthew G. Knepley   PetscFE          *fe    = fem->fe;
4769a559087SMatthew G. Knepley   PetscFE          *feAux = fem->feAux;
477f1ea0e2fSMatthew G. Knepley   PetscFE          *feBd  = fem->feBd;
478a0845e3aSMatthew G. Knepley   const char       *name  = "Residual";
4799a559087SMatthew G. Knepley   DM                dmAux;
4809a559087SMatthew G. Knepley   Vec               A;
481a0845e3aSMatthew G. Knepley   PetscQuadrature   q;
482a0845e3aSMatthew G. Knepley   PetscCellGeometry geom;
4839a559087SMatthew G. Knepley   PetscSection      section, sectionAux;
484a0845e3aSMatthew G. Knepley   PetscReal        *v0, *J, *invJ, *detJ;
4859a559087SMatthew G. Knepley   PetscScalar      *elemVec, *u, *a;
4869a559087SMatthew G. Knepley   PetscInt          dim, Nf, NfAux = 0, f, numCells, cStart, cEnd, c;
487a0845e3aSMatthew G. Knepley   PetscInt          cellDof = 0, numComponents = 0;
4889a559087SMatthew G. Knepley   PetscInt          cellDofAux = 0, numComponentsAux = 0;
489a0845e3aSMatthew G. Knepley   PetscErrorCode    ierr;
490a0845e3aSMatthew G. Knepley 
491a0845e3aSMatthew G. Knepley   PetscFunctionBegin;
492a0845e3aSMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_ResidualFEM,dm,0,0,0);CHKERRQ(ierr);
493a0845e3aSMatthew G. Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
494a0845e3aSMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
4959a559087SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
496a0845e3aSMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
497a0845e3aSMatthew G. Knepley   numCells = cEnd - cStart;
4989a559087SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
499a0845e3aSMatthew G. Knepley     PetscInt Nb, Nc;
500a0845e3aSMatthew G. Knepley 
501a0845e3aSMatthew G. Knepley     ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr);
502a0845e3aSMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe[f], &Nc);CHKERRQ(ierr);
503a0845e3aSMatthew G. Knepley     cellDof       += Nb*Nc;
504a0845e3aSMatthew G. Knepley     numComponents += Nc;
505a0845e3aSMatthew G. Knepley   }
5069a559087SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
5079a559087SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr);
5089a559087SMatthew G. Knepley   if (dmAux) {
5099a559087SMatthew G. Knepley     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
5109a559087SMatthew G. Knepley     ierr = PetscSectionGetNumFields(sectionAux, &NfAux);CHKERRQ(ierr);
5119a559087SMatthew G. Knepley   }
5129a559087SMatthew G. Knepley   for (f = 0; f < NfAux; ++f) {
5139a559087SMatthew G. Knepley     PetscInt Nb, Nc;
5149a559087SMatthew G. Knepley 
5159a559087SMatthew G. Knepley     ierr = PetscFEGetDimension(feAux[f], &Nb);CHKERRQ(ierr);
5169a559087SMatthew G. Knepley     ierr = PetscFEGetNumComponents(feAux[f], &Nc);CHKERRQ(ierr);
5179a559087SMatthew G. Knepley     cellDofAux       += Nb*Nc;
5189a559087SMatthew G. Knepley     numComponentsAux += Nc;
5199a559087SMatthew G. Knepley   }
520*72f94c41SMatthew G. Knepley   ierr = DMPlexProjectFunctionLocal(dm, fe, fem->bcFuncs, INSERT_BC_VALUES, X);CHKERRQ(ierr);
521a0845e3aSMatthew G. Knepley   ierr = VecSet(F, 0.0);CHKERRQ(ierr);
522a0845e3aSMatthew G. Knepley   ierr = PetscMalloc6(numCells*cellDof,PetscScalar,&u,numCells*dim,PetscReal,&v0,numCells*dim*dim,PetscReal,&J,numCells*dim*dim,PetscReal,&invJ,numCells,PetscReal,&detJ,numCells*cellDof,PetscScalar,&elemVec);CHKERRQ(ierr);
5239a559087SMatthew G. Knepley   if (dmAux) {ierr = PetscMalloc(numCells*cellDofAux * sizeof(PetscScalar), &a);CHKERRQ(ierr);}
524a0845e3aSMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
525a0845e3aSMatthew G. Knepley     PetscScalar *x = NULL;
526a0845e3aSMatthew G. Knepley     PetscInt     i;
527a0845e3aSMatthew G. Knepley 
528a0845e3aSMatthew G. Knepley     ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr);
529a0845e3aSMatthew G. Knepley     if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c);
530a0845e3aSMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr);
531a0845e3aSMatthew G. Knepley     for (i = 0; i < cellDof; ++i) u[c*cellDof+i] = x[i];
532a0845e3aSMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr);
5339a559087SMatthew G. Knepley     if (dmAux) {
5349a559087SMatthew G. Knepley       ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
5359a559087SMatthew G. Knepley       for (i = 0; i < cellDofAux; ++i) a[c*cellDofAux+i] = x[i];
5369a559087SMatthew G. Knepley       ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
537a0845e3aSMatthew G. Knepley     }
5389a559087SMatthew G. Knepley   }
5399a559087SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
540c012ea0aSMatthew G. Knepley     void   (*f0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->f0Funcs[f];
541c012ea0aSMatthew G. Knepley     void   (*f1)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->f1Funcs[f];
542a0845e3aSMatthew G. Knepley     PetscInt Nb;
543a0845e3aSMatthew G. Knepley     /* Conforming batches */
544f30c5766SMatthew G. Knepley     PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
545a0845e3aSMatthew G. Knepley     /* Remainder */
546a0845e3aSMatthew G. Knepley     PetscInt Nr, offset;
547a0845e3aSMatthew G. Knepley 
548a0845e3aSMatthew G. Knepley     ierr = PetscFEGetQuadrature(fe[f], &q);CHKERRQ(ierr);
549a0845e3aSMatthew G. Knepley     ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr);
550f30c5766SMatthew G. Knepley     ierr = PetscFEGetTileSizes(fe[f], NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
551a0845e3aSMatthew G. Knepley     blockSize = Nb*q.numQuadPoints;
552a0845e3aSMatthew G. Knepley     batchSize = numBlocks * blockSize;
553f30c5766SMatthew G. Knepley     ierr =  PetscFESetTileSizes(fe[f], blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
554a0845e3aSMatthew G. Knepley     numChunks = numCells / (numBatches*batchSize);
555a0845e3aSMatthew G. Knepley     Ne        = numChunks*numBatches*batchSize;
556a0845e3aSMatthew G. Knepley     Nr        = numCells % (numBatches*batchSize);
557a0845e3aSMatthew G. Knepley     offset    = numCells - Nr;
558a0845e3aSMatthew G. Knepley     geom.v0   = v0;
559a0845e3aSMatthew G. Knepley     geom.J    = J;
560a0845e3aSMatthew G. Knepley     geom.invJ = invJ;
561a0845e3aSMatthew G. Knepley     geom.detJ = detJ;
5629a559087SMatthew G. Knepley     ierr = PetscFEIntegrateResidual(fe[f], Ne, Nf, fe, f, geom, u, NfAux, feAux, a, f0, f1, elemVec);CHKERRQ(ierr);
563a0845e3aSMatthew G. Knepley     geom.v0   = &v0[offset*dim];
564a0845e3aSMatthew G. Knepley     geom.J    = &J[offset*dim*dim];
565a0845e3aSMatthew G. Knepley     geom.invJ = &invJ[offset*dim*dim];
566a0845e3aSMatthew G. Knepley     geom.detJ = &detJ[offset];
5679a559087SMatthew 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);
568a0845e3aSMatthew G. Knepley   }
569a0845e3aSMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
570a0845e3aSMatthew G. Knepley     if (mesh->printFEM > 1) {ierr = DMPrintCellVector(c, name, cellDof, &elemVec[c*cellDof]);CHKERRQ(ierr);}
571a0845e3aSMatthew G. Knepley     ierr = DMPlexVecSetClosure(dm, section, F, c, &elemVec[c*cellDof], ADD_VALUES);CHKERRQ(ierr);
572a0845e3aSMatthew G. Knepley   }
573a0845e3aSMatthew G. Knepley   ierr = PetscFree6(u,v0,J,invJ,detJ,elemVec);CHKERRQ(ierr);
5749a559087SMatthew G. Knepley   if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);}
575f1ea0e2fSMatthew G. Knepley   if (feBd) {
576f1ea0e2fSMatthew G. Knepley     DMLabel         label;
577f1ea0e2fSMatthew G. Knepley     IS              pointIS;
578f1ea0e2fSMatthew G. Knepley     const PetscInt *points;
579f1ea0e2fSMatthew G. Knepley     PetscInt        numPoints, p;
580f1ea0e2fSMatthew G. Knepley     PetscReal      *n;
581f1ea0e2fSMatthew G. Knepley 
582f1ea0e2fSMatthew G. Knepley     ierr = DMPlexGetLabel(dm, "boundary", &label);CHKERRQ(ierr);
583f1ea0e2fSMatthew G. Knepley     ierr = DMLabelGetStratumSize(label, 1, &numPoints);CHKERRQ(ierr);
584f1ea0e2fSMatthew G. Knepley     ierr = DMLabelGetStratumIS(label, 1, &pointIS);CHKERRQ(ierr);
585f1ea0e2fSMatthew G. Knepley     ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
586f1ea0e2fSMatthew G. Knepley     for (f = 0, cellDof = 0, numComponents = 0; f < Nf; ++f) {
587f1ea0e2fSMatthew G. Knepley       PetscInt Nb, Nc;
588f1ea0e2fSMatthew G. Knepley 
589f1ea0e2fSMatthew G. Knepley       ierr = PetscFEGetDimension(feBd[f], &Nb);CHKERRQ(ierr);
590f1ea0e2fSMatthew G. Knepley       ierr = PetscFEGetNumComponents(feBd[f], &Nc);CHKERRQ(ierr);
591f1ea0e2fSMatthew G. Knepley       cellDof       += Nb*Nc;
592f1ea0e2fSMatthew G. Knepley       numComponents += Nc;
593f1ea0e2fSMatthew G. Knepley     }
594f1ea0e2fSMatthew G. Knepley     ierr = PetscMalloc7(numPoints*cellDof,PetscScalar,&u,numPoints*dim,PetscReal,&v0,numPoints*dim,PetscReal,&n,numPoints*dim*dim,PetscReal,&J,numPoints*dim*dim,PetscReal,&invJ,numPoints,PetscReal,&detJ,numPoints*cellDof,PetscScalar,&elemVec);CHKERRQ(ierr);
595f1ea0e2fSMatthew G. Knepley     for (p = 0; p < numPoints; ++p) {
596f1ea0e2fSMatthew G. Knepley       const PetscInt point = points[p];
597f1ea0e2fSMatthew G. Knepley       PetscScalar   *x     = NULL;
598f1ea0e2fSMatthew G. Knepley       PetscInt       i;
599f1ea0e2fSMatthew G. Knepley 
600f1ea0e2fSMatthew G. Knepley       /* TODO: Add normal determination here */
601f1ea0e2fSMatthew G. Knepley       ierr = DMPlexComputeCellGeometry(dm, point, &v0[p*dim], &J[p*dim*dim], &invJ[p*dim*dim], &detJ[p]);CHKERRQ(ierr);
602f1ea0e2fSMatthew G. Knepley       if (detJ[p] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for face %d", detJ[p], point);
603f1ea0e2fSMatthew G. Knepley       ierr = DMPlexVecGetClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr);
604f1ea0e2fSMatthew G. Knepley       for (i = 0; i < cellDof; ++i) u[p*cellDof+i] = x[i];
605f1ea0e2fSMatthew G. Knepley       ierr = DMPlexVecRestoreClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr);
606f1ea0e2fSMatthew G. Knepley     }
607f1ea0e2fSMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
608f1ea0e2fSMatthew G. Knepley       void   (*f0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->f0BdFuncs[f];
609f1ea0e2fSMatthew G. Knepley       void   (*f1)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->f1BdFuncs[f];
610f1ea0e2fSMatthew G. Knepley       PetscInt Nb;
611f1ea0e2fSMatthew G. Knepley       /* Conforming batches */
612f1ea0e2fSMatthew G. Knepley       PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
613f1ea0e2fSMatthew G. Knepley       /* Remainder */
614f1ea0e2fSMatthew G. Knepley       PetscInt Nr, offset;
615f1ea0e2fSMatthew G. Knepley 
616f1ea0e2fSMatthew G. Knepley       ierr = PetscFEGetQuadrature(feBd[f], &q);CHKERRQ(ierr);
617f1ea0e2fSMatthew G. Knepley       ierr = PetscFEGetDimension(feBd[f], &Nb);CHKERRQ(ierr);
618f1ea0e2fSMatthew G. Knepley       ierr = PetscFEGetTileSizes(feBd[f], NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
619f1ea0e2fSMatthew G. Knepley       blockSize = Nb*q.numQuadPoints;
620f1ea0e2fSMatthew G. Knepley       batchSize = numBlocks * blockSize;
621f1ea0e2fSMatthew G. Knepley       ierr =  PetscFESetTileSizes(feBd[f], blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
622f1ea0e2fSMatthew G. Knepley       numChunks = numPoints / (numBatches*batchSize);
623f1ea0e2fSMatthew G. Knepley       Ne        = numChunks*numBatches*batchSize;
624f1ea0e2fSMatthew G. Knepley       Nr        = numPoints % (numBatches*batchSize);
625f1ea0e2fSMatthew G. Knepley       offset    = numPoints - Nr;
626f1ea0e2fSMatthew G. Knepley       geom.v0   = v0;
627f1ea0e2fSMatthew G. Knepley       geom.n    = n;
628f1ea0e2fSMatthew G. Knepley       geom.J    = J;
629f1ea0e2fSMatthew G. Knepley       geom.invJ = invJ;
630f1ea0e2fSMatthew G. Knepley       geom.detJ = detJ;
631f1ea0e2fSMatthew G. Knepley       ierr = PetscFEIntegrateBdResidual(feBd[f], Ne, Nf, feBd, f, geom, u, 0, NULL, NULL, f0, f1, elemVec);CHKERRQ(ierr);
632f1ea0e2fSMatthew G. Knepley       geom.v0   = &v0[offset*dim];
633f1ea0e2fSMatthew G. Knepley       geom.n    = &n[offset*dim];
634f1ea0e2fSMatthew G. Knepley       geom.J    = &J[offset*dim*dim];
635f1ea0e2fSMatthew G. Knepley       geom.invJ = &invJ[offset*dim*dim];
636f1ea0e2fSMatthew G. Knepley       geom.detJ = &detJ[offset];
637f1ea0e2fSMatthew G. Knepley       ierr = PetscFEIntegrateBdResidual(feBd[f], Nr, Nf, feBd, f, geom, &u[offset*cellDof], 0, NULL, NULL, f0, f1, &elemVec[offset*cellDof]);CHKERRQ(ierr);
638f1ea0e2fSMatthew G. Knepley     }
639f1ea0e2fSMatthew G. Knepley     for (p = 0; p < numPoints; ++p) {
640f1ea0e2fSMatthew G. Knepley       const PetscInt point = points[p];
641f1ea0e2fSMatthew G. Knepley 
642f1ea0e2fSMatthew G. Knepley       if (mesh->printFEM > 1) {ierr = DMPrintCellVector(point, "BdResidual", cellDof, &elemVec[p*cellDof]);CHKERRQ(ierr);}
643f1ea0e2fSMatthew G. Knepley       ierr = DMPlexVecSetClosure(dm, NULL, F, point, &elemVec[p*cellDof], ADD_VALUES);CHKERRQ(ierr);
644f1ea0e2fSMatthew G. Knepley     }
645f1ea0e2fSMatthew G. Knepley     ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
646f1ea0e2fSMatthew G. Knepley     ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
647f1ea0e2fSMatthew G. Knepley     ierr = PetscFree7(u,v0,n,J,invJ,detJ,elemVec);CHKERRQ(ierr);
648f1ea0e2fSMatthew G. Knepley   }
6496113b454SMatthew G. Knepley   if (mesh->printFEM) {ierr = DMPrintLocalVec(dm, name, mesh->printTol, F);CHKERRQ(ierr);}
650a0845e3aSMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_ResidualFEM,dm,0,0,0);CHKERRQ(ierr);
651a0845e3aSMatthew G. Knepley   PetscFunctionReturn(0);
652a0845e3aSMatthew G. Knepley }
653a0845e3aSMatthew G. Knepley 
654a0845e3aSMatthew G. Knepley #endif
655a0845e3aSMatthew G. Knepley 
656cb1e1211SMatthew G Knepley #undef __FUNCT__
657cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexComputeJacobianActionFEM"
658cb1e1211SMatthew G Knepley /*@C
659cb1e1211SMatthew G Knepley   DMPlexComputeJacobianActionFEM - Form the local action of Jacobian J(u) on the local input X using pointwise functions specified by the user
660cb1e1211SMatthew G Knepley 
661cb1e1211SMatthew G Knepley   Input Parameters:
662cb1e1211SMatthew G Knepley + dm - The mesh
663cb1e1211SMatthew G Knepley . J  - The Jacobian shell matrix
664cb1e1211SMatthew G Knepley . X  - Local input vector
665cb1e1211SMatthew G Knepley - user - The user context
666cb1e1211SMatthew G Knepley 
667cb1e1211SMatthew G Knepley   Output Parameter:
668cb1e1211SMatthew G Knepley . F  - Local output vector
669cb1e1211SMatthew G Knepley 
670cb1e1211SMatthew G Knepley   Note:
671cb1e1211SMatthew G Knepley   The second member of the user context must be an FEMContext.
672cb1e1211SMatthew G Knepley 
673cb1e1211SMatthew G Knepley   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
674cb1e1211SMatthew G Knepley   like a GPU, or vectorize on a multicore machine.
675cb1e1211SMatthew G Knepley 
6760059ad2aSSatish Balay   Level: developer
6770059ad2aSSatish Balay 
678cb1e1211SMatthew G Knepley .seealso: DMPlexComputeResidualFEM()
679878cb397SSatish Balay @*/
680cb1e1211SMatthew G Knepley PetscErrorCode DMPlexComputeJacobianActionFEM(DM dm, Mat Jac, Vec X, Vec F, void *user)
681cb1e1211SMatthew G Knepley {
682cb1e1211SMatthew G Knepley   DM_Plex          *mesh = (DM_Plex *) dm->data;
6839a559087SMatthew G. Knepley   PetscFEM         *fem  = (PetscFEM *) user;
6840483ade4SMatthew G. Knepley   PetscFE          *fe   = fem->fe;
6850483ade4SMatthew G. Knepley   PetscQuadrature   quad;
6860483ade4SMatthew G. Knepley   PetscCellGeometry geom;
687cb1e1211SMatthew G Knepley   PetscSection      section;
688cb1e1211SMatthew G Knepley   JacActionCtx     *jctx;
689cb1e1211SMatthew G Knepley   PetscReal        *v0, *J, *invJ, *detJ;
690cb1e1211SMatthew G Knepley   PetscScalar      *elemVec, *u, *a;
6910483ade4SMatthew G. Knepley   PetscInt          dim, numFields, field, numCells, cStart, cEnd, c;
692cb1e1211SMatthew G Knepley   PetscInt          cellDof = 0;
693cb1e1211SMatthew G Knepley   PetscErrorCode    ierr;
694cb1e1211SMatthew G Knepley 
695cb1e1211SMatthew G Knepley   PetscFunctionBegin;
6960483ade4SMatthew G. Knepley   /* ierr = PetscLogEventBegin(DMPLEX_JacobianActionFEM,dm,0,0,0);CHKERRQ(ierr); */
697cb1e1211SMatthew G Knepley   ierr = MatShellGetContext(Jac, &jctx);CHKERRQ(ierr);
698cb1e1211SMatthew G Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
699cb1e1211SMatthew G Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
700cb1e1211SMatthew G Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
701cb1e1211SMatthew G Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
702cb1e1211SMatthew G Knepley   numCells = cEnd - cStart;
703cb1e1211SMatthew G Knepley   for (field = 0; field < numFields; ++field) {
7040483ade4SMatthew G. Knepley     PetscInt Nb, Nc;
7050483ade4SMatthew G. Knepley 
7060483ade4SMatthew G. Knepley     ierr = PetscFEGetDimension(fe[field], &Nb);CHKERRQ(ierr);
7070483ade4SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe[field], &Nc);CHKERRQ(ierr);
7080483ade4SMatthew G. Knepley     cellDof += Nb*Nc;
709cb1e1211SMatthew G Knepley   }
710cb1e1211SMatthew G Knepley   ierr = VecSet(F, 0.0);CHKERRQ(ierr);
711cb1e1211SMatthew G Knepley   ierr = PetscMalloc7(numCells*cellDof,PetscScalar,&u,numCells*cellDof,PetscScalar,&a,numCells*dim,PetscReal,&v0,numCells*dim*dim,PetscReal,&J,numCells*dim*dim,PetscReal,&invJ,numCells,PetscReal,&detJ,numCells*cellDof,PetscScalar,&elemVec);CHKERRQ(ierr);
712cb1e1211SMatthew G Knepley   for (c = cStart; c < cEnd; ++c) {
713a1e44745SMatthew G. Knepley     PetscScalar *x = NULL;
714cb1e1211SMatthew G Knepley     PetscInt     i;
715cb1e1211SMatthew G Knepley 
716cb1e1211SMatthew G Knepley     ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr);
717cb1e1211SMatthew G Knepley     if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c);
718cb1e1211SMatthew G Knepley     ierr = DMPlexVecGetClosure(dm, NULL, jctx->u, c, NULL, &x);CHKERRQ(ierr);
719cb1e1211SMatthew G Knepley     for (i = 0; i < cellDof; ++i) u[c*cellDof+i] = x[i];
720cb1e1211SMatthew G Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, jctx->u, c, NULL, &x);CHKERRQ(ierr);
721cb1e1211SMatthew G Knepley     ierr = DMPlexVecGetClosure(dm, NULL, X, c, NULL, &x);CHKERRQ(ierr);
722cb1e1211SMatthew G Knepley     for (i = 0; i < cellDof; ++i) a[c*cellDof+i] = x[i];
723cb1e1211SMatthew G Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, X, c, NULL, &x);CHKERRQ(ierr);
724cb1e1211SMatthew G Knepley   }
725cb1e1211SMatthew G Knepley   for (field = 0; field < numFields; ++field) {
7260483ade4SMatthew G. Knepley     PetscInt Nb;
727cb1e1211SMatthew G Knepley     /* Conforming batches */
728cb1e1211SMatthew G Knepley     PetscInt numBlocks  = 1;
7290483ade4SMatthew G. Knepley     PetscInt numBatches = 1;
7300483ade4SMatthew G. Knepley     PetscInt numChunks, Ne, blockSize, batchSize;
731cb1e1211SMatthew G Knepley     /* Remainder */
7320483ade4SMatthew G. Knepley     PetscInt Nr, offset;
733cb1e1211SMatthew G Knepley 
7340483ade4SMatthew G. Knepley     ierr = PetscFEGetQuadrature(fe[field], &quad);CHKERRQ(ierr);
7350483ade4SMatthew G. Knepley     ierr = PetscFEGetDimension(fe[field], &Nb);CHKERRQ(ierr);
7360483ade4SMatthew G. Knepley     blockSize = Nb*quad.numQuadPoints;
7370483ade4SMatthew G. Knepley     batchSize = numBlocks * blockSize;
7380483ade4SMatthew G. Knepley     numChunks = numCells / (numBatches*batchSize);
7390483ade4SMatthew G. Knepley     Ne        = numChunks*numBatches*batchSize;
7400483ade4SMatthew G. Knepley     Nr        = numCells % (numBatches*batchSize);
7410483ade4SMatthew G. Knepley     offset    = numCells - Nr;
7420483ade4SMatthew G. Knepley     geom.v0   = v0;
7430483ade4SMatthew G. Knepley     geom.J    = J;
7440483ade4SMatthew G. Knepley     geom.invJ = invJ;
7450483ade4SMatthew G. Knepley     geom.detJ = detJ;
7460483ade4SMatthew G. Knepley     ierr = PetscFEIntegrateJacobianAction(fe[field], Ne, numFields, fe, field, geom, u, a, fem->g0Funcs, fem->g1Funcs, fem->g2Funcs, fem->g3Funcs, elemVec);CHKERRQ(ierr);
7470483ade4SMatthew G. Knepley     geom.v0   = &v0[offset*dim];
7480483ade4SMatthew G. Knepley     geom.J    = &J[offset*dim*dim];
7490483ade4SMatthew G. Knepley     geom.invJ = &invJ[offset*dim*dim];
7500483ade4SMatthew G. Knepley     geom.detJ = &detJ[offset];
7510483ade4SMatthew G. Knepley     ierr = PetscFEIntegrateJacobianAction(fe[field], Nr, numFields, fe, field, geom, &u[offset*cellDof], &a[offset*cellDof],
752cb1e1211SMatthew G Knepley                                           fem->g0Funcs, fem->g1Funcs, fem->g2Funcs, fem->g3Funcs, &elemVec[offset*cellDof]);CHKERRQ(ierr);
753cb1e1211SMatthew G Knepley   }
754cb1e1211SMatthew G Knepley   for (c = cStart; c < cEnd; ++c) {
755cb1e1211SMatthew G Knepley     if (mesh->printFEM > 1) {ierr = DMPrintCellVector(c, "Jacobian Action", cellDof, &elemVec[c*cellDof]);CHKERRQ(ierr);}
756cb1e1211SMatthew G Knepley     ierr = DMPlexVecSetClosure(dm, NULL, F, c, &elemVec[c*cellDof], ADD_VALUES);CHKERRQ(ierr);
757cb1e1211SMatthew G Knepley   }
758cb1e1211SMatthew G Knepley   ierr = PetscFree7(u,a,v0,J,invJ,detJ,elemVec);CHKERRQ(ierr);
759cb1e1211SMatthew G Knepley   if (mesh->printFEM) {
760cb1e1211SMatthew G Knepley     PetscMPIInt rank, numProcs;
761cb1e1211SMatthew G Knepley     PetscInt    p;
762cb1e1211SMatthew G Knepley 
763cb1e1211SMatthew G Knepley     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr);
764cb1e1211SMatthew G Knepley     ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm), &numProcs);CHKERRQ(ierr);
76586a74ee0SMatthew G. Knepley     ierr = PetscPrintf(PetscObjectComm((PetscObject)dm), "Jacobian Action:\n");CHKERRQ(ierr);
766cb1e1211SMatthew G Knepley     for (p = 0; p < numProcs; ++p) {
767cb1e1211SMatthew G Knepley       if (p == rank) {ierr = VecView(F, PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);}
768cb1e1211SMatthew G Knepley       ierr = PetscBarrier((PetscObject) dm);CHKERRQ(ierr);
769cb1e1211SMatthew G Knepley     }
770cb1e1211SMatthew G Knepley   }
7710483ade4SMatthew G. Knepley   /* ierr = PetscLogEventEnd(DMPLEX_JacobianActionFEM,dm,0,0,0);CHKERRQ(ierr); */
772cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
773cb1e1211SMatthew G Knepley }
774cb1e1211SMatthew G Knepley 
775cb1e1211SMatthew G Knepley #undef __FUNCT__
776cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexComputeJacobianFEM"
777cb1e1211SMatthew G Knepley /*@
778cb1e1211SMatthew G Knepley   DMPlexComputeJacobianFEM - Form the local portion of the Jacobian matrix J at the local solution X using pointwise functions specified by the user.
779cb1e1211SMatthew G Knepley 
780cb1e1211SMatthew G Knepley   Input Parameters:
781cb1e1211SMatthew G Knepley + dm - The mesh
782cb1e1211SMatthew G Knepley . X  - Local input vector
783cb1e1211SMatthew G Knepley - user - The user context
784cb1e1211SMatthew G Knepley 
785cb1e1211SMatthew G Knepley   Output Parameter:
786cb1e1211SMatthew G Knepley . Jac  - Jacobian matrix
787cb1e1211SMatthew G Knepley 
788cb1e1211SMatthew G Knepley   Note:
789cb1e1211SMatthew G Knepley   The second member of the user context must be an FEMContext.
790cb1e1211SMatthew G Knepley 
791cb1e1211SMatthew G Knepley   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
792cb1e1211SMatthew G Knepley   like a GPU, or vectorize on a multicore machine.
793cb1e1211SMatthew G Knepley 
7940059ad2aSSatish Balay   Level: developer
7950059ad2aSSatish Balay 
796cb1e1211SMatthew G Knepley .seealso: FormFunctionLocal()
797878cb397SSatish Balay @*/
798cb1e1211SMatthew G Knepley PetscErrorCode DMPlexComputeJacobianFEM(DM dm, Vec X, Mat Jac, Mat JacP, MatStructure *str,void *user)
799cb1e1211SMatthew G Knepley {
800cb1e1211SMatthew G Knepley   DM_Plex          *mesh  = (DM_Plex *) dm->data;
8019a559087SMatthew G. Knepley   PetscFEM         *fem   = (PetscFEM *) user;
802a319912fSMatthew G. Knepley   PetscFE          *fe    = fem->fe;
803754551f4SMatthew G. Knepley   PetscFE          *feAux = fem->feAux;
804a319912fSMatthew G. Knepley   const char       *name  = "Jacobian";
805754551f4SMatthew G. Knepley   DM                dmAux;
806754551f4SMatthew G. Knepley   Vec               A;
807a319912fSMatthew G. Knepley   PetscQuadrature   quad;
808a319912fSMatthew G. Knepley   PetscCellGeometry geom;
809754551f4SMatthew G. Knepley   PetscSection      section, globalSection, sectionAux;
810cb1e1211SMatthew G Knepley   PetscReal        *v0, *J, *invJ, *detJ;
811754551f4SMatthew G. Knepley   PetscScalar      *elemMat, *u, *a;
812754551f4SMatthew G. Knepley   PetscInt          dim, Nf, NfAux = 0, f, fieldI, fieldJ, numCells, cStart, cEnd, c;
813cb1e1211SMatthew G Knepley   PetscInt          cellDof = 0, numComponents = 0;
814754551f4SMatthew G. Knepley   PetscInt          cellDofAux = 0, numComponentsAux = 0;
815cb1e1211SMatthew G Knepley   PetscBool         isShell;
816cb1e1211SMatthew G Knepley   PetscErrorCode    ierr;
817cb1e1211SMatthew G Knepley 
818cb1e1211SMatthew G Knepley   PetscFunctionBegin;
819a319912fSMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_JacobianFEM,dm,0,0,0);CHKERRQ(ierr);
820cb1e1211SMatthew G Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
821cb1e1211SMatthew G Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
822a319912fSMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dm, &globalSection);CHKERRQ(ierr);
823754551f4SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
824cb1e1211SMatthew G Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
825cb1e1211SMatthew G Knepley   numCells = cEnd - cStart;
826754551f4SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
827a319912fSMatthew G. Knepley     PetscInt Nb, Nc;
828a319912fSMatthew G. Knepley 
829a319912fSMatthew G. Knepley     ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr);
830a319912fSMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe[f], &Nc);CHKERRQ(ierr);
831a319912fSMatthew G. Knepley     cellDof       += Nb*Nc;
832a319912fSMatthew G. Knepley     numComponents += Nc;
833cb1e1211SMatthew G Knepley   }
834754551f4SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
835754551f4SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr);
836754551f4SMatthew G. Knepley   if (dmAux) {
837754551f4SMatthew G. Knepley     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
838754551f4SMatthew G. Knepley     ierr = PetscSectionGetNumFields(sectionAux, &NfAux);CHKERRQ(ierr);
839754551f4SMatthew G. Knepley   }
840754551f4SMatthew G. Knepley   for (f = 0; f < NfAux; ++f) {
841754551f4SMatthew G. Knepley     PetscInt Nb, Nc;
842754551f4SMatthew G. Knepley 
843754551f4SMatthew G. Knepley     ierr = PetscFEGetDimension(feAux[f], &Nb);CHKERRQ(ierr);
844754551f4SMatthew G. Knepley     ierr = PetscFEGetNumComponents(feAux[f], &Nc);CHKERRQ(ierr);
845754551f4SMatthew G. Knepley     cellDofAux       += Nb*Nc;
846754551f4SMatthew G. Knepley     numComponentsAux += Nc;
847754551f4SMatthew G. Knepley   }
848*72f94c41SMatthew G. Knepley   ierr = DMPlexProjectFunctionLocal(dm, fe, fem->bcFuncs, INSERT_BC_VALUES, X);CHKERRQ(ierr);
849cb1e1211SMatthew G Knepley   ierr = MatZeroEntries(JacP);CHKERRQ(ierr);
850cb1e1211SMatthew G Knepley   ierr = PetscMalloc6(numCells*cellDof,PetscScalar,&u,numCells*dim,PetscReal,&v0,numCells*dim*dim,PetscReal,&J,numCells*dim*dim,PetscReal,&invJ,numCells,PetscReal,&detJ,numCells*cellDof*cellDof,PetscScalar,&elemMat);CHKERRQ(ierr);
851754551f4SMatthew G. Knepley   if (dmAux) {ierr = PetscMalloc(numCells*cellDofAux * sizeof(PetscScalar), &a);CHKERRQ(ierr);}
852cb1e1211SMatthew G Knepley   for (c = cStart; c < cEnd; ++c) {
853a1e44745SMatthew G. Knepley     PetscScalar *x = NULL;
854cb1e1211SMatthew G Knepley     PetscInt     i;
855cb1e1211SMatthew G Knepley 
856cb1e1211SMatthew G Knepley     ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr);
857cb1e1211SMatthew G Knepley     if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c);
858a319912fSMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr);
859cb1e1211SMatthew G Knepley     for (i = 0; i < cellDof; ++i) u[c*cellDof+i] = x[i];
860a319912fSMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr);
861754551f4SMatthew G. Knepley     if (dmAux) {
862754551f4SMatthew G. Knepley       ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
863754551f4SMatthew G. Knepley       for (i = 0; i < cellDofAux; ++i) a[c*cellDofAux+i] = x[i];
864754551f4SMatthew G. Knepley       ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
865754551f4SMatthew G. Knepley     }
866cb1e1211SMatthew G Knepley   }
867cb1e1211SMatthew G Knepley   ierr = PetscMemzero(elemMat, numCells*cellDof*cellDof * sizeof(PetscScalar));CHKERRQ(ierr);
868754551f4SMatthew G. Knepley   for (fieldI = 0; fieldI < Nf; ++fieldI) {
869a319912fSMatthew G. Knepley     PetscInt Nb;
870cb1e1211SMatthew G Knepley     /* Conforming batches */
871754551f4SMatthew G. Knepley     PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
872cb1e1211SMatthew G Knepley     /* Remainder */
873a319912fSMatthew G. Knepley     PetscInt Nr, offset;
874cb1e1211SMatthew G Knepley 
875754551f4SMatthew G. Knepley     ierr = PetscFEGetQuadrature(fe[fieldI], &quad);CHKERRQ(ierr);
876754551f4SMatthew G. Knepley     ierr = PetscFEGetDimension(fe[fieldI], &Nb);CHKERRQ(ierr);
877754551f4SMatthew G. Knepley     ierr = PetscFEGetTileSizes(fe[fieldI], NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
878a319912fSMatthew G. Knepley     blockSize = Nb*quad.numQuadPoints;
879a319912fSMatthew G. Knepley     batchSize = numBlocks * blockSize;
880754551f4SMatthew G. Knepley     ierr = PetscFESetTileSizes(fe[fieldI], blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
881a319912fSMatthew G. Knepley     numChunks = numCells / (numBatches*batchSize);
882a319912fSMatthew G. Knepley     Ne        = numChunks*numBatches*batchSize;
883a319912fSMatthew G. Knepley     Nr        = numCells % (numBatches*batchSize);
884a319912fSMatthew G. Knepley     offset    = numCells - Nr;
885754551f4SMatthew G. Knepley     for (fieldJ = 0; fieldJ < Nf; ++fieldJ) {
886754551f4SMatthew G. Knepley       void   (*g0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g0Funcs[fieldI*Nf+fieldJ];
887754551f4SMatthew G. Knepley       void   (*g1)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g1Funcs[fieldI*Nf+fieldJ];
888754551f4SMatthew G. Knepley       void   (*g2)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g2Funcs[fieldI*Nf+fieldJ];
889754551f4SMatthew G. Knepley       void   (*g3)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g3Funcs[fieldI*Nf+fieldJ];
890754551f4SMatthew G. Knepley 
891a319912fSMatthew G. Knepley       geom.v0   = v0;
892a319912fSMatthew G. Knepley       geom.J    = J;
893a319912fSMatthew G. Knepley       geom.invJ = invJ;
894a319912fSMatthew G. Knepley       geom.detJ = detJ;
895754551f4SMatthew G. Knepley       ierr = PetscFEIntegrateJacobian(fe[fieldI], Ne, Nf, fe, fieldI, fieldJ, geom, u, NfAux, feAux, a, g0, g1, g2, g3, elemMat);CHKERRQ(ierr);
896a319912fSMatthew G. Knepley       geom.v0   = &v0[offset*dim];
897a319912fSMatthew G. Knepley       geom.J    = &J[offset*dim*dim];
898a319912fSMatthew G. Knepley       geom.invJ = &invJ[offset*dim*dim];
899a319912fSMatthew G. Knepley       geom.detJ = &detJ[offset];
900754551f4SMatthew 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);
901cb1e1211SMatthew G Knepley     }
902cb1e1211SMatthew G Knepley   }
903cb1e1211SMatthew G Knepley   for (c = cStart; c < cEnd; ++c) {
904a319912fSMatthew G. Knepley     if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(c, name, cellDof, cellDof, &elemMat[c*cellDof*cellDof]);CHKERRQ(ierr);}
905a319912fSMatthew G. Knepley     ierr = DMPlexMatSetClosure(dm, section, globalSection, JacP, c, &elemMat[c*cellDof*cellDof], ADD_VALUES);CHKERRQ(ierr);
906cb1e1211SMatthew G Knepley   }
907cb1e1211SMatthew G Knepley   ierr = PetscFree6(u,v0,J,invJ,detJ,elemMat);CHKERRQ(ierr);
908754551f4SMatthew G. Knepley   if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);}
909cb1e1211SMatthew G Knepley   ierr = MatAssemblyBegin(JacP, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
910cb1e1211SMatthew G Knepley   ierr = MatAssemblyEnd(JacP, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
911cb1e1211SMatthew G Knepley   if (mesh->printFEM) {
912a319912fSMatthew G. Knepley     ierr = PetscPrintf(PETSC_COMM_WORLD, "%s:\n", name);CHKERRQ(ierr);
913cb1e1211SMatthew G Knepley     ierr = MatChop(JacP, 1.0e-10);CHKERRQ(ierr);
914cb1e1211SMatthew G Knepley     ierr = MatView(JacP, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
915cb1e1211SMatthew G Knepley   }
916a319912fSMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_JacobianFEM,dm,0,0,0);CHKERRQ(ierr);
917cb1e1211SMatthew G Knepley   ierr = PetscObjectTypeCompare((PetscObject) Jac, MATSHELL, &isShell);CHKERRQ(ierr);
918cb1e1211SMatthew G Knepley   if (isShell) {
919cb1e1211SMatthew G Knepley     JacActionCtx *jctx;
920cb1e1211SMatthew G Knepley 
921cb1e1211SMatthew G Knepley     ierr = MatShellGetContext(Jac, &jctx);CHKERRQ(ierr);
922cb1e1211SMatthew G Knepley     ierr = VecCopy(X, jctx->u);CHKERRQ(ierr);
923cb1e1211SMatthew G Knepley   }
924cb1e1211SMatthew G Knepley   *str = SAME_NONZERO_PATTERN;
925cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
926cb1e1211SMatthew G Knepley }
927