xref: /petsc/src/dm/impls/plex/plexfem.c (revision a1d24da52eb909373f67ff3d47d8bc95ef57cc3d)
1cb1e1211SMatthew G Knepley #include <petsc-private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
2cb1e1211SMatthew G Knepley 
3cb1e1211SMatthew G Knepley #undef __FUNCT__
4cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexGetScale"
5cb1e1211SMatthew G Knepley PetscErrorCode DMPlexGetScale(DM dm, PetscUnit unit, PetscReal *scale)
6cb1e1211SMatthew G Knepley {
7cb1e1211SMatthew G Knepley   DM_Plex *mesh = (DM_Plex*) dm->data;
8cb1e1211SMatthew G Knepley 
9cb1e1211SMatthew G Knepley   PetscFunctionBegin;
10cb1e1211SMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
11cb1e1211SMatthew G Knepley   PetscValidPointer(scale, 3);
12cb1e1211SMatthew G Knepley   *scale = mesh->scale[unit];
13cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
14cb1e1211SMatthew G Knepley }
15cb1e1211SMatthew G Knepley 
16cb1e1211SMatthew G Knepley #undef __FUNCT__
17cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexSetScale"
18cb1e1211SMatthew G Knepley PetscErrorCode DMPlexSetScale(DM dm, PetscUnit unit, PetscReal scale)
19cb1e1211SMatthew G Knepley {
20cb1e1211SMatthew G Knepley   DM_Plex *mesh = (DM_Plex*) dm->data;
21cb1e1211SMatthew G Knepley 
22cb1e1211SMatthew G Knepley   PetscFunctionBegin;
23cb1e1211SMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
24cb1e1211SMatthew G Knepley   mesh->scale[unit] = scale;
25cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
26cb1e1211SMatthew G Knepley }
27cb1e1211SMatthew G Knepley 
28cb1e1211SMatthew G Knepley PETSC_STATIC_INLINE PetscInt epsilon(PetscInt i, PetscInt j, PetscInt k)
29cb1e1211SMatthew G Knepley {
30cb1e1211SMatthew G Knepley   switch (i) {
31cb1e1211SMatthew G Knepley   case 0:
32cb1e1211SMatthew G Knepley     switch (j) {
33cb1e1211SMatthew G Knepley     case 0: return 0;
34cb1e1211SMatthew G Knepley     case 1:
35cb1e1211SMatthew G Knepley       switch (k) {
36cb1e1211SMatthew G Knepley       case 0: return 0;
37cb1e1211SMatthew G Knepley       case 1: return 0;
38cb1e1211SMatthew G Knepley       case 2: return 1;
39cb1e1211SMatthew G Knepley       }
40cb1e1211SMatthew G Knepley     case 2:
41cb1e1211SMatthew G Knepley       switch (k) {
42cb1e1211SMatthew G Knepley       case 0: return 0;
43cb1e1211SMatthew G Knepley       case 1: return -1;
44cb1e1211SMatthew G Knepley       case 2: return 0;
45cb1e1211SMatthew G Knepley       }
46cb1e1211SMatthew G Knepley     }
47cb1e1211SMatthew G Knepley   case 1:
48cb1e1211SMatthew G Knepley     switch (j) {
49cb1e1211SMatthew G Knepley     case 0:
50cb1e1211SMatthew G Knepley       switch (k) {
51cb1e1211SMatthew G Knepley       case 0: return 0;
52cb1e1211SMatthew G Knepley       case 1: return 0;
53cb1e1211SMatthew G Knepley       case 2: return -1;
54cb1e1211SMatthew G Knepley       }
55cb1e1211SMatthew G Knepley     case 1: return 0;
56cb1e1211SMatthew G Knepley     case 2:
57cb1e1211SMatthew G Knepley       switch (k) {
58cb1e1211SMatthew G Knepley       case 0: return 1;
59cb1e1211SMatthew G Knepley       case 1: return 0;
60cb1e1211SMatthew G Knepley       case 2: return 0;
61cb1e1211SMatthew G Knepley       }
62cb1e1211SMatthew G Knepley     }
63cb1e1211SMatthew G Knepley   case 2:
64cb1e1211SMatthew G Knepley     switch (j) {
65cb1e1211SMatthew G Knepley     case 0:
66cb1e1211SMatthew G Knepley       switch (k) {
67cb1e1211SMatthew G Knepley       case 0: return 0;
68cb1e1211SMatthew G Knepley       case 1: return 1;
69cb1e1211SMatthew G Knepley       case 2: return 0;
70cb1e1211SMatthew G Knepley       }
71cb1e1211SMatthew G Knepley     case 1:
72cb1e1211SMatthew G Knepley       switch (k) {
73cb1e1211SMatthew G Knepley       case 0: return -1;
74cb1e1211SMatthew G Knepley       case 1: return 0;
75cb1e1211SMatthew G Knepley       case 2: return 0;
76cb1e1211SMatthew G Knepley       }
77cb1e1211SMatthew G Knepley     case 2: return 0;
78cb1e1211SMatthew G Knepley     }
79cb1e1211SMatthew G Knepley   }
80cb1e1211SMatthew G Knepley   return 0;
81cb1e1211SMatthew G Knepley }
82cb1e1211SMatthew G Knepley 
83cb1e1211SMatthew G Knepley #undef __FUNCT__
84cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexCreateRigidBody"
85cb1e1211SMatthew G Knepley /*@C
86cb1e1211SMatthew G Knepley   DMPlexCreateRigidBody - create rigid body modes from coordinates
87cb1e1211SMatthew G Knepley 
88cb1e1211SMatthew G Knepley   Collective on DM
89cb1e1211SMatthew G Knepley 
90cb1e1211SMatthew G Knepley   Input Arguments:
91cb1e1211SMatthew G Knepley + dm - the DM
92cb1e1211SMatthew G Knepley . section - the local section associated with the rigid field, or NULL for the default section
93cb1e1211SMatthew G Knepley - globalSection - the global section associated with the rigid field, or NULL for the default section
94cb1e1211SMatthew G Knepley 
95cb1e1211SMatthew G Knepley   Output Argument:
96cb1e1211SMatthew G Knepley . sp - the null space
97cb1e1211SMatthew G Knepley 
98cb1e1211SMatthew G Knepley   Note: This is necessary to take account of Dirichlet conditions on the displacements
99cb1e1211SMatthew G Knepley 
100cb1e1211SMatthew G Knepley   Level: advanced
101cb1e1211SMatthew G Knepley 
102cb1e1211SMatthew G Knepley .seealso: MatNullSpaceCreate()
103cb1e1211SMatthew G Knepley @*/
104cb1e1211SMatthew G Knepley PetscErrorCode DMPlexCreateRigidBody(DM dm, PetscSection section, PetscSection globalSection, MatNullSpace *sp)
105cb1e1211SMatthew G Knepley {
106cb1e1211SMatthew G Knepley   MPI_Comm       comm;
107cb1e1211SMatthew G Knepley   Vec            coordinates, localMode, mode[6];
108cb1e1211SMatthew G Knepley   PetscSection   coordSection;
109cb1e1211SMatthew G Knepley   PetscScalar   *coords;
110cb1e1211SMatthew G Knepley   PetscInt       dim, vStart, vEnd, v, n, m, d, i, j;
111cb1e1211SMatthew G Knepley   PetscErrorCode ierr;
112cb1e1211SMatthew G Knepley 
113cb1e1211SMatthew G Knepley   PetscFunctionBegin;
114cb1e1211SMatthew G Knepley   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
115cb1e1211SMatthew G Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
116cb1e1211SMatthew G Knepley   if (dim == 1) {
117cb1e1211SMatthew G Knepley     ierr = MatNullSpaceCreate(comm, PETSC_TRUE, 0, NULL, sp);CHKERRQ(ierr);
118cb1e1211SMatthew G Knepley     PetscFunctionReturn(0);
119cb1e1211SMatthew G Knepley   }
120cb1e1211SMatthew G Knepley   if (!section)       {ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);}
121cb1e1211SMatthew G Knepley   if (!globalSection) {ierr = DMGetDefaultGlobalSection(dm, &globalSection);CHKERRQ(ierr);}
122cb1e1211SMatthew G Knepley   ierr = PetscSectionGetConstrainedStorageSize(globalSection, &n);CHKERRQ(ierr);
123cb1e1211SMatthew G Knepley   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
124cb1e1211SMatthew G Knepley   ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
125cb1e1211SMatthew G Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
126cb1e1211SMatthew G Knepley   m    = (dim*(dim+1))/2;
127cb1e1211SMatthew G Knepley   ierr = VecCreate(comm, &mode[0]);CHKERRQ(ierr);
128cb1e1211SMatthew G Knepley   ierr = VecSetSizes(mode[0], n, PETSC_DETERMINE);CHKERRQ(ierr);
129cb1e1211SMatthew G Knepley   ierr = VecSetUp(mode[0]);CHKERRQ(ierr);
130cb1e1211SMatthew G Knepley   for (i = 1; i < m; ++i) {ierr = VecDuplicate(mode[0], &mode[i]);CHKERRQ(ierr);}
131cb1e1211SMatthew G Knepley   /* Assume P1 */
132cb1e1211SMatthew G Knepley   ierr = DMGetLocalVector(dm, &localMode);CHKERRQ(ierr);
133cb1e1211SMatthew G Knepley   for (d = 0; d < dim; ++d) {
134cb1e1211SMatthew G Knepley     PetscScalar values[3] = {0.0, 0.0, 0.0};
135cb1e1211SMatthew G Knepley 
136cb1e1211SMatthew G Knepley     values[d] = 1.0;
137cb1e1211SMatthew G Knepley     ierr      = VecSet(localMode, 0.0);CHKERRQ(ierr);
138cb1e1211SMatthew G Knepley     for (v = vStart; v < vEnd; ++v) {
139cb1e1211SMatthew G Knepley       ierr = DMPlexVecSetClosure(dm, section, localMode, v, values, INSERT_VALUES);CHKERRQ(ierr);
140cb1e1211SMatthew G Knepley     }
141cb1e1211SMatthew G Knepley     ierr = DMLocalToGlobalBegin(dm, localMode, INSERT_VALUES, mode[d]);CHKERRQ(ierr);
142cb1e1211SMatthew G Knepley     ierr = DMLocalToGlobalEnd(dm, localMode, INSERT_VALUES, mode[d]);CHKERRQ(ierr);
143cb1e1211SMatthew G Knepley   }
144cb1e1211SMatthew G Knepley   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
145cb1e1211SMatthew G Knepley   for (d = dim; d < dim*(dim+1)/2; ++d) {
146cb1e1211SMatthew G Knepley     PetscInt i, j, k = dim > 2 ? d - dim : d;
147cb1e1211SMatthew G Knepley 
148cb1e1211SMatthew G Knepley     ierr = VecSet(localMode, 0.0);CHKERRQ(ierr);
149cb1e1211SMatthew G Knepley     for (v = vStart; v < vEnd; ++v) {
150cb1e1211SMatthew G Knepley       PetscScalar values[3] = {0.0, 0.0, 0.0};
151cb1e1211SMatthew G Knepley       PetscInt    off;
152cb1e1211SMatthew G Knepley 
153cb1e1211SMatthew G Knepley       ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr);
154cb1e1211SMatthew G Knepley       for (i = 0; i < dim; ++i) {
155cb1e1211SMatthew G Knepley         for (j = 0; j < dim; ++j) {
156cb1e1211SMatthew G Knepley           values[j] += epsilon(i, j, k)*PetscRealPart(coords[off+i]);
157cb1e1211SMatthew G Knepley         }
158cb1e1211SMatthew G Knepley       }
159cb1e1211SMatthew G Knepley       ierr = DMPlexVecSetClosure(dm, section, localMode, v, values, INSERT_VALUES);CHKERRQ(ierr);
160cb1e1211SMatthew G Knepley     }
161cb1e1211SMatthew G Knepley     ierr = DMLocalToGlobalBegin(dm, localMode, INSERT_VALUES, mode[d]);CHKERRQ(ierr);
162cb1e1211SMatthew G Knepley     ierr = DMLocalToGlobalEnd(dm, localMode, INSERT_VALUES, mode[d]);CHKERRQ(ierr);
163cb1e1211SMatthew G Knepley   }
164cb1e1211SMatthew G Knepley   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
165cb1e1211SMatthew G Knepley   ierr = DMRestoreLocalVector(dm, &localMode);CHKERRQ(ierr);
166cb1e1211SMatthew G Knepley   for (i = 0; i < dim; ++i) {ierr = VecNormalize(mode[i], NULL);CHKERRQ(ierr);}
167cb1e1211SMatthew G Knepley   /* Orthonormalize system */
168cb1e1211SMatthew G Knepley   for (i = dim; i < m; ++i) {
169cb1e1211SMatthew G Knepley     PetscScalar dots[6];
170cb1e1211SMatthew G Knepley 
171cb1e1211SMatthew G Knepley     ierr = VecMDot(mode[i], i, mode, dots);CHKERRQ(ierr);
172cb1e1211SMatthew G Knepley     for (j = 0; j < i; ++j) dots[j] *= -1.0;
173cb1e1211SMatthew G Knepley     ierr = VecMAXPY(mode[i], i, dots, mode);CHKERRQ(ierr);
174cb1e1211SMatthew G Knepley     ierr = VecNormalize(mode[i], NULL);CHKERRQ(ierr);
175cb1e1211SMatthew G Knepley   }
176cb1e1211SMatthew G Knepley   ierr = MatNullSpaceCreate(comm, PETSC_FALSE, m, mode, sp);CHKERRQ(ierr);
177cb1e1211SMatthew G Knepley   for (i = 0; i< m; ++i) {ierr = VecDestroy(&mode[i]);CHKERRQ(ierr);}
178cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
179cb1e1211SMatthew G Knepley }
180cb1e1211SMatthew G Knepley /*******************************************************************************
181cb1e1211SMatthew G Knepley This should be in a separate Discretization object, but I am not sure how to lay
182cb1e1211SMatthew G Knepley it out yet, so I am stuffing things here while I experiment.
183cb1e1211SMatthew G Knepley *******************************************************************************/
184cb1e1211SMatthew G Knepley #undef __FUNCT__
185cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexSetFEMIntegration"
186cb1e1211SMatthew G Knepley PetscErrorCode DMPlexSetFEMIntegration(DM dm,
187cb1e1211SMatthew G Knepley                                           PetscErrorCode (*integrateResidualFEM)(PetscInt, PetscInt, PetscInt, PetscQuadrature[], const PetscScalar[],
188cb1e1211SMatthew G Knepley                                                                                  const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[],
189cb1e1211SMatthew G Knepley                                                                                  void (*)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]),
190cb1e1211SMatthew G Knepley                                                                                  void (*)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), PetscScalar[]),
19102a80efeSMatthew G. Knepley                                           PetscErrorCode (*integrateBdResidualFEM)(PetscInt, PetscInt, PetscInt, PetscQuadrature[], const PetscScalar[],
19202a80efeSMatthew G. Knepley                                                                                    const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[],
19302a80efeSMatthew G. Knepley                                                                                    void (*)(const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]),
19402a80efeSMatthew G. Knepley                                                                                    void (*)(const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]), PetscScalar[]),
195cb1e1211SMatthew G Knepley                                           PetscErrorCode (*integrateJacobianActionFEM)(PetscInt, PetscInt, PetscInt, PetscQuadrature[], const PetscScalar[], const PetscScalar[],
196cb1e1211SMatthew G Knepley                                                                                        const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[],
197cb1e1211SMatthew G Knepley                                                                                        void (**)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]),
198cb1e1211SMatthew G Knepley                                                                                        void (**)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]),
199cb1e1211SMatthew G Knepley                                                                                        void (**)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]),
200cb1e1211SMatthew G Knepley                                                                                        void (**)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), PetscScalar[]),
201cb1e1211SMatthew G Knepley                                           PetscErrorCode (*integrateJacobianFEM)(PetscInt, PetscInt, PetscInt, PetscInt, PetscQuadrature[], const PetscScalar[],
202cb1e1211SMatthew G Knepley                                                                                  const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[],
203cb1e1211SMatthew G Knepley                                                                                  void (*)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]),
204cb1e1211SMatthew G Knepley                                                                                  void (*)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]),
205cb1e1211SMatthew G Knepley                                                                                  void (*)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]),
206cb1e1211SMatthew G Knepley                                                                                  void (*)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), PetscScalar[]))
207cb1e1211SMatthew G Knepley {
208cb1e1211SMatthew G Knepley   DM_Plex *mesh = (DM_Plex*) dm->data;
209cb1e1211SMatthew G Knepley 
210cb1e1211SMatthew G Knepley   PetscFunctionBegin;
211cb1e1211SMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
212cb1e1211SMatthew G Knepley   mesh->integrateResidualFEM       = integrateResidualFEM;
21302a80efeSMatthew G. Knepley   mesh->integrateBdResidualFEM     = integrateBdResidualFEM;
214cb1e1211SMatthew G Knepley   mesh->integrateJacobianActionFEM = integrateJacobianActionFEM;
215cb1e1211SMatthew G Knepley   mesh->integrateJacobianFEM       = integrateJacobianFEM;
216cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
217cb1e1211SMatthew G Knepley }
218cb1e1211SMatthew G Knepley 
219cb1e1211SMatthew G Knepley #undef __FUNCT__
220cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexProjectFunctionLocal"
221*a1d24da5SMatthew G. Knepley PetscErrorCode DMPlexProjectFunctionLocal(DM dm, PetscInt numComp, void (**funcs)(const PetscReal [], PetscScalar *), InsertMode mode, Vec localX)
222cb1e1211SMatthew G Knepley {
223cb1e1211SMatthew G Knepley   Vec            coordinates;
224cb1e1211SMatthew G Knepley   PetscSection   section, cSection;
225cb1e1211SMatthew G Knepley   PetscInt       dim, vStart, vEnd, v, c, d;
226cb1e1211SMatthew G Knepley   PetscScalar   *values, *cArray;
227cb1e1211SMatthew G Knepley   PetscReal     *coords;
228cb1e1211SMatthew G Knepley   PetscErrorCode ierr;
229cb1e1211SMatthew G Knepley 
230cb1e1211SMatthew G Knepley   PetscFunctionBegin;
231cb1e1211SMatthew G Knepley   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
232cb1e1211SMatthew G Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
233cb1e1211SMatthew G Knepley   ierr = DMPlexGetCoordinateSection(dm, &cSection);CHKERRQ(ierr);
234cb1e1211SMatthew G Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
235cb1e1211SMatthew G Knepley   ierr = PetscMalloc(numComp * sizeof(PetscScalar), &values);CHKERRQ(ierr);
236cb1e1211SMatthew G Knepley   ierr = VecGetArray(coordinates, &cArray);CHKERRQ(ierr);
237cb1e1211SMatthew G Knepley   ierr = PetscSectionGetDof(cSection, vStart, &dim);CHKERRQ(ierr);
238cb1e1211SMatthew G Knepley   ierr = PetscMalloc(dim * sizeof(PetscReal),&coords);CHKERRQ(ierr);
239cb1e1211SMatthew G Knepley   for (v = vStart; v < vEnd; ++v) {
240cb1e1211SMatthew G Knepley     PetscInt dof, off;
241cb1e1211SMatthew G Knepley 
242cb1e1211SMatthew G Knepley     ierr = PetscSectionGetDof(cSection, v, &dof);CHKERRQ(ierr);
243cb1e1211SMatthew G Knepley     ierr = PetscSectionGetOffset(cSection, v, &off);CHKERRQ(ierr);
244cb1e1211SMatthew G Knepley     if (dof > dim) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Cannot have more coordinates %d then dimensions %d", dof, dim);
245cb1e1211SMatthew G Knepley     for (d = 0; d < dof; ++d) coords[d] = PetscRealPart(cArray[off+d]);
246*a1d24da5SMatthew G. Knepley     for (c = 0; c < numComp; ++c) (*funcs[c])(coords, &values[c]);
247cb1e1211SMatthew G Knepley     ierr = VecSetValuesSection(localX, section, v, values, mode);CHKERRQ(ierr);
248cb1e1211SMatthew G Knepley   }
249cb1e1211SMatthew G Knepley   ierr = VecRestoreArray(coordinates, &cArray);CHKERRQ(ierr);
250cb1e1211SMatthew G Knepley   /* Temporary, must be replaced by a projection on the finite element basis */
251cb1e1211SMatthew G Knepley   {
252cb1e1211SMatthew G Knepley     PetscInt eStart = 0, eEnd = 0, e, depth;
253cb1e1211SMatthew G Knepley 
254cb1e1211SMatthew G Knepley     ierr = DMPlexGetLabelSize(dm, "depth", &depth);CHKERRQ(ierr);
255cb1e1211SMatthew G Knepley     --depth;
256cb1e1211SMatthew G Knepley     if (depth > 1) {ierr = DMPlexGetDepthStratum(dm, 1, &eStart, &eEnd);CHKERRQ(ierr);}
257cb1e1211SMatthew G Knepley     for (e = eStart; e < eEnd; ++e) {
258cb1e1211SMatthew G Knepley       const PetscInt *cone = NULL;
259cb1e1211SMatthew G Knepley       PetscInt        coneSize, d;
260cb1e1211SMatthew G Knepley       PetscScalar    *coordsA, *coordsB;
261cb1e1211SMatthew G Knepley 
262cb1e1211SMatthew G Knepley       ierr = DMPlexGetConeSize(dm, e, &coneSize);CHKERRQ(ierr);
263cb1e1211SMatthew G Knepley       ierr = DMPlexGetCone(dm, e, &cone);CHKERRQ(ierr);
264cb1e1211SMatthew G Knepley       if (coneSize != 2) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_SIZ, "Cone size %d for point %d should be 2", coneSize, e);
265cb1e1211SMatthew G Knepley       ierr = VecGetValuesSection(coordinates, cSection, cone[0], &coordsA);CHKERRQ(ierr);
266cb1e1211SMatthew G Knepley       ierr = VecGetValuesSection(coordinates, cSection, cone[1], &coordsB);CHKERRQ(ierr);
267cb1e1211SMatthew G Knepley       for (d = 0; d < dim; ++d) {
268cb1e1211SMatthew G Knepley         coords[d] = 0.5*(PetscRealPart(coordsA[d]) + PetscRealPart(coordsB[d]));
269cb1e1211SMatthew G Knepley       }
270*a1d24da5SMatthew G. Knepley       for (c = 0; c < numComp; ++c) (*funcs[c])(coords, &values[c]);
271cb1e1211SMatthew G Knepley       ierr = VecSetValuesSection(localX, section, e, values, mode);CHKERRQ(ierr);
272cb1e1211SMatthew G Knepley     }
273cb1e1211SMatthew G Knepley   }
274cb1e1211SMatthew G Knepley 
275cb1e1211SMatthew G Knepley   ierr = PetscFree(coords);CHKERRQ(ierr);
276cb1e1211SMatthew G Knepley   ierr = PetscFree(values);CHKERRQ(ierr);
277cb1e1211SMatthew G Knepley #if 0
278cb1e1211SMatthew G Knepley   const PetscInt localDof = this->_mesh->sizeWithBC(s, *cells->begin());
279cb1e1211SMatthew G Knepley   PetscReal      detJ;
280cb1e1211SMatthew G Knepley 
281cb1e1211SMatthew G Knepley   ierr = PetscMalloc(localDof * sizeof(PetscScalar), &values);CHKERRQ(ierr);
282cb1e1211SMatthew G Knepley   ierr = PetscMalloc2(dim,PetscReal,&v0,dim*dim,PetscReal,&J);CHKERRQ(ierr);
283cb1e1211SMatthew G Knepley   ALE::ISieveVisitor::PointRetriever<PETSC_MESH_TYPE::sieve_type> pV(PetscPowInt(this->_mesh->getSieve()->getMaxConeSize(),dim+1), true);
284cb1e1211SMatthew G Knepley 
285cb1e1211SMatthew G Knepley   for (PetscInt c = cStart; c < cEnd; ++c) {
286cb1e1211SMatthew G Knepley     ALE::ISieveTraversal<PETSC_MESH_TYPE::sieve_type>::orientedClosure(*this->_mesh->getSieve(), c, pV);
287cb1e1211SMatthew G Knepley     const PETSC_MESH_TYPE::point_type *oPoints = pV.getPoints();
288cb1e1211SMatthew G Knepley     const int                          oSize   = pV.getSize();
289cb1e1211SMatthew G Knepley     int                                v       = 0;
290cb1e1211SMatthew G Knepley 
291cb1e1211SMatthew G Knepley     ierr = DMPlexComputeCellGeometry(dm, c, v0, J, NULL, &detJ);CHKERRQ(ierr);
292cb1e1211SMatthew G Knepley     for (PetscInt cl = 0; cl < oSize; ++cl) {
293cb1e1211SMatthew G Knepley       const PetscInt fDim;
294cb1e1211SMatthew G Knepley 
295cb1e1211SMatthew G Knepley       ierr = PetscSectionGetDof(oPoints[cl], &fDim);CHKERRQ(ierr);
296cb1e1211SMatthew G Knepley       if (pointDim) {
297cb1e1211SMatthew G Knepley         for (PetscInt d = 0; d < fDim; ++d, ++v) {
298cb1e1211SMatthew G Knepley           values[v] = (*this->_options.integrate)(v0, J, v, initFunc);
299cb1e1211SMatthew G Knepley         }
300cb1e1211SMatthew G Knepley       }
301cb1e1211SMatthew G Knepley     }
302cb1e1211SMatthew G Knepley     ierr = DMPlexVecSetClosure(dm, NULL, localX, c, values);CHKERRQ(ierr);
303cb1e1211SMatthew G Knepley     pV.clear();
304cb1e1211SMatthew G Knepley   }
305cb1e1211SMatthew G Knepley   ierr = PetscFree2(v0,J);CHKERRQ(ierr);
306cb1e1211SMatthew G Knepley   ierr = PetscFree(values);CHKERRQ(ierr);
307cb1e1211SMatthew G Knepley #endif
308cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
309cb1e1211SMatthew G Knepley }
310cb1e1211SMatthew G Knepley 
311cb1e1211SMatthew G Knepley #undef __FUNCT__
312cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexProjectFunction"
313cb1e1211SMatthew G Knepley /*@C
314cb1e1211SMatthew G Knepley   DMPlexProjectFunction - This projects the given function into the function space provided.
315cb1e1211SMatthew G Knepley 
316cb1e1211SMatthew G Knepley   Input Parameters:
317cb1e1211SMatthew G Knepley + dm      - The DM
318cb1e1211SMatthew G Knepley . numComp - The number of components (functions)
319cb1e1211SMatthew G Knepley . funcs   - The coordinate functions to evaluate
320cb1e1211SMatthew G Knepley - mode    - The insertion mode for values
321cb1e1211SMatthew G Knepley 
322cb1e1211SMatthew G Knepley   Output Parameter:
323cb1e1211SMatthew G Knepley . X - vector
324cb1e1211SMatthew G Knepley 
325cb1e1211SMatthew G Knepley   Level: developer
326cb1e1211SMatthew G Knepley 
327cb1e1211SMatthew G Knepley   Note:
328cb1e1211SMatthew G Knepley   This currently just calls the function with the coordinates of each vertex and edge midpoint, and stores the result in a vector.
329cb1e1211SMatthew G Knepley   We will eventually fix it.
330cb1e1211SMatthew G Knepley 
331878cb397SSatish Balay .seealso: DMPlexComputeL2Diff()
332878cb397SSatish Balay @*/
333*a1d24da5SMatthew G. Knepley PetscErrorCode DMPlexProjectFunction(DM dm, PetscInt numComp, void (**funcs)(const PetscReal [], PetscScalar *), InsertMode mode, Vec X)
334cb1e1211SMatthew G Knepley {
335cb1e1211SMatthew G Knepley   Vec            localX;
336cb1e1211SMatthew G Knepley   PetscErrorCode ierr;
337cb1e1211SMatthew G Knepley 
338cb1e1211SMatthew G Knepley   PetscFunctionBegin;
339cb1e1211SMatthew G Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
340cb1e1211SMatthew G Knepley   ierr = DMPlexProjectFunctionLocal(dm, numComp, funcs, mode, localX);CHKERRQ(ierr);
341cb1e1211SMatthew G Knepley   ierr = DMLocalToGlobalBegin(dm, localX, mode, X);CHKERRQ(ierr);
342cb1e1211SMatthew G Knepley   ierr = DMLocalToGlobalEnd(dm, localX, mode, X);CHKERRQ(ierr);
343cb1e1211SMatthew G Knepley   ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
344cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
345cb1e1211SMatthew G Knepley }
346cb1e1211SMatthew G Knepley 
347cb1e1211SMatthew G Knepley #undef __FUNCT__
348cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexComputeL2Diff"
349cb1e1211SMatthew G Knepley /*@C
350cb1e1211SMatthew G Knepley   DMPlexComputeL2Diff - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h.
351cb1e1211SMatthew G Knepley 
352cb1e1211SMatthew G Knepley   Input Parameters:
353cb1e1211SMatthew G Knepley + dm    - The DM
354cb1e1211SMatthew G Knepley . quad  - The PetscQuadrature object for each field
355cb1e1211SMatthew G Knepley . funcs - The functions to evaluate for each field component
356cb1e1211SMatthew G Knepley - X     - The coefficient vector u_h
357cb1e1211SMatthew G Knepley 
358cb1e1211SMatthew G Knepley   Output Parameter:
359cb1e1211SMatthew G Knepley . diff - The diff ||u - u_h||_2
360cb1e1211SMatthew G Knepley 
361cb1e1211SMatthew G Knepley   Level: developer
362cb1e1211SMatthew G Knepley 
363cb1e1211SMatthew G Knepley .seealso: DMPlexProjectFunction()
364878cb397SSatish Balay @*/
365*a1d24da5SMatthew G. Knepley PetscErrorCode DMPlexComputeL2Diff(DM dm, PetscQuadrature quad[], void (**funcs)(const PetscReal [], PetscScalar *), Vec X, PetscReal *diff)
366cb1e1211SMatthew G Knepley {
367cb1e1211SMatthew G Knepley   const PetscInt debug = 0;
368cb1e1211SMatthew G Knepley   PetscSection   section;
369cb1e1211SMatthew G Knepley   Vec            localX;
370cb1e1211SMatthew G Knepley   PetscReal     *coords, *v0, *J, *invJ, detJ;
371cb1e1211SMatthew G Knepley   PetscReal      localDiff = 0.0;
372cb1e1211SMatthew G Knepley   PetscInt       dim, numFields, numComponents = 0, cStart, cEnd, c, field, fieldOffset, comp;
373cb1e1211SMatthew G Knepley   PetscErrorCode ierr;
374cb1e1211SMatthew G Knepley 
375cb1e1211SMatthew G Knepley   PetscFunctionBegin;
376cb1e1211SMatthew G Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
377cb1e1211SMatthew G Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
378cb1e1211SMatthew G Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
379cb1e1211SMatthew G Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
380cb1e1211SMatthew G Knepley   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
381cb1e1211SMatthew G Knepley   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
382cb1e1211SMatthew G Knepley   for (field = 0; field < numFields; ++field) {
383cb1e1211SMatthew G Knepley     numComponents += quad[field].numComponents;
384cb1e1211SMatthew G Knepley   }
385cb1e1211SMatthew G Knepley   ierr = DMPlexProjectFunctionLocal(dm, numComponents, funcs, INSERT_BC_VALUES, localX);CHKERRQ(ierr);
386cb1e1211SMatthew G Knepley   ierr = PetscMalloc4(dim,PetscReal,&coords,dim,PetscReal,&v0,dim*dim,PetscReal,&J,dim*dim,PetscReal,&invJ);CHKERRQ(ierr);
387cb1e1211SMatthew G Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
388cb1e1211SMatthew G Knepley   for (c = cStart; c < cEnd; ++c) {
389cb1e1211SMatthew G Knepley     PetscScalar *x;
390cb1e1211SMatthew G Knepley     PetscReal    elemDiff = 0.0;
391cb1e1211SMatthew G Knepley 
392cb1e1211SMatthew G Knepley     ierr = DMPlexComputeCellGeometry(dm, c, v0, J, invJ, &detJ);CHKERRQ(ierr);
393cb1e1211SMatthew G Knepley     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
394cb1e1211SMatthew G Knepley     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
395cb1e1211SMatthew G Knepley 
396cb1e1211SMatthew G Knepley     for (field = 0, comp = 0, fieldOffset = 0; field < numFields; ++field) {
397cb1e1211SMatthew G Knepley       const PetscInt   numQuadPoints = quad[field].numQuadPoints;
398cb1e1211SMatthew G Knepley       const PetscReal *quadPoints    = quad[field].quadPoints;
399cb1e1211SMatthew G Knepley       const PetscReal *quadWeights   = quad[field].quadWeights;
400cb1e1211SMatthew G Knepley       const PetscInt   numBasisFuncs = quad[field].numBasisFuncs;
401cb1e1211SMatthew G Knepley       const PetscInt   numBasisComps = quad[field].numComponents;
402cb1e1211SMatthew G Knepley       const PetscReal *basis         = quad[field].basis;
403cb1e1211SMatthew G Knepley       PetscInt         q, d, e, fc, f;
404cb1e1211SMatthew G Knepley 
405cb1e1211SMatthew G Knepley       if (debug) {
406cb1e1211SMatthew G Knepley         char title[1024];
407cb1e1211SMatthew G Knepley         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
408cb1e1211SMatthew G Knepley         ierr = DMPrintCellVector(c, title, numBasisFuncs*numBasisComps, &x[fieldOffset]);CHKERRQ(ierr);
409cb1e1211SMatthew G Knepley       }
410cb1e1211SMatthew G Knepley       for (q = 0; q < numQuadPoints; ++q) {
411cb1e1211SMatthew G Knepley         for (d = 0; d < dim; d++) {
412cb1e1211SMatthew G Knepley           coords[d] = v0[d];
413cb1e1211SMatthew G Knepley           for (e = 0; e < dim; e++) {
414cb1e1211SMatthew G Knepley             coords[d] += J[d*dim+e]*(quadPoints[q*dim+e] + 1.0);
415cb1e1211SMatthew G Knepley           }
416cb1e1211SMatthew G Knepley         }
417cb1e1211SMatthew G Knepley         for (fc = 0; fc < numBasisComps; ++fc) {
418*a1d24da5SMatthew G. Knepley           PetscScalar funcVal;
419*a1d24da5SMatthew G. Knepley           PetscScalar interpolant = 0.0;
420*a1d24da5SMatthew G. Knepley 
421*a1d24da5SMatthew G. Knepley           (*funcs[comp+fc])(coords, &funcVal);
422cb1e1211SMatthew G Knepley           for (f = 0; f < numBasisFuncs; ++f) {
423cb1e1211SMatthew G Knepley             const PetscInt fidx = f*numBasisComps+fc;
424*a1d24da5SMatthew G. Knepley             interpolant += x[fieldOffset+fidx]*basis[q*numBasisFuncs*numBasisComps+fidx];
425cb1e1211SMatthew G Knepley           }
426*a1d24da5SMatthew G. Knepley           if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "    elem %d field %d diff %g\n", c, field, PetscSqr(PetscRealPart(interpolant - funcVal))*quadWeights[q]*detJ);CHKERRQ(ierr);}
427*a1d24da5SMatthew G. Knepley           elemDiff += PetscSqr(PetscRealPart(interpolant - funcVal))*quadWeights[q]*detJ;
428cb1e1211SMatthew G Knepley         }
429cb1e1211SMatthew G Knepley       }
430cb1e1211SMatthew G Knepley       comp        += numBasisComps;
431cb1e1211SMatthew G Knepley       fieldOffset += numBasisFuncs*numBasisComps;
432cb1e1211SMatthew G Knepley     }
433cb1e1211SMatthew G Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
434cb1e1211SMatthew G Knepley     if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);}
435cb1e1211SMatthew G Knepley     localDiff += elemDiff;
436cb1e1211SMatthew G Knepley   }
437cb1e1211SMatthew G Knepley   ierr  = PetscFree4(coords,v0,J,invJ);CHKERRQ(ierr);
438cb1e1211SMatthew G Knepley   ierr  = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
439cb1e1211SMatthew G Knepley   ierr  = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPI_SUM, PETSC_COMM_WORLD);CHKERRQ(ierr);
440cb1e1211SMatthew G Knepley   *diff = PetscSqrtReal(*diff);
441cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
442cb1e1211SMatthew G Knepley }
443cb1e1211SMatthew G Knepley 
444cb1e1211SMatthew G Knepley #undef __FUNCT__
445cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexComputeResidualFEM"
446cb1e1211SMatthew G Knepley /*@
447cb1e1211SMatthew G Knepley   DMPlexComputeResidualFEM - Form the local residual F from the local input X using pointwise functions specified by the user
448cb1e1211SMatthew G Knepley 
449cb1e1211SMatthew G Knepley   Input Parameters:
450cb1e1211SMatthew G Knepley + dm - The mesh
451cb1e1211SMatthew G Knepley . X  - Local input vector
452cb1e1211SMatthew G Knepley - user - The user context
453cb1e1211SMatthew G Knepley 
454cb1e1211SMatthew G Knepley   Output Parameter:
455cb1e1211SMatthew G Knepley . F  - Local output vector
456cb1e1211SMatthew G Knepley 
457cb1e1211SMatthew G Knepley   Note:
458cb1e1211SMatthew G Knepley   The second member of the user context must be an FEMContext.
459cb1e1211SMatthew G Knepley 
460cb1e1211SMatthew G Knepley   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
461cb1e1211SMatthew G Knepley   like a GPU, or vectorize on a multicore machine.
462cb1e1211SMatthew G Knepley 
4630059ad2aSSatish Balay   Level: developer
4640059ad2aSSatish Balay 
465cb1e1211SMatthew G Knepley .seealso: DMPlexComputeJacobianActionFEM()
466878cb397SSatish Balay @*/
467cb1e1211SMatthew G Knepley PetscErrorCode DMPlexComputeResidualFEM(DM dm, Vec X, Vec F, void *user)
468cb1e1211SMatthew G Knepley {
469cb1e1211SMatthew G Knepley   DM_Plex         *mesh   = (DM_Plex*) dm->data;
470cb1e1211SMatthew G Knepley   PetscFEM        *fem    = (PetscFEM*) &((DM*) user)[1];
471cb1e1211SMatthew G Knepley   PetscQuadrature *quad   = fem->quad;
472652b88e8SMatthew G. Knepley   PetscQuadrature *quadBd = fem->quadBd;
473cb1e1211SMatthew G Knepley   PetscSection     section;
47402a80efeSMatthew G. Knepley   PetscReal       *v0, *n, *J, *invJ, *detJ;
475cb1e1211SMatthew G Knepley   PetscScalar     *elemVec, *u;
476cb1e1211SMatthew G Knepley   PetscInt         dim, numFields, field, numBatchesTmp = 1, numCells, cStart, cEnd, c;
477652b88e8SMatthew G. Knepley   PetscInt         cellDof, numComponents;
478652b88e8SMatthew G. Knepley   PetscBool        has;
479cb1e1211SMatthew G Knepley   PetscErrorCode   ierr;
480cb1e1211SMatthew G Knepley 
481cb1e1211SMatthew G Knepley   PetscFunctionBegin;
482cb1e1211SMatthew G Knepley   /* ierr = PetscLogEventBegin(ResidualFEMEvent,0,0,0,0);CHKERRQ(ierr); */
483cb1e1211SMatthew G Knepley   ierr     = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
484cb1e1211SMatthew G Knepley   ierr     = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
485cb1e1211SMatthew G Knepley   ierr     = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
486cb1e1211SMatthew G Knepley   ierr     = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
487cb1e1211SMatthew G Knepley   numCells = cEnd - cStart;
488652b88e8SMatthew G. Knepley   for (field = 0, cellDof = 0, numComponents = 0; field < numFields; ++field) {
489cb1e1211SMatthew G Knepley     cellDof       += quad[field].numBasisFuncs*quad[field].numComponents;
490cb1e1211SMatthew G Knepley     numComponents += quad[field].numComponents;
491cb1e1211SMatthew G Knepley   }
492cb1e1211SMatthew G Knepley   ierr = DMPlexProjectFunctionLocal(dm, numComponents, fem->bcFuncs, INSERT_BC_VALUES, X);CHKERRQ(ierr);
493cb1e1211SMatthew G Knepley   ierr = VecSet(F, 0.0);CHKERRQ(ierr);
494cb1e1211SMatthew 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);
495cb1e1211SMatthew G Knepley   for (c = cStart; c < cEnd; ++c) {
496cb1e1211SMatthew G Knepley     PetscScalar *x;
497cb1e1211SMatthew G Knepley     PetscInt     i;
498cb1e1211SMatthew G Knepley 
499cb1e1211SMatthew G Knepley     ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr);
500cb1e1211SMatthew G Knepley     if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c);
501cb1e1211SMatthew G Knepley     ierr = DMPlexVecGetClosure(dm, NULL, X, c, NULL, &x);CHKERRQ(ierr);
502cb1e1211SMatthew G Knepley 
503cb1e1211SMatthew G Knepley     for (i = 0; i < cellDof; ++i) u[c*cellDof+i] = x[i];
504cb1e1211SMatthew G Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, X, c, NULL, &x);CHKERRQ(ierr);
505cb1e1211SMatthew G Knepley   }
506cb1e1211SMatthew G Knepley   for (field = 0; field < numFields; ++field) {
507cb1e1211SMatthew G Knepley     const PetscInt numQuadPoints = quad[field].numQuadPoints;
508cb1e1211SMatthew G Knepley     const PetscInt numBasisFuncs = quad[field].numBasisFuncs;
509cb1e1211SMatthew G Knepley     void           (*f0)(const PetscScalar u[], const PetscScalar gradU[], const PetscReal x[], PetscScalar f0[]) = fem->f0Funcs[field];
510cb1e1211SMatthew G Knepley     void           (*f1)(const PetscScalar u[], const PetscScalar gradU[], const PetscReal x[], PetscScalar f1[]) = fem->f1Funcs[field];
511cb1e1211SMatthew G Knepley     /* Conforming batches */
512cb1e1211SMatthew G Knepley     PetscInt blockSize  = numBasisFuncs*numQuadPoints;
513cb1e1211SMatthew G Knepley     PetscInt numBlocks  = 1;
514cb1e1211SMatthew G Knepley     PetscInt batchSize  = numBlocks * blockSize;
515cb1e1211SMatthew G Knepley     PetscInt numBatches = numBatchesTmp;
516cb1e1211SMatthew G Knepley     PetscInt numChunks  = numCells / (numBatches*batchSize);
517cb1e1211SMatthew G Knepley     /* Remainder */
518cb1e1211SMatthew G Knepley     PetscInt numRemainder = numCells % (numBatches * batchSize);
519cb1e1211SMatthew G Knepley     PetscInt offset       = numCells - numRemainder;
520cb1e1211SMatthew G Knepley 
521cb1e1211SMatthew G Knepley     ierr = (*mesh->integrateResidualFEM)(numChunks*numBatches*batchSize, numFields, field, quad, u, v0, J, invJ, detJ, f0, f1, elemVec);CHKERRQ(ierr);
522cb1e1211SMatthew G Knepley     ierr = (*mesh->integrateResidualFEM)(numRemainder, numFields, field, quad, &u[offset*cellDof], &v0[offset*dim], &J[offset*dim*dim], &invJ[offset*dim*dim], &detJ[offset],
523cb1e1211SMatthew G Knepley                                          f0, f1, &elemVec[offset*cellDof]);CHKERRQ(ierr);
524cb1e1211SMatthew G Knepley   }
525cb1e1211SMatthew G Knepley   for (c = cStart; c < cEnd; ++c) {
526cb1e1211SMatthew G Knepley     if (mesh->printFEM > 1) {ierr = DMPrintCellVector(c, "Residual", cellDof, &elemVec[c*cellDof]);CHKERRQ(ierr);}
527cb1e1211SMatthew G Knepley     ierr = DMPlexVecSetClosure(dm, NULL, F, c, &elemVec[c*cellDof], ADD_VALUES);CHKERRQ(ierr);
528cb1e1211SMatthew G Knepley   }
529cb1e1211SMatthew G Knepley   ierr = PetscFree6(u,v0,J,invJ,detJ,elemVec);CHKERRQ(ierr);
530652b88e8SMatthew G. Knepley   /* Integration over the boundary:
531652b88e8SMatthew G. Knepley      - This can probably be generalized to integration over a set of labels, however
532652b88e8SMatthew G. Knepley        the idea here is to do integration where we need the cell normal
533652b88e8SMatthew G. Knepley      - We can replace hardcoding with a registration process, and this is how we hook
534652b88e8SMatthew G. Knepley        up the system to something like FEniCS
535652b88e8SMatthew G. Knepley   */
536652b88e8SMatthew G. Knepley   ierr = DMPlexHasLabel(dm, "boundary", &has);CHKERRQ(ierr);
537652b88e8SMatthew G. Knepley   if (has && quadBd) {
538652b88e8SMatthew G. Knepley     DMLabel         label;
539652b88e8SMatthew G. Knepley     IS              pointIS;
540652b88e8SMatthew G. Knepley     const PetscInt *points;
541652b88e8SMatthew G. Knepley     PetscInt        numPoints, p;
542652b88e8SMatthew G. Knepley 
543652b88e8SMatthew G. Knepley     ierr = DMPlexGetLabel(dm, "boundary", &label);CHKERRQ(ierr);
544652b88e8SMatthew G. Knepley     ierr = DMLabelGetStratumSize(label, 1, &numPoints);CHKERRQ(ierr);
545652b88e8SMatthew G. Knepley     ierr = DMLabelGetStratumIS(label, 1, &pointIS);CHKERRQ(ierr);
546652b88e8SMatthew G. Knepley     ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
547652b88e8SMatthew G. Knepley     for (field = 0, cellDof = 0, numComponents = 0; field < numFields; ++field) {
548652b88e8SMatthew G. Knepley       cellDof       += quadBd[field].numBasisFuncs*quadBd[field].numComponents;
549652b88e8SMatthew G. Knepley       numComponents += quadBd[field].numComponents;
550652b88e8SMatthew G. Knepley     }
55102a80efeSMatthew 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);
552652b88e8SMatthew G. Knepley     for (p = 0; p < numPoints; ++p) {
553652b88e8SMatthew G. Knepley       const PetscInt point = points[p];
554652b88e8SMatthew G. Knepley       PetscScalar   *x;
555652b88e8SMatthew G. Knepley       PetscInt       i;
556652b88e8SMatthew G. Knepley 
55702a80efeSMatthew G. Knepley       /* TODO: Add normal determination here */
558652b88e8SMatthew G. Knepley       ierr = DMPlexComputeCellGeometry(dm, point, &v0[p*dim], &J[p*dim*dim], &invJ[p*dim*dim], &detJ[p]);CHKERRQ(ierr);
559652b88e8SMatthew G. Knepley       if (detJ[p] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[p], p);
560652b88e8SMatthew G. Knepley       ierr = DMPlexVecGetClosure(dm, NULL, X, point, NULL, &x);CHKERRQ(ierr);
561652b88e8SMatthew G. Knepley 
562652b88e8SMatthew G. Knepley       for (i = 0; i < cellDof; ++i) u[p*cellDof+i] = x[i];
563652b88e8SMatthew G. Knepley       ierr = DMPlexVecRestoreClosure(dm, NULL, X, point, NULL, &x);CHKERRQ(ierr);
564652b88e8SMatthew G. Knepley     }
565652b88e8SMatthew G. Knepley     for (field = 0; field < numFields; ++field) {
566652b88e8SMatthew G. Knepley       const PetscInt numQuadPoints = quadBd[field].numQuadPoints;
567652b88e8SMatthew G. Knepley       const PetscInt numBasisFuncs = quadBd[field].numBasisFuncs;
56802a80efeSMatthew G. Knepley       void           (*f0)(const PetscScalar u[], const PetscScalar gradU[], const PetscReal x[], const PetscReal n[], PetscScalar f0[]) = fem->f0BdFuncs[field];
56902a80efeSMatthew G. Knepley       void           (*f1)(const PetscScalar u[], const PetscScalar gradU[], const PetscReal x[], const PetscReal n[], PetscScalar f1[]) = fem->f1BdFuncs[field];
570652b88e8SMatthew G. Knepley       /* Conforming batches */
571652b88e8SMatthew G. Knepley       PetscInt blockSize  = numBasisFuncs*numQuadPoints;
572652b88e8SMatthew G. Knepley       PetscInt numBlocks  = 1;
573652b88e8SMatthew G. Knepley       PetscInt batchSize  = numBlocks * blockSize;
574652b88e8SMatthew G. Knepley       PetscInt numBatches = numBatchesTmp;
575652b88e8SMatthew G. Knepley       PetscInt numChunks  = numPoints / (numBatches*batchSize);
576652b88e8SMatthew G. Knepley       /* Remainder */
577652b88e8SMatthew G. Knepley       PetscInt numRemainder = numPoints % (numBatches * batchSize);
578652b88e8SMatthew G. Knepley       PetscInt offset       = numPoints - numRemainder;
579652b88e8SMatthew G. Knepley 
58002a80efeSMatthew G. Knepley       ierr = (*mesh->integrateBdResidualFEM)(numChunks*numBatches*batchSize, numFields, field, quadBd, u, v0, n, J, invJ, detJ, f0, f1, elemVec);CHKERRQ(ierr);
58102a80efeSMatthew 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],
582652b88e8SMatthew G. Knepley                                              f0, f1, &elemVec[offset*cellDof]);CHKERRQ(ierr);
583652b88e8SMatthew G. Knepley     }
584652b88e8SMatthew G. Knepley     for (p = 0; p < numPoints; ++p) {
585652b88e8SMatthew G. Knepley       const PetscInt point = points[p];
586652b88e8SMatthew G. Knepley 
587652b88e8SMatthew G. Knepley       if (mesh->printFEM > 1) {ierr = DMPrintCellVector(point, "Residual", cellDof, &elemVec[p*cellDof]);CHKERRQ(ierr);}
588652b88e8SMatthew G. Knepley       ierr = DMPlexVecSetClosure(dm, NULL, F, point, &elemVec[p*cellDof], ADD_VALUES);CHKERRQ(ierr);
589652b88e8SMatthew G. Knepley     }
590652b88e8SMatthew G. Knepley     ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
591652b88e8SMatthew G. Knepley     ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
59202a80efeSMatthew G. Knepley     ierr = PetscFree7(u,v0,n,J,invJ,detJ,elemVec);CHKERRQ(ierr);
593652b88e8SMatthew G. Knepley   }
594cb1e1211SMatthew G Knepley   if (mesh->printFEM) {
595cb1e1211SMatthew G Knepley     PetscMPIInt rank, numProcs;
596cb1e1211SMatthew G Knepley     PetscInt    p;
597cb1e1211SMatthew G Knepley 
598cb1e1211SMatthew G Knepley     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr);
599cb1e1211SMatthew G Knepley     ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm), &numProcs);CHKERRQ(ierr);
600cb1e1211SMatthew G Knepley     ierr = PetscPrintf(PETSC_COMM_WORLD, "Residual:\n");CHKERRQ(ierr);
601cb1e1211SMatthew G Knepley     for (p = 0; p < numProcs; ++p) {
602cb1e1211SMatthew G Knepley       if (p == rank) {
603cb1e1211SMatthew G Knepley         Vec f;
604cb1e1211SMatthew G Knepley 
605cb1e1211SMatthew G Knepley         ierr = VecDuplicate(F, &f);CHKERRQ(ierr);
606cb1e1211SMatthew G Knepley         ierr = VecCopy(F, f);CHKERRQ(ierr);
607cb1e1211SMatthew G Knepley         ierr = VecChop(f, 1.0e-10);CHKERRQ(ierr);
608cb1e1211SMatthew G Knepley         ierr = VecView(f, PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
609cb1e1211SMatthew G Knepley         ierr = VecDestroy(&f);CHKERRQ(ierr);
610cb1e1211SMatthew G Knepley         ierr = PetscViewerFlush(PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
611cb1e1211SMatthew G Knepley       }
612cb1e1211SMatthew G Knepley       ierr = PetscBarrier((PetscObject) dm);CHKERRQ(ierr);
613cb1e1211SMatthew G Knepley     }
614cb1e1211SMatthew G Knepley   }
615cb1e1211SMatthew G Knepley   /* ierr = PetscLogEventEnd(ResidualFEMEvent,0,0,0,0);CHKERRQ(ierr); */
616cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
617cb1e1211SMatthew G Knepley }
618cb1e1211SMatthew G Knepley 
619cb1e1211SMatthew G Knepley #undef __FUNCT__
620cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexComputeJacobianActionFEM"
621cb1e1211SMatthew G Knepley /*@C
622cb1e1211SMatthew G Knepley   DMPlexComputeJacobianActionFEM - Form the local action of Jacobian J(u) on the local input X using pointwise functions specified by the user
623cb1e1211SMatthew G Knepley 
624cb1e1211SMatthew G Knepley   Input Parameters:
625cb1e1211SMatthew G Knepley + dm - The mesh
626cb1e1211SMatthew G Knepley . J  - The Jacobian shell matrix
627cb1e1211SMatthew G Knepley . X  - Local input vector
628cb1e1211SMatthew G Knepley - user - The user context
629cb1e1211SMatthew G Knepley 
630cb1e1211SMatthew G Knepley   Output Parameter:
631cb1e1211SMatthew G Knepley . F  - Local output vector
632cb1e1211SMatthew G Knepley 
633cb1e1211SMatthew G Knepley   Note:
634cb1e1211SMatthew G Knepley   The second member of the user context must be an FEMContext.
635cb1e1211SMatthew G Knepley 
636cb1e1211SMatthew G Knepley   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
637cb1e1211SMatthew G Knepley   like a GPU, or vectorize on a multicore machine.
638cb1e1211SMatthew G Knepley 
6390059ad2aSSatish Balay   Level: developer
6400059ad2aSSatish Balay 
641cb1e1211SMatthew G Knepley .seealso: DMPlexComputeResidualFEM()
642878cb397SSatish Balay @*/
643cb1e1211SMatthew G Knepley PetscErrorCode DMPlexComputeJacobianActionFEM(DM dm, Mat Jac, Vec X, Vec F, void *user)
644cb1e1211SMatthew G Knepley {
645cb1e1211SMatthew G Knepley   DM_Plex         *mesh = (DM_Plex*) dm->data;
646cb1e1211SMatthew G Knepley   PetscFEM        *fem  = (PetscFEM*) &((DM*) user)[1];
647cb1e1211SMatthew G Knepley   PetscQuadrature *quad = fem->quad;
648cb1e1211SMatthew G Knepley   PetscSection     section;
649cb1e1211SMatthew G Knepley   JacActionCtx    *jctx;
650cb1e1211SMatthew G Knepley   PetscReal       *v0, *J, *invJ, *detJ;
651cb1e1211SMatthew G Knepley   PetscScalar     *elemVec, *u, *a;
652cb1e1211SMatthew G Knepley   PetscInt         dim, numFields, field, numBatchesTmp = 1, numCells, cStart, cEnd, c;
653cb1e1211SMatthew G Knepley   PetscInt         cellDof = 0;
654cb1e1211SMatthew G Knepley   PetscErrorCode   ierr;
655cb1e1211SMatthew G Knepley 
656cb1e1211SMatthew G Knepley   PetscFunctionBegin;
657cb1e1211SMatthew G Knepley   /* ierr = PetscLogEventBegin(JacobianActionFEMEvent,0,0,0,0);CHKERRQ(ierr); */
658cb1e1211SMatthew G Knepley   ierr     = MatShellGetContext(Jac, &jctx);CHKERRQ(ierr);
659cb1e1211SMatthew G Knepley   ierr     = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
660cb1e1211SMatthew G Knepley   ierr     = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
661cb1e1211SMatthew G Knepley   ierr     = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
662cb1e1211SMatthew G Knepley   ierr     = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
663cb1e1211SMatthew G Knepley   numCells = cEnd - cStart;
664cb1e1211SMatthew G Knepley   for (field = 0; field < numFields; ++field) {
665cb1e1211SMatthew G Knepley     cellDof += quad[field].numBasisFuncs*quad[field].numComponents;
666cb1e1211SMatthew G Knepley   }
667cb1e1211SMatthew G Knepley   ierr = VecSet(F, 0.0);CHKERRQ(ierr);
668cb1e1211SMatthew 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);
669cb1e1211SMatthew G Knepley   for (c = cStart; c < cEnd; ++c) {
670cb1e1211SMatthew G Knepley     PetscScalar *x;
671cb1e1211SMatthew G Knepley     PetscInt     i;
672cb1e1211SMatthew G Knepley 
673cb1e1211SMatthew G Knepley     ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr);
674cb1e1211SMatthew G Knepley     if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c);
675cb1e1211SMatthew G Knepley     ierr = DMPlexVecGetClosure(dm, NULL, jctx->u, c, NULL, &x);CHKERRQ(ierr);
676cb1e1211SMatthew G Knepley     for (i = 0; i < cellDof; ++i) u[c*cellDof+i] = x[i];
677cb1e1211SMatthew G Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, jctx->u, c, NULL, &x);CHKERRQ(ierr);
678cb1e1211SMatthew G Knepley     ierr = DMPlexVecGetClosure(dm, NULL, X, c, NULL, &x);CHKERRQ(ierr);
679cb1e1211SMatthew G Knepley     for (i = 0; i < cellDof; ++i) a[c*cellDof+i] = x[i];
680cb1e1211SMatthew G Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, X, c, NULL, &x);CHKERRQ(ierr);
681cb1e1211SMatthew G Knepley   }
682cb1e1211SMatthew G Knepley   for (field = 0; field < numFields; ++field) {
683cb1e1211SMatthew G Knepley     const PetscInt numQuadPoints = quad[field].numQuadPoints;
684cb1e1211SMatthew G Knepley     const PetscInt numBasisFuncs = quad[field].numBasisFuncs;
685cb1e1211SMatthew G Knepley     /* Conforming batches */
686cb1e1211SMatthew G Knepley     PetscInt blockSize  = numBasisFuncs*numQuadPoints;
687cb1e1211SMatthew G Knepley     PetscInt numBlocks  = 1;
688cb1e1211SMatthew G Knepley     PetscInt batchSize  = numBlocks * blockSize;
689cb1e1211SMatthew G Knepley     PetscInt numBatches = numBatchesTmp;
690cb1e1211SMatthew G Knepley     PetscInt numChunks  = numCells / (numBatches*batchSize);
691cb1e1211SMatthew G Knepley     /* Remainder */
692cb1e1211SMatthew G Knepley     PetscInt numRemainder = numCells % (numBatches * batchSize);
693cb1e1211SMatthew G Knepley     PetscInt offset       = numCells - numRemainder;
694cb1e1211SMatthew G Knepley 
695cb1e1211SMatthew G Knepley     ierr = (*mesh->integrateJacobianActionFEM)(numChunks*numBatches*batchSize, numFields, field, quad, u, a, v0, J, invJ, detJ, fem->g0Funcs, fem->g1Funcs, fem->g2Funcs, fem->g3Funcs, elemVec);CHKERRQ(ierr);
696cb1e1211SMatthew G Knepley     ierr = (*mesh->integrateJacobianActionFEM)(numRemainder, numFields, field, quad, &u[offset*cellDof], &a[offset*cellDof], &v0[offset*dim], &J[offset*dim*dim], &invJ[offset*dim*dim], &detJ[offset],
697cb1e1211SMatthew G Knepley                                                fem->g0Funcs, fem->g1Funcs, fem->g2Funcs, fem->g3Funcs, &elemVec[offset*cellDof]);CHKERRQ(ierr);
698cb1e1211SMatthew G Knepley   }
699cb1e1211SMatthew G Knepley   for (c = cStart; c < cEnd; ++c) {
700cb1e1211SMatthew G Knepley     if (mesh->printFEM > 1) {ierr = DMPrintCellVector(c, "Jacobian Action", cellDof, &elemVec[c*cellDof]);CHKERRQ(ierr);}
701cb1e1211SMatthew G Knepley     ierr = DMPlexVecSetClosure(dm, NULL, F, c, &elemVec[c*cellDof], ADD_VALUES);CHKERRQ(ierr);
702cb1e1211SMatthew G Knepley   }
703cb1e1211SMatthew G Knepley   ierr = PetscFree7(u,a,v0,J,invJ,detJ,elemVec);CHKERRQ(ierr);
704cb1e1211SMatthew G Knepley   if (mesh->printFEM) {
705cb1e1211SMatthew G Knepley     PetscMPIInt rank, numProcs;
706cb1e1211SMatthew G Knepley     PetscInt    p;
707cb1e1211SMatthew G Knepley 
708cb1e1211SMatthew G Knepley     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr);
709cb1e1211SMatthew G Knepley     ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm), &numProcs);CHKERRQ(ierr);
710cb1e1211SMatthew G Knepley     ierr = PetscPrintf(PETSC_COMM_WORLD, "Jacobian Action:\n");CHKERRQ(ierr);
711cb1e1211SMatthew G Knepley     for (p = 0; p < numProcs; ++p) {
712cb1e1211SMatthew G Knepley       if (p == rank) {ierr = VecView(F, PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);}
713cb1e1211SMatthew G Knepley       ierr = PetscBarrier((PetscObject) dm);CHKERRQ(ierr);
714cb1e1211SMatthew G Knepley     }
715cb1e1211SMatthew G Knepley   }
716cb1e1211SMatthew G Knepley   /* ierr = PetscLogEventEnd(JacobianActionFEMEvent,0,0,0,0);CHKERRQ(ierr); */
717cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
718cb1e1211SMatthew G Knepley }
719cb1e1211SMatthew G Knepley 
720cb1e1211SMatthew G Knepley #undef __FUNCT__
721cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexComputeJacobianFEM"
722cb1e1211SMatthew G Knepley /*@
723cb1e1211SMatthew G Knepley   DMPlexComputeJacobianFEM - Form the local portion of the Jacobian matrix J at the local solution X using pointwise functions specified by the user.
724cb1e1211SMatthew G Knepley 
725cb1e1211SMatthew G Knepley   Input Parameters:
726cb1e1211SMatthew G Knepley + dm - The mesh
727cb1e1211SMatthew G Knepley . X  - Local input vector
728cb1e1211SMatthew G Knepley - user - The user context
729cb1e1211SMatthew G Knepley 
730cb1e1211SMatthew G Knepley   Output Parameter:
731cb1e1211SMatthew G Knepley . Jac  - Jacobian matrix
732cb1e1211SMatthew G Knepley 
733cb1e1211SMatthew G Knepley   Note:
734cb1e1211SMatthew G Knepley   The second member of the user context must be an FEMContext.
735cb1e1211SMatthew G Knepley 
736cb1e1211SMatthew G Knepley   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
737cb1e1211SMatthew G Knepley   like a GPU, or vectorize on a multicore machine.
738cb1e1211SMatthew G Knepley 
7390059ad2aSSatish Balay   Level: developer
7400059ad2aSSatish Balay 
741cb1e1211SMatthew G Knepley .seealso: FormFunctionLocal()
742878cb397SSatish Balay @*/
743cb1e1211SMatthew G Knepley PetscErrorCode DMPlexComputeJacobianFEM(DM dm, Vec X, Mat Jac, Mat JacP, MatStructure *str,void *user)
744cb1e1211SMatthew G Knepley {
745cb1e1211SMatthew G Knepley   DM_Plex         *mesh = (DM_Plex*) dm->data;
746cb1e1211SMatthew G Knepley   PetscFEM        *fem  = (PetscFEM*) &((DM*) user)[1];
747cb1e1211SMatthew G Knepley   PetscQuadrature *quad = fem->quad;
748cb1e1211SMatthew G Knepley   PetscSection     section;
749cb1e1211SMatthew G Knepley   PetscReal       *v0, *J, *invJ, *detJ;
750cb1e1211SMatthew G Knepley   PetscScalar     *elemMat, *u;
751cb1e1211SMatthew G Knepley   PetscInt         dim, numFields, field, fieldI, numBatchesTmp = 1, numCells, cStart, cEnd, c;
752cb1e1211SMatthew G Knepley   PetscInt         cellDof = 0, numComponents = 0;
753cb1e1211SMatthew G Knepley   PetscBool        isShell;
754cb1e1211SMatthew G Knepley   PetscErrorCode   ierr;
755cb1e1211SMatthew G Knepley 
756cb1e1211SMatthew G Knepley   PetscFunctionBegin;
757cb1e1211SMatthew G Knepley   /* ierr = PetscLogEventBegin(JacobianFEMEvent,0,0,0,0);CHKERRQ(ierr); */
758cb1e1211SMatthew G Knepley   ierr     = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
759cb1e1211SMatthew G Knepley   ierr     = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
760cb1e1211SMatthew G Knepley   ierr     = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
761cb1e1211SMatthew G Knepley   ierr     = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
762cb1e1211SMatthew G Knepley   numCells = cEnd - cStart;
763cb1e1211SMatthew G Knepley   for (field = 0; field < numFields; ++field) {
764cb1e1211SMatthew G Knepley     cellDof       += quad[field].numBasisFuncs*quad[field].numComponents;
765cb1e1211SMatthew G Knepley     numComponents += quad[field].numComponents;
766cb1e1211SMatthew G Knepley   }
767cb1e1211SMatthew G Knepley   ierr = DMPlexProjectFunctionLocal(dm, numComponents, fem->bcFuncs, INSERT_BC_VALUES, X);CHKERRQ(ierr);
768cb1e1211SMatthew G Knepley   ierr = MatZeroEntries(JacP);CHKERRQ(ierr);
769cb1e1211SMatthew 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);
770cb1e1211SMatthew G Knepley   for (c = cStart; c < cEnd; ++c) {
771cb1e1211SMatthew G Knepley     PetscScalar *x;
772cb1e1211SMatthew G Knepley     PetscInt     i;
773cb1e1211SMatthew G Knepley 
774cb1e1211SMatthew G Knepley     ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr);
775cb1e1211SMatthew G Knepley     if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c);
776cb1e1211SMatthew G Knepley     ierr = DMPlexVecGetClosure(dm, NULL, X, c, NULL, &x);CHKERRQ(ierr);
777cb1e1211SMatthew G Knepley 
778cb1e1211SMatthew G Knepley     for (i = 0; i < cellDof; ++i) u[c*cellDof+i] = x[i];
779cb1e1211SMatthew G Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, X, c, NULL, &x);CHKERRQ(ierr);
780cb1e1211SMatthew G Knepley   }
781cb1e1211SMatthew G Knepley   ierr = PetscMemzero(elemMat, numCells*cellDof*cellDof * sizeof(PetscScalar));CHKERRQ(ierr);
782cb1e1211SMatthew G Knepley   for (fieldI = 0; fieldI < numFields; ++fieldI) {
783cb1e1211SMatthew G Knepley     const PetscInt numQuadPoints = quad[fieldI].numQuadPoints;
784cb1e1211SMatthew G Knepley     const PetscInt numBasisFuncs = quad[fieldI].numBasisFuncs;
785cb1e1211SMatthew G Knepley     PetscInt       fieldJ;
786cb1e1211SMatthew G Knepley 
787cb1e1211SMatthew G Knepley     for (fieldJ = 0; fieldJ < numFields; ++fieldJ) {
788cb1e1211SMatthew G Knepley       void (*g0)(const PetscScalar u[], const PetscScalar gradU[], const PetscReal x[], PetscScalar g0[]) = fem->g0Funcs[fieldI*numFields+fieldJ];
789cb1e1211SMatthew G Knepley       void (*g1)(const PetscScalar u[], const PetscScalar gradU[], const PetscReal x[], PetscScalar g1[]) = fem->g1Funcs[fieldI*numFields+fieldJ];
790cb1e1211SMatthew G Knepley       void (*g2)(const PetscScalar u[], const PetscScalar gradU[], const PetscReal x[], PetscScalar g2[]) = fem->g2Funcs[fieldI*numFields+fieldJ];
791cb1e1211SMatthew G Knepley       void (*g3)(const PetscScalar u[], const PetscScalar gradU[], const PetscReal x[], PetscScalar g3[]) = fem->g3Funcs[fieldI*numFields+fieldJ];
792cb1e1211SMatthew G Knepley       /* Conforming batches */
793cb1e1211SMatthew G Knepley       PetscInt blockSize  = numBasisFuncs*numQuadPoints;
794cb1e1211SMatthew G Knepley       PetscInt numBlocks  = 1;
795cb1e1211SMatthew G Knepley       PetscInt batchSize  = numBlocks * blockSize;
796cb1e1211SMatthew G Knepley       PetscInt numBatches = numBatchesTmp;
797cb1e1211SMatthew G Knepley       PetscInt numChunks  = numCells / (numBatches*batchSize);
798cb1e1211SMatthew G Knepley       /* Remainder */
799cb1e1211SMatthew G Knepley       PetscInt numRemainder = numCells % (numBatches * batchSize);
800cb1e1211SMatthew G Knepley       PetscInt offset       = numCells - numRemainder;
801cb1e1211SMatthew G Knepley 
802cb1e1211SMatthew G Knepley       ierr = (*mesh->integrateJacobianFEM)(numChunks*numBatches*batchSize, numFields, fieldI, fieldJ, quad, u, v0, J, invJ, detJ, g0, g1, g2, g3, elemMat);CHKERRQ(ierr);
803cb1e1211SMatthew G Knepley       ierr = (*mesh->integrateJacobianFEM)(numRemainder, numFields, fieldI, fieldJ, quad, &u[offset*cellDof], &v0[offset*dim], &J[offset*dim*dim], &invJ[offset*dim*dim], &detJ[offset],
804cb1e1211SMatthew G Knepley                                            g0, g1, g2, g3, &elemMat[offset*cellDof*cellDof]);CHKERRQ(ierr);
805cb1e1211SMatthew G Knepley     }
806cb1e1211SMatthew G Knepley   }
807cb1e1211SMatthew G Knepley   for (c = cStart; c < cEnd; ++c) {
808cb1e1211SMatthew G Knepley     if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(c, "Jacobian", cellDof, cellDof, &elemMat[c*cellDof*cellDof]);CHKERRQ(ierr);}
809cb1e1211SMatthew G Knepley     ierr = DMPlexMatSetClosure(dm, NULL, NULL, JacP, c, &elemMat[c*cellDof*cellDof], ADD_VALUES);CHKERRQ(ierr);
810cb1e1211SMatthew G Knepley   }
811cb1e1211SMatthew G Knepley   ierr = PetscFree6(u,v0,J,invJ,detJ,elemMat);CHKERRQ(ierr);
812cb1e1211SMatthew G Knepley 
813cb1e1211SMatthew G Knepley   /* Assemble matrix, using the 2-step process:
814cb1e1211SMatthew G Knepley        MatAssemblyBegin(), MatAssemblyEnd(). */
815cb1e1211SMatthew G Knepley   ierr = MatAssemblyBegin(JacP, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
816cb1e1211SMatthew G Knepley   ierr = MatAssemblyEnd(JacP, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
817cb1e1211SMatthew G Knepley 
818cb1e1211SMatthew G Knepley   if (mesh->printFEM) {
819cb1e1211SMatthew G Knepley     ierr = PetscPrintf(PETSC_COMM_WORLD, "Jacobian:\n");CHKERRQ(ierr);
820cb1e1211SMatthew G Knepley     ierr = MatChop(JacP, 1.0e-10);CHKERRQ(ierr);
821cb1e1211SMatthew G Knepley     ierr = MatView(JacP, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
822cb1e1211SMatthew G Knepley   }
823cb1e1211SMatthew G Knepley   /* ierr = PetscLogEventEnd(JacobianFEMEvent,0,0,0,0);CHKERRQ(ierr); */
824cb1e1211SMatthew G Knepley   ierr = PetscObjectTypeCompare((PetscObject)Jac, MATSHELL, &isShell);CHKERRQ(ierr);
825cb1e1211SMatthew G Knepley   if (isShell) {
826cb1e1211SMatthew G Knepley     JacActionCtx *jctx;
827cb1e1211SMatthew G Knepley 
828cb1e1211SMatthew G Knepley     ierr = MatShellGetContext(Jac, &jctx);CHKERRQ(ierr);
829cb1e1211SMatthew G Knepley     ierr = VecCopy(X, jctx->u);CHKERRQ(ierr);
830cb1e1211SMatthew G Knepley   }
831cb1e1211SMatthew G Knepley   *str = SAME_NONZERO_PATTERN;
832cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
833cb1e1211SMatthew G Knepley }
834