xref: /petsc/src/dm/impls/plex/plexfem.c (revision c012ea0a8af393ff14026d82ebbcc54675353cb4)
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"
185a1d24da5SMatthew G. Knepley PetscErrorCode DMPlexProjectFunctionLocal(DM dm, PetscInt numComp, void (**funcs)(const PetscReal [], PetscScalar *), InsertMode mode, Vec localX)
186cb1e1211SMatthew G Knepley {
187cb1e1211SMatthew G Knepley   Vec            coordinates;
188cb1e1211SMatthew G Knepley   PetscSection   section, cSection;
189cb1e1211SMatthew G Knepley   PetscInt       dim, vStart, vEnd, v, c, d;
190cb1e1211SMatthew G Knepley   PetscScalar   *values, *cArray;
191cb1e1211SMatthew G Knepley   PetscReal     *coords;
192cb1e1211SMatthew G Knepley   PetscErrorCode ierr;
193cb1e1211SMatthew G Knepley 
194cb1e1211SMatthew G Knepley   PetscFunctionBegin;
195cb1e1211SMatthew G Knepley   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
196cb1e1211SMatthew G Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
197cb1e1211SMatthew G Knepley   ierr = DMPlexGetCoordinateSection(dm, &cSection);CHKERRQ(ierr);
198cb1e1211SMatthew G Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
199cb1e1211SMatthew G Knepley   ierr = PetscMalloc(numComp * sizeof(PetscScalar), &values);CHKERRQ(ierr);
200cb1e1211SMatthew G Knepley   ierr = VecGetArray(coordinates, &cArray);CHKERRQ(ierr);
201cb1e1211SMatthew G Knepley   ierr = PetscSectionGetDof(cSection, vStart, &dim);CHKERRQ(ierr);
202cb1e1211SMatthew G Knepley   ierr = PetscMalloc(dim * sizeof(PetscReal),&coords);CHKERRQ(ierr);
203cb1e1211SMatthew G Knepley   for (v = vStart; v < vEnd; ++v) {
204cb1e1211SMatthew G Knepley     PetscInt dof, off;
205cb1e1211SMatthew G Knepley 
206cb1e1211SMatthew G Knepley     ierr = PetscSectionGetDof(cSection, v, &dof);CHKERRQ(ierr);
207cb1e1211SMatthew G Knepley     ierr = PetscSectionGetOffset(cSection, v, &off);CHKERRQ(ierr);
208cb1e1211SMatthew G Knepley     if (dof > dim) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Cannot have more coordinates %d then dimensions %d", dof, dim);
209cb1e1211SMatthew G Knepley     for (d = 0; d < dof; ++d) coords[d] = PetscRealPart(cArray[off+d]);
210a1d24da5SMatthew G. Knepley     for (c = 0; c < numComp; ++c) (*funcs[c])(coords, &values[c]);
211cb1e1211SMatthew G Knepley     ierr = VecSetValuesSection(localX, section, v, values, mode);CHKERRQ(ierr);
212cb1e1211SMatthew G Knepley   }
213cb1e1211SMatthew G Knepley   ierr = VecRestoreArray(coordinates, &cArray);CHKERRQ(ierr);
214cb1e1211SMatthew G Knepley   /* Temporary, must be replaced by a projection on the finite element basis */
215cb1e1211SMatthew G Knepley   {
216cb1e1211SMatthew G Knepley     PetscInt eStart = 0, eEnd = 0, e, depth;
217cb1e1211SMatthew G Knepley 
218cb1e1211SMatthew G Knepley     ierr = DMPlexGetLabelSize(dm, "depth", &depth);CHKERRQ(ierr);
219cb1e1211SMatthew G Knepley     --depth;
220cb1e1211SMatthew G Knepley     if (depth > 1) {ierr = DMPlexGetDepthStratum(dm, 1, &eStart, &eEnd);CHKERRQ(ierr);}
221cb1e1211SMatthew G Knepley     for (e = eStart; e < eEnd; ++e) {
222cb1e1211SMatthew G Knepley       const PetscInt *cone = NULL;
223cb1e1211SMatthew G Knepley       PetscInt        coneSize, d;
224cb1e1211SMatthew G Knepley       PetscScalar    *coordsA, *coordsB;
225cb1e1211SMatthew G Knepley 
226cb1e1211SMatthew G Knepley       ierr = DMPlexGetConeSize(dm, e, &coneSize);CHKERRQ(ierr);
227cb1e1211SMatthew G Knepley       ierr = DMPlexGetCone(dm, e, &cone);CHKERRQ(ierr);
228cb1e1211SMatthew G Knepley       if (coneSize != 2) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_SIZ, "Cone size %d for point %d should be 2", coneSize, e);
229cb1e1211SMatthew G Knepley       ierr = VecGetValuesSection(coordinates, cSection, cone[0], &coordsA);CHKERRQ(ierr);
230cb1e1211SMatthew G Knepley       ierr = VecGetValuesSection(coordinates, cSection, cone[1], &coordsB);CHKERRQ(ierr);
231cb1e1211SMatthew G Knepley       for (d = 0; d < dim; ++d) {
232cb1e1211SMatthew G Knepley         coords[d] = 0.5*(PetscRealPart(coordsA[d]) + PetscRealPart(coordsB[d]));
233cb1e1211SMatthew G Knepley       }
234a1d24da5SMatthew G. Knepley       for (c = 0; c < numComp; ++c) (*funcs[c])(coords, &values[c]);
235cb1e1211SMatthew G Knepley       ierr = VecSetValuesSection(localX, section, e, values, mode);CHKERRQ(ierr);
236cb1e1211SMatthew G Knepley     }
237cb1e1211SMatthew G Knepley   }
238cb1e1211SMatthew G Knepley 
239cb1e1211SMatthew G Knepley   ierr = PetscFree(coords);CHKERRQ(ierr);
240cb1e1211SMatthew G Knepley   ierr = PetscFree(values);CHKERRQ(ierr);
241cb1e1211SMatthew G Knepley #if 0
242cb1e1211SMatthew G Knepley   const PetscInt localDof = this->_mesh->sizeWithBC(s, *cells->begin());
243cb1e1211SMatthew G Knepley   PetscReal      detJ;
244cb1e1211SMatthew G Knepley 
245cb1e1211SMatthew G Knepley   ierr = PetscMalloc(localDof * sizeof(PetscScalar), &values);CHKERRQ(ierr);
246cb1e1211SMatthew G Knepley   ierr = PetscMalloc2(dim,PetscReal,&v0,dim*dim,PetscReal,&J);CHKERRQ(ierr);
247cb1e1211SMatthew G Knepley   ALE::ISieveVisitor::PointRetriever<PETSC_MESH_TYPE::sieve_type> pV(PetscPowInt(this->_mesh->getSieve()->getMaxConeSize(),dim+1), true);
248cb1e1211SMatthew G Knepley 
249cb1e1211SMatthew G Knepley   for (PetscInt c = cStart; c < cEnd; ++c) {
250cb1e1211SMatthew G Knepley     ALE::ISieveTraversal<PETSC_MESH_TYPE::sieve_type>::orientedClosure(*this->_mesh->getSieve(), c, pV);
251cb1e1211SMatthew G Knepley     const PETSC_MESH_TYPE::point_type *oPoints = pV.getPoints();
252cb1e1211SMatthew G Knepley     const int                          oSize   = pV.getSize();
253cb1e1211SMatthew G Knepley     int                                v       = 0;
254cb1e1211SMatthew G Knepley 
255cb1e1211SMatthew G Knepley     ierr = DMPlexComputeCellGeometry(dm, c, v0, J, NULL, &detJ);CHKERRQ(ierr);
256cb1e1211SMatthew G Knepley     for (PetscInt cl = 0; cl < oSize; ++cl) {
257cb1e1211SMatthew G Knepley       const PetscInt fDim;
258cb1e1211SMatthew G Knepley 
259cb1e1211SMatthew G Knepley       ierr = PetscSectionGetDof(oPoints[cl], &fDim);CHKERRQ(ierr);
260cb1e1211SMatthew G Knepley       if (pointDim) {
261cb1e1211SMatthew G Knepley         for (PetscInt d = 0; d < fDim; ++d, ++v) {
262cb1e1211SMatthew G Knepley           values[v] = (*this->_options.integrate)(v0, J, v, initFunc);
263cb1e1211SMatthew G Knepley         }
264cb1e1211SMatthew G Knepley       }
265cb1e1211SMatthew G Knepley     }
266cb1e1211SMatthew G Knepley     ierr = DMPlexVecSetClosure(dm, NULL, localX, c, values);CHKERRQ(ierr);
267cb1e1211SMatthew G Knepley     pV.clear();
268cb1e1211SMatthew G Knepley   }
269cb1e1211SMatthew G Knepley   ierr = PetscFree2(v0,J);CHKERRQ(ierr);
270cb1e1211SMatthew G Knepley   ierr = PetscFree(values);CHKERRQ(ierr);
271cb1e1211SMatthew G Knepley #endif
272cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
273cb1e1211SMatthew G Knepley }
274cb1e1211SMatthew G Knepley 
275cb1e1211SMatthew G Knepley #undef __FUNCT__
276cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexProjectFunction"
277cb1e1211SMatthew G Knepley /*@C
278cb1e1211SMatthew G Knepley   DMPlexProjectFunction - This projects the given function into the function space provided.
279cb1e1211SMatthew G Knepley 
280cb1e1211SMatthew G Knepley   Input Parameters:
281cb1e1211SMatthew G Knepley + dm      - The DM
282cb1e1211SMatthew G Knepley . numComp - The number of components (functions)
283cb1e1211SMatthew G Knepley . funcs   - The coordinate functions to evaluate
284cb1e1211SMatthew G Knepley - mode    - The insertion mode for values
285cb1e1211SMatthew G Knepley 
286cb1e1211SMatthew G Knepley   Output Parameter:
287cb1e1211SMatthew G Knepley . X - vector
288cb1e1211SMatthew G Knepley 
289cb1e1211SMatthew G Knepley   Level: developer
290cb1e1211SMatthew G Knepley 
291cb1e1211SMatthew G Knepley   Note:
292cb1e1211SMatthew G Knepley   This currently just calls the function with the coordinates of each vertex and edge midpoint, and stores the result in a vector.
293cb1e1211SMatthew G Knepley   We will eventually fix it.
294cb1e1211SMatthew G Knepley 
295878cb397SSatish Balay .seealso: DMPlexComputeL2Diff()
296878cb397SSatish Balay @*/
297a1d24da5SMatthew G. Knepley PetscErrorCode DMPlexProjectFunction(DM dm, PetscInt numComp, void (**funcs)(const PetscReal [], PetscScalar *), InsertMode mode, Vec X)
298cb1e1211SMatthew G Knepley {
299cb1e1211SMatthew G Knepley   Vec            localX;
300cb1e1211SMatthew G Knepley   PetscErrorCode ierr;
301cb1e1211SMatthew G Knepley 
302cb1e1211SMatthew G Knepley   PetscFunctionBegin;
303cb1e1211SMatthew G Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
304cb1e1211SMatthew G Knepley   ierr = DMPlexProjectFunctionLocal(dm, numComp, funcs, mode, localX);CHKERRQ(ierr);
305cb1e1211SMatthew G Knepley   ierr = DMLocalToGlobalBegin(dm, localX, mode, X);CHKERRQ(ierr);
306cb1e1211SMatthew G Knepley   ierr = DMLocalToGlobalEnd(dm, localX, mode, X);CHKERRQ(ierr);
307cb1e1211SMatthew G Knepley   ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
308cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
309cb1e1211SMatthew G Knepley }
310cb1e1211SMatthew G Knepley 
311cb1e1211SMatthew G Knepley #undef __FUNCT__
312cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexComputeL2Diff"
313cb1e1211SMatthew G Knepley /*@C
314cb1e1211SMatthew G Knepley   DMPlexComputeL2Diff - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h.
315cb1e1211SMatthew G Knepley 
316cb1e1211SMatthew G Knepley   Input Parameters:
317cb1e1211SMatthew G Knepley + dm    - The DM
318c5bbbd5bSMatthew G. Knepley . fe    - The PetscFE object for each field
319cb1e1211SMatthew G Knepley . funcs - The functions to evaluate for each field component
320cb1e1211SMatthew G Knepley - X     - The coefficient vector u_h
321cb1e1211SMatthew G Knepley 
322cb1e1211SMatthew G Knepley   Output Parameter:
323cb1e1211SMatthew G Knepley . diff - The diff ||u - u_h||_2
324cb1e1211SMatthew G Knepley 
325cb1e1211SMatthew G Knepley   Level: developer
326cb1e1211SMatthew G Knepley 
327cb1e1211SMatthew G Knepley .seealso: DMPlexProjectFunction()
328878cb397SSatish Balay @*/
329c5bbbd5bSMatthew G. Knepley PetscErrorCode DMPlexComputeL2Diff(DM dm, PetscFE fe[], void (**funcs)(const PetscReal [], PetscScalar *), Vec X, PetscReal *diff)
330cb1e1211SMatthew G Knepley {
331cb1e1211SMatthew G Knepley   const PetscInt  debug = 0;
332cb1e1211SMatthew G Knepley   PetscSection    section;
333c5bbbd5bSMatthew G. Knepley   PetscQuadrature quad;
334cb1e1211SMatthew G Knepley   Vec             localX;
335cb1e1211SMatthew G Knepley   PetscReal      *coords, *v0, *J, *invJ, detJ;
336cb1e1211SMatthew G Knepley   PetscReal       localDiff = 0.0;
337cb1e1211SMatthew G Knepley   PetscInt        dim, numFields, numComponents = 0, cStart, cEnd, c, field, fieldOffset, comp;
338cb1e1211SMatthew G Knepley   PetscErrorCode  ierr;
339cb1e1211SMatthew G Knepley 
340cb1e1211SMatthew G Knepley   PetscFunctionBegin;
341cb1e1211SMatthew G Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
342cb1e1211SMatthew G Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
343cb1e1211SMatthew G Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
344cb1e1211SMatthew G Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
345cb1e1211SMatthew G Knepley   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
346cb1e1211SMatthew G Knepley   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
347cb1e1211SMatthew G Knepley   for (field = 0; field < numFields; ++field) {
348c5bbbd5bSMatthew G. Knepley     PetscInt Nc;
349c5bbbd5bSMatthew G. Knepley 
350c5bbbd5bSMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe[field], &Nc);CHKERRQ(ierr);
351c5bbbd5bSMatthew G. Knepley     numComponents += Nc;
352cb1e1211SMatthew G Knepley   }
353cb1e1211SMatthew G Knepley   ierr = DMPlexProjectFunctionLocal(dm, numComponents, funcs, INSERT_BC_VALUES, localX);CHKERRQ(ierr);
354cb1e1211SMatthew G Knepley   ierr = PetscMalloc4(dim,PetscReal,&coords,dim,PetscReal,&v0,dim*dim,PetscReal,&J,dim*dim,PetscReal,&invJ);CHKERRQ(ierr);
355cb1e1211SMatthew G Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
356c5bbbd5bSMatthew G. Knepley   ierr = PetscFEGetQuadrature(fe[0], &quad);CHKERRQ(ierr);
357cb1e1211SMatthew G Knepley   for (c = cStart; c < cEnd; ++c) {
358a1e44745SMatthew G. Knepley     PetscScalar *x = NULL;
359cb1e1211SMatthew G Knepley     PetscReal    elemDiff = 0.0;
360cb1e1211SMatthew G Knepley 
361cb1e1211SMatthew G Knepley     ierr = DMPlexComputeCellGeometry(dm, c, v0, J, invJ, &detJ);CHKERRQ(ierr);
362cb1e1211SMatthew G Knepley     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
363cb1e1211SMatthew G Knepley     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
364cb1e1211SMatthew G Knepley 
365cb1e1211SMatthew G Knepley     for (field = 0, comp = 0, fieldOffset = 0; field < numFields; ++field) {
366c5bbbd5bSMatthew G. Knepley       const PetscInt   numQuadPoints = quad.numQuadPoints;
367c5bbbd5bSMatthew G. Knepley       const PetscReal *quadPoints    = quad.quadPoints;
368c5bbbd5bSMatthew G. Knepley       const PetscReal *quadWeights   = quad.quadWeights;
369c5bbbd5bSMatthew G. Knepley       PetscReal       *basis;
370c5bbbd5bSMatthew G. Knepley       PetscInt         numBasisFuncs, numBasisComps, q, d, e, fc, f;
371cb1e1211SMatthew G Knepley 
372c5bbbd5bSMatthew G. Knepley       ierr = PetscFEGetDimension(fe[field], &numBasisFuncs);CHKERRQ(ierr);
373c5bbbd5bSMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe[field], &numBasisComps);CHKERRQ(ierr);
374c5bbbd5bSMatthew G. Knepley       ierr = PetscFEGetDefaultTabulation(fe[field], &basis, NULL, NULL);CHKERRQ(ierr);
375cb1e1211SMatthew G Knepley       if (debug) {
376cb1e1211SMatthew G Knepley         char title[1024];
377cb1e1211SMatthew G Knepley         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
378cb1e1211SMatthew G Knepley         ierr = DMPrintCellVector(c, title, numBasisFuncs*numBasisComps, &x[fieldOffset]);CHKERRQ(ierr);
379cb1e1211SMatthew G Knepley       }
380cb1e1211SMatthew G Knepley       for (q = 0; q < numQuadPoints; ++q) {
381cb1e1211SMatthew G Knepley         for (d = 0; d < dim; d++) {
382cb1e1211SMatthew G Knepley           coords[d] = v0[d];
383cb1e1211SMatthew G Knepley           for (e = 0; e < dim; e++) {
384cb1e1211SMatthew G Knepley             coords[d] += J[d*dim+e]*(quadPoints[q*dim+e] + 1.0);
385cb1e1211SMatthew G Knepley           }
386cb1e1211SMatthew G Knepley         }
387cb1e1211SMatthew G Knepley         for (fc = 0; fc < numBasisComps; ++fc) {
388a1d24da5SMatthew G. Knepley           PetscScalar funcVal;
389a1d24da5SMatthew G. Knepley           PetscScalar interpolant = 0.0;
390a1d24da5SMatthew G. Knepley 
391a1d24da5SMatthew G. Knepley           (*funcs[comp+fc])(coords, &funcVal);
392cb1e1211SMatthew G Knepley           for (f = 0; f < numBasisFuncs; ++f) {
393cb1e1211SMatthew G Knepley             const PetscInt fidx = f*numBasisComps+fc;
394a1d24da5SMatthew G. Knepley             interpolant += x[fieldOffset+fidx]*basis[q*numBasisFuncs*numBasisComps+fidx];
395cb1e1211SMatthew G Knepley           }
396a1d24da5SMatthew 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);}
397a1d24da5SMatthew G. Knepley           elemDiff += PetscSqr(PetscRealPart(interpolant - funcVal))*quadWeights[q]*detJ;
398cb1e1211SMatthew G Knepley         }
399cb1e1211SMatthew G Knepley       }
400cb1e1211SMatthew G Knepley       comp        += numBasisComps;
401cb1e1211SMatthew G Knepley       fieldOffset += numBasisFuncs*numBasisComps;
402cb1e1211SMatthew G Knepley     }
403cb1e1211SMatthew G Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
404cb1e1211SMatthew G Knepley     if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);}
405cb1e1211SMatthew G Knepley     localDiff += elemDiff;
406cb1e1211SMatthew G Knepley   }
407cb1e1211SMatthew G Knepley   ierr  = PetscFree4(coords,v0,J,invJ);CHKERRQ(ierr);
408cb1e1211SMatthew G Knepley   ierr  = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
40986a74ee0SMatthew G. Knepley   ierr  = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
410cb1e1211SMatthew G Knepley   *diff = PetscSqrtReal(*diff);
411cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
412cb1e1211SMatthew G Knepley }
413cb1e1211SMatthew G Knepley 
414a0845e3aSMatthew G. Knepley #if 0
415a0845e3aSMatthew G. Knepley 
416cb1e1211SMatthew G Knepley #undef __FUNCT__
417cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexComputeResidualFEM"
418cb1e1211SMatthew G Knepley PetscErrorCode DMPlexComputeResidualFEM(DM dm, Vec X, Vec F, void *user)
419cb1e1211SMatthew G Knepley {
420cb1e1211SMatthew G Knepley   DM_Plex         *mesh   = (DM_Plex *) dm->data;
4219a559087SMatthew G. Knepley   PetscFEM        *fem    = (PetscFEM *) user;
422cb1e1211SMatthew G Knepley   PetscQuadrature *quad   = fem->quad;
423652b88e8SMatthew G. Knepley   PetscQuadrature *quadBd = fem->quadBd;
424cb1e1211SMatthew G Knepley   PetscSection     section;
42502a80efeSMatthew G. Knepley   PetscReal       *v0, *n, *J, *invJ, *detJ;
426cb1e1211SMatthew G Knepley   PetscScalar     *elemVec, *u;
427cb1e1211SMatthew G Knepley   PetscInt         dim, numFields, field, numBatchesTmp = 1, numCells, cStart, cEnd, c;
428652b88e8SMatthew G. Knepley   PetscInt         cellDof, numComponents;
429652b88e8SMatthew G. Knepley   PetscBool        has;
430cb1e1211SMatthew G Knepley   PetscErrorCode   ierr;
431cb1e1211SMatthew G Knepley 
432cb1e1211SMatthew G Knepley   PetscFunctionBegin;
433652b88e8SMatthew G. Knepley   if (has && quadBd) {
434652b88e8SMatthew G. Knepley     DMLabel         label;
435652b88e8SMatthew G. Knepley     IS              pointIS;
436652b88e8SMatthew G. Knepley     const PetscInt *points;
437652b88e8SMatthew G. Knepley     PetscInt        numPoints, p;
438652b88e8SMatthew G. Knepley 
439652b88e8SMatthew G. Knepley     ierr = DMPlexGetLabel(dm, "boundary", &label);CHKERRQ(ierr);
440652b88e8SMatthew G. Knepley     ierr = DMLabelGetStratumSize(label, 1, &numPoints);CHKERRQ(ierr);
441652b88e8SMatthew G. Knepley     ierr = DMLabelGetStratumIS(label, 1, &pointIS);CHKERRQ(ierr);
442652b88e8SMatthew G. Knepley     ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
443652b88e8SMatthew G. Knepley     for (field = 0, cellDof = 0, numComponents = 0; field < numFields; ++field) {
444652b88e8SMatthew G. Knepley       cellDof       += quadBd[field].numBasisFuncs*quadBd[field].numComponents;
445652b88e8SMatthew G. Knepley       numComponents += quadBd[field].numComponents;
446652b88e8SMatthew G. Knepley     }
44702a80efeSMatthew 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);
448652b88e8SMatthew G. Knepley     for (p = 0; p < numPoints; ++p) {
449652b88e8SMatthew G. Knepley       const PetscInt point = points[p];
450652b88e8SMatthew G. Knepley       PetscScalar   *x;
451652b88e8SMatthew G. Knepley       PetscInt       i;
452652b88e8SMatthew G. Knepley 
45302a80efeSMatthew G. Knepley       /* TODO: Add normal determination here */
454652b88e8SMatthew G. Knepley       ierr = DMPlexComputeCellGeometry(dm, point, &v0[p*dim], &J[p*dim*dim], &invJ[p*dim*dim], &detJ[p]);CHKERRQ(ierr);
4551d930511SMatthew G. Knepley       if (detJ[p] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for face %d", detJ[p], point);
456652b88e8SMatthew G. Knepley       ierr = DMPlexVecGetClosure(dm, NULL, X, point, NULL, &x);CHKERRQ(ierr);
457652b88e8SMatthew G. Knepley 
458652b88e8SMatthew G. Knepley       for (i = 0; i < cellDof; ++i) u[p*cellDof+i] = x[i];
459652b88e8SMatthew G. Knepley       ierr = DMPlexVecRestoreClosure(dm, NULL, X, point, NULL, &x);CHKERRQ(ierr);
460652b88e8SMatthew G. Knepley     }
461652b88e8SMatthew G. Knepley     for (field = 0; field < numFields; ++field) {
462652b88e8SMatthew G. Knepley       const PetscInt numQuadPoints = quadBd[field].numQuadPoints;
463652b88e8SMatthew G. Knepley       const PetscInt numBasisFuncs = quadBd[field].numBasisFuncs;
464a9dc2124SMatthew G. Knepley       void           (*f0)(const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->f0BdFuncs[field];
465a9dc2124SMatthew G. Knepley       void           (*f1)(const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->f1BdFuncs[field];
466652b88e8SMatthew G. Knepley       /* Conforming batches */
467652b88e8SMatthew G. Knepley       PetscInt blockSize  = numBasisFuncs*numQuadPoints;
468652b88e8SMatthew G. Knepley       PetscInt numBlocks  = 1;
469652b88e8SMatthew G. Knepley       PetscInt batchSize  = numBlocks * blockSize;
470652b88e8SMatthew G. Knepley       PetscInt numBatches = numBatchesTmp;
471652b88e8SMatthew G. Knepley       PetscInt numChunks  = numPoints / (numBatches*batchSize);
472652b88e8SMatthew G. Knepley       /* Remainder */
473652b88e8SMatthew G. Knepley       PetscInt numRemainder = numPoints % (numBatches * batchSize);
474652b88e8SMatthew G. Knepley       PetscInt offset       = numPoints - numRemainder;
475652b88e8SMatthew G. Knepley 
47602a80efeSMatthew G. Knepley       ierr = (*mesh->integrateBdResidualFEM)(numChunks*numBatches*batchSize, numFields, field, quadBd, u, v0, n, J, invJ, detJ, f0, f1, elemVec);CHKERRQ(ierr);
47702a80efeSMatthew 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],
478652b88e8SMatthew G. Knepley                                              f0, f1, &elemVec[offset*cellDof]);CHKERRQ(ierr);
479652b88e8SMatthew G. Knepley     }
480652b88e8SMatthew G. Knepley     for (p = 0; p < numPoints; ++p) {
481652b88e8SMatthew G. Knepley       const PetscInt point = points[p];
482652b88e8SMatthew G. Knepley 
483652b88e8SMatthew G. Knepley       if (mesh->printFEM > 1) {ierr = DMPrintCellVector(point, "Residual", cellDof, &elemVec[p*cellDof]);CHKERRQ(ierr);}
484652b88e8SMatthew G. Knepley       ierr = DMPlexVecSetClosure(dm, NULL, F, point, &elemVec[p*cellDof], ADD_VALUES);CHKERRQ(ierr);
485652b88e8SMatthew G. Knepley     }
486652b88e8SMatthew G. Knepley     ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
487652b88e8SMatthew G. Knepley     ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
48802a80efeSMatthew G. Knepley     ierr = PetscFree7(u,v0,n,J,invJ,detJ,elemVec);CHKERRQ(ierr);
489652b88e8SMatthew G. Knepley   }
490cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
491cb1e1211SMatthew G Knepley }
492cb1e1211SMatthew G Knepley 
493a0845e3aSMatthew G. Knepley #else
494a0845e3aSMatthew G. Knepley 
495a0845e3aSMatthew G. Knepley #undef __FUNCT__
496a0845e3aSMatthew G. Knepley #define __FUNCT__ "DMPlexComputeResidualFEM"
497a0845e3aSMatthew G. Knepley /*@
498a0845e3aSMatthew G. Knepley   DMPlexComputeResidualFEM - Form the local residual F from the local input X using pointwise functions specified by the user
499a0845e3aSMatthew G. Knepley 
500a0845e3aSMatthew G. Knepley   Input Parameters:
501a0845e3aSMatthew G. Knepley + dm - The mesh
502a0845e3aSMatthew G. Knepley . X  - Local input vector
503a0845e3aSMatthew G. Knepley - user - The user context
504a0845e3aSMatthew G. Knepley 
505a0845e3aSMatthew G. Knepley   Output Parameter:
506a0845e3aSMatthew G. Knepley . F  - Local output vector
507a0845e3aSMatthew G. Knepley 
508a0845e3aSMatthew G. Knepley   Note:
509a0845e3aSMatthew G. Knepley   The second member of the user context must be an FEMContext.
510a0845e3aSMatthew G. Knepley 
511a0845e3aSMatthew G. Knepley   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
512a0845e3aSMatthew G. Knepley   like a GPU, or vectorize on a multicore machine.
513a0845e3aSMatthew G. Knepley 
514a0845e3aSMatthew G. Knepley   Level: developer
515a0845e3aSMatthew G. Knepley 
516a0845e3aSMatthew G. Knepley .seealso: DMPlexComputeJacobianActionFEM()
517a0845e3aSMatthew G. Knepley @*/
518a0845e3aSMatthew G. Knepley PetscErrorCode DMPlexComputeResidualFEM(DM dm, Vec X, Vec F, void *user)
519a0845e3aSMatthew G. Knepley {
520a0845e3aSMatthew G. Knepley   DM_Plex          *mesh  = (DM_Plex *) dm->data;
5219a559087SMatthew G. Knepley   PetscFEM         *fem   = (PetscFEM *) user;
522a0845e3aSMatthew G. Knepley   PetscFE          *fe    = fem->fe;
5239a559087SMatthew G. Knepley   PetscFE          *feAux = fem->feAux;
524f1ea0e2fSMatthew G. Knepley   PetscFE          *feBd  = fem->feBd;
525a0845e3aSMatthew G. Knepley   const char       *name  = "Residual";
5269a559087SMatthew G. Knepley   DM                dmAux;
5279a559087SMatthew G. Knepley   Vec               A;
528a0845e3aSMatthew G. Knepley   PetscQuadrature   q;
529a0845e3aSMatthew G. Knepley   PetscCellGeometry geom;
5309a559087SMatthew G. Knepley   PetscSection      section, sectionAux;
531a0845e3aSMatthew G. Knepley   PetscReal        *v0, *J, *invJ, *detJ;
5329a559087SMatthew G. Knepley   PetscScalar      *elemVec, *u, *a;
5339a559087SMatthew G. Knepley   PetscInt          dim, Nf, NfAux = 0, f, numCells, cStart, cEnd, c;
534a0845e3aSMatthew G. Knepley   PetscInt          cellDof = 0, numComponents = 0;
5359a559087SMatthew G. Knepley   PetscInt          cellDofAux = 0, numComponentsAux = 0;
536a0845e3aSMatthew G. Knepley   PetscErrorCode    ierr;
537a0845e3aSMatthew G. Knepley 
538a0845e3aSMatthew G. Knepley   PetscFunctionBegin;
539a0845e3aSMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_ResidualFEM,dm,0,0,0);CHKERRQ(ierr);
540a0845e3aSMatthew G. Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
541a0845e3aSMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
5429a559087SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
543a0845e3aSMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
544a0845e3aSMatthew G. Knepley   numCells = cEnd - cStart;
5459a559087SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
546a0845e3aSMatthew G. Knepley     PetscInt Nb, Nc;
547a0845e3aSMatthew G. Knepley 
548a0845e3aSMatthew G. Knepley     ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr);
549a0845e3aSMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe[f], &Nc);CHKERRQ(ierr);
550a0845e3aSMatthew G. Knepley     cellDof       += Nb*Nc;
551a0845e3aSMatthew G. Knepley     numComponents += Nc;
552a0845e3aSMatthew G. Knepley   }
5539a559087SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
5549a559087SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr);
5559a559087SMatthew G. Knepley   if (dmAux) {
5569a559087SMatthew G. Knepley     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
5579a559087SMatthew G. Knepley     ierr = PetscSectionGetNumFields(sectionAux, &NfAux);CHKERRQ(ierr);
5589a559087SMatthew G. Knepley   }
5599a559087SMatthew G. Knepley   for (f = 0; f < NfAux; ++f) {
5609a559087SMatthew G. Knepley     PetscInt Nb, Nc;
5619a559087SMatthew G. Knepley 
5629a559087SMatthew G. Knepley     ierr = PetscFEGetDimension(feAux[f], &Nb);CHKERRQ(ierr);
5639a559087SMatthew G. Knepley     ierr = PetscFEGetNumComponents(feAux[f], &Nc);CHKERRQ(ierr);
5649a559087SMatthew G. Knepley     cellDofAux       += Nb*Nc;
5659a559087SMatthew G. Knepley     numComponentsAux += Nc;
5669a559087SMatthew G. Knepley   }
567a0845e3aSMatthew G. Knepley   ierr = DMPlexProjectFunctionLocal(dm, numComponents, fem->bcFuncs, INSERT_BC_VALUES, X);CHKERRQ(ierr);
568a0845e3aSMatthew G. Knepley   ierr = VecSet(F, 0.0);CHKERRQ(ierr);
569a0845e3aSMatthew 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);
5709a559087SMatthew G. Knepley   if (dmAux) {ierr = PetscMalloc(numCells*cellDofAux * sizeof(PetscScalar), &a);CHKERRQ(ierr);}
571a0845e3aSMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
572a0845e3aSMatthew G. Knepley     PetscScalar *x = NULL;
573a0845e3aSMatthew G. Knepley     PetscInt     i;
574a0845e3aSMatthew G. Knepley 
575a0845e3aSMatthew G. Knepley     ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr);
576a0845e3aSMatthew G. Knepley     if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c);
577a0845e3aSMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr);
578a0845e3aSMatthew G. Knepley     for (i = 0; i < cellDof; ++i) u[c*cellDof+i] = x[i];
579a0845e3aSMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr);
5809a559087SMatthew G. Knepley     if (dmAux) {
5819a559087SMatthew G. Knepley       ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
5829a559087SMatthew G. Knepley       for (i = 0; i < cellDofAux; ++i) a[c*cellDofAux+i] = x[i];
5839a559087SMatthew G. Knepley       ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
584a0845e3aSMatthew G. Knepley     }
5859a559087SMatthew G. Knepley   }
5869a559087SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
587*c012ea0aSMatthew G. Knepley     void   (*f0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->f0Funcs[f];
588*c012ea0aSMatthew G. Knepley     void   (*f1)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->f1Funcs[f];
589a0845e3aSMatthew G. Knepley     PetscInt Nb;
590a0845e3aSMatthew G. Knepley     /* Conforming batches */
591f30c5766SMatthew G. Knepley     PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
592a0845e3aSMatthew G. Knepley     /* Remainder */
593a0845e3aSMatthew G. Knepley     PetscInt Nr, offset;
594a0845e3aSMatthew G. Knepley 
595a0845e3aSMatthew G. Knepley     ierr = PetscFEGetQuadrature(fe[f], &q);CHKERRQ(ierr);
596a0845e3aSMatthew G. Knepley     ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr);
597f30c5766SMatthew G. Knepley     ierr = PetscFEGetTileSizes(fe[f], NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
598a0845e3aSMatthew G. Knepley     blockSize = Nb*q.numQuadPoints;
599a0845e3aSMatthew G. Knepley     batchSize = numBlocks * blockSize;
600f30c5766SMatthew G. Knepley     ierr =  PetscFESetTileSizes(fe[f], blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
601a0845e3aSMatthew G. Knepley     numChunks = numCells / (numBatches*batchSize);
602a0845e3aSMatthew G. Knepley     Ne        = numChunks*numBatches*batchSize;
603a0845e3aSMatthew G. Knepley     Nr        = numCells % (numBatches*batchSize);
604a0845e3aSMatthew G. Knepley     offset    = numCells - Nr;
605a0845e3aSMatthew G. Knepley     geom.v0   = v0;
606a0845e3aSMatthew G. Knepley     geom.J    = J;
607a0845e3aSMatthew G. Knepley     geom.invJ = invJ;
608a0845e3aSMatthew G. Knepley     geom.detJ = detJ;
6099a559087SMatthew G. Knepley     ierr = PetscFEIntegrateResidual(fe[f], Ne, Nf, fe, f, geom, u, NfAux, feAux, a, f0, f1, elemVec);CHKERRQ(ierr);
610a0845e3aSMatthew G. Knepley     geom.v0   = &v0[offset*dim];
611a0845e3aSMatthew G. Knepley     geom.J    = &J[offset*dim*dim];
612a0845e3aSMatthew G. Knepley     geom.invJ = &invJ[offset*dim*dim];
613a0845e3aSMatthew G. Knepley     geom.detJ = &detJ[offset];
6149a559087SMatthew 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);
615a0845e3aSMatthew G. Knepley   }
616a0845e3aSMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
617a0845e3aSMatthew G. Knepley     if (mesh->printFEM > 1) {ierr = DMPrintCellVector(c, name, cellDof, &elemVec[c*cellDof]);CHKERRQ(ierr);}
618a0845e3aSMatthew G. Knepley     ierr = DMPlexVecSetClosure(dm, section, F, c, &elemVec[c*cellDof], ADD_VALUES);CHKERRQ(ierr);
619a0845e3aSMatthew G. Knepley   }
620a0845e3aSMatthew G. Knepley   ierr = PetscFree6(u,v0,J,invJ,detJ,elemVec);CHKERRQ(ierr);
6219a559087SMatthew G. Knepley   if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);}
622f1ea0e2fSMatthew G. Knepley   if (feBd) {
623f1ea0e2fSMatthew G. Knepley     DMLabel         label;
624f1ea0e2fSMatthew G. Knepley     IS              pointIS;
625f1ea0e2fSMatthew G. Knepley     const PetscInt *points;
626f1ea0e2fSMatthew G. Knepley     PetscInt        numPoints, p;
627f1ea0e2fSMatthew G. Knepley     PetscReal      *n;
628f1ea0e2fSMatthew G. Knepley 
629f1ea0e2fSMatthew G. Knepley     ierr = DMPlexGetLabel(dm, "boundary", &label);CHKERRQ(ierr);
630f1ea0e2fSMatthew G. Knepley     ierr = DMLabelGetStratumSize(label, 1, &numPoints);CHKERRQ(ierr);
631f1ea0e2fSMatthew G. Knepley     ierr = DMLabelGetStratumIS(label, 1, &pointIS);CHKERRQ(ierr);
632f1ea0e2fSMatthew G. Knepley     ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
633f1ea0e2fSMatthew G. Knepley     for (f = 0, cellDof = 0, numComponents = 0; f < Nf; ++f) {
634f1ea0e2fSMatthew G. Knepley       PetscInt Nb, Nc;
635f1ea0e2fSMatthew G. Knepley 
636f1ea0e2fSMatthew G. Knepley       ierr = PetscFEGetDimension(feBd[f], &Nb);CHKERRQ(ierr);
637f1ea0e2fSMatthew G. Knepley       ierr = PetscFEGetNumComponents(feBd[f], &Nc);CHKERRQ(ierr);
638f1ea0e2fSMatthew G. Knepley       cellDof       += Nb*Nc;
639f1ea0e2fSMatthew G. Knepley       numComponents += Nc;
640f1ea0e2fSMatthew G. Knepley     }
641f1ea0e2fSMatthew 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);
642f1ea0e2fSMatthew G. Knepley     for (p = 0; p < numPoints; ++p) {
643f1ea0e2fSMatthew G. Knepley       const PetscInt point = points[p];
644f1ea0e2fSMatthew G. Knepley       PetscScalar   *x     = NULL;
645f1ea0e2fSMatthew G. Knepley       PetscInt       i;
646f1ea0e2fSMatthew G. Knepley 
647f1ea0e2fSMatthew G. Knepley       /* TODO: Add normal determination here */
648f1ea0e2fSMatthew G. Knepley       ierr = DMPlexComputeCellGeometry(dm, point, &v0[p*dim], &J[p*dim*dim], &invJ[p*dim*dim], &detJ[p]);CHKERRQ(ierr);
649f1ea0e2fSMatthew G. Knepley       if (detJ[p] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for face %d", detJ[p], point);
650f1ea0e2fSMatthew G. Knepley       ierr = DMPlexVecGetClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr);
651f1ea0e2fSMatthew G. Knepley       for (i = 0; i < cellDof; ++i) u[p*cellDof+i] = x[i];
652f1ea0e2fSMatthew G. Knepley       ierr = DMPlexVecRestoreClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr);
653f1ea0e2fSMatthew G. Knepley     }
654f1ea0e2fSMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
655f1ea0e2fSMatthew G. Knepley       void   (*f0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->f0BdFuncs[f];
656f1ea0e2fSMatthew G. Knepley       void   (*f1)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->f1BdFuncs[f];
657f1ea0e2fSMatthew G. Knepley       PetscInt Nb;
658f1ea0e2fSMatthew G. Knepley       /* Conforming batches */
659f1ea0e2fSMatthew G. Knepley       PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
660f1ea0e2fSMatthew G. Knepley       /* Remainder */
661f1ea0e2fSMatthew G. Knepley       PetscInt Nr, offset;
662f1ea0e2fSMatthew G. Knepley 
663f1ea0e2fSMatthew G. Knepley       ierr = PetscFEGetQuadrature(feBd[f], &q);CHKERRQ(ierr);
664f1ea0e2fSMatthew G. Knepley       ierr = PetscFEGetDimension(feBd[f], &Nb);CHKERRQ(ierr);
665f1ea0e2fSMatthew G. Knepley       ierr = PetscFEGetTileSizes(feBd[f], NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
666f1ea0e2fSMatthew G. Knepley       blockSize = Nb*q.numQuadPoints;
667f1ea0e2fSMatthew G. Knepley       batchSize = numBlocks * blockSize;
668f1ea0e2fSMatthew G. Knepley       ierr =  PetscFESetTileSizes(feBd[f], blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
669f1ea0e2fSMatthew G. Knepley       numChunks = numPoints / (numBatches*batchSize);
670f1ea0e2fSMatthew G. Knepley       Ne        = numChunks*numBatches*batchSize;
671f1ea0e2fSMatthew G. Knepley       Nr        = numPoints % (numBatches*batchSize);
672f1ea0e2fSMatthew G. Knepley       offset    = numPoints - Nr;
673f1ea0e2fSMatthew G. Knepley       geom.v0   = v0;
674f1ea0e2fSMatthew G. Knepley       geom.n    = n;
675f1ea0e2fSMatthew G. Knepley       geom.J    = J;
676f1ea0e2fSMatthew G. Knepley       geom.invJ = invJ;
677f1ea0e2fSMatthew G. Knepley       geom.detJ = detJ;
678f1ea0e2fSMatthew G. Knepley       ierr = PetscFEIntegrateBdResidual(feBd[f], Ne, Nf, feBd, f, geom, u, 0, NULL, NULL, f0, f1, elemVec);CHKERRQ(ierr);
679f1ea0e2fSMatthew G. Knepley       geom.v0   = &v0[offset*dim];
680f1ea0e2fSMatthew G. Knepley       geom.n    = &n[offset*dim];
681f1ea0e2fSMatthew G. Knepley       geom.J    = &J[offset*dim*dim];
682f1ea0e2fSMatthew G. Knepley       geom.invJ = &invJ[offset*dim*dim];
683f1ea0e2fSMatthew G. Knepley       geom.detJ = &detJ[offset];
684f1ea0e2fSMatthew G. Knepley       ierr = PetscFEIntegrateBdResidual(feBd[f], Nr, Nf, feBd, f, geom, &u[offset*cellDof], 0, NULL, NULL, f0, f1, &elemVec[offset*cellDof]);CHKERRQ(ierr);
685f1ea0e2fSMatthew G. Knepley     }
686f1ea0e2fSMatthew G. Knepley     for (p = 0; p < numPoints; ++p) {
687f1ea0e2fSMatthew G. Knepley       const PetscInt point = points[p];
688f1ea0e2fSMatthew G. Knepley 
689f1ea0e2fSMatthew G. Knepley       if (mesh->printFEM > 1) {ierr = DMPrintCellVector(point, "BdResidual", cellDof, &elemVec[p*cellDof]);CHKERRQ(ierr);}
690f1ea0e2fSMatthew G. Knepley       ierr = DMPlexVecSetClosure(dm, NULL, F, point, &elemVec[p*cellDof], ADD_VALUES);CHKERRQ(ierr);
691f1ea0e2fSMatthew G. Knepley     }
692f1ea0e2fSMatthew G. Knepley     ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
693f1ea0e2fSMatthew G. Knepley     ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
694f1ea0e2fSMatthew G. Knepley     ierr = PetscFree7(u,v0,n,J,invJ,detJ,elemVec);CHKERRQ(ierr);
695f1ea0e2fSMatthew G. Knepley   }
6966113b454SMatthew G. Knepley   if (mesh->printFEM) {ierr = DMPrintLocalVec(dm, name, mesh->printTol, F);CHKERRQ(ierr);}
697a0845e3aSMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_ResidualFEM,dm,0,0,0);CHKERRQ(ierr);
698a0845e3aSMatthew G. Knepley   PetscFunctionReturn(0);
699a0845e3aSMatthew G. Knepley }
700a0845e3aSMatthew G. Knepley 
701a0845e3aSMatthew G. Knepley #endif
702a0845e3aSMatthew G. Knepley 
703cb1e1211SMatthew G Knepley #undef __FUNCT__
704cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexComputeJacobianActionFEM"
705cb1e1211SMatthew G Knepley /*@C
706cb1e1211SMatthew G Knepley   DMPlexComputeJacobianActionFEM - Form the local action of Jacobian J(u) on the local input X using pointwise functions specified by the user
707cb1e1211SMatthew G Knepley 
708cb1e1211SMatthew G Knepley   Input Parameters:
709cb1e1211SMatthew G Knepley + dm - The mesh
710cb1e1211SMatthew G Knepley . J  - The Jacobian shell matrix
711cb1e1211SMatthew G Knepley . X  - Local input vector
712cb1e1211SMatthew G Knepley - user - The user context
713cb1e1211SMatthew G Knepley 
714cb1e1211SMatthew G Knepley   Output Parameter:
715cb1e1211SMatthew G Knepley . F  - Local output vector
716cb1e1211SMatthew G Knepley 
717cb1e1211SMatthew G Knepley   Note:
718cb1e1211SMatthew G Knepley   The second member of the user context must be an FEMContext.
719cb1e1211SMatthew G Knepley 
720cb1e1211SMatthew G Knepley   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
721cb1e1211SMatthew G Knepley   like a GPU, or vectorize on a multicore machine.
722cb1e1211SMatthew G Knepley 
7230059ad2aSSatish Balay   Level: developer
7240059ad2aSSatish Balay 
725cb1e1211SMatthew G Knepley .seealso: DMPlexComputeResidualFEM()
726878cb397SSatish Balay @*/
727cb1e1211SMatthew G Knepley PetscErrorCode DMPlexComputeJacobianActionFEM(DM dm, Mat Jac, Vec X, Vec F, void *user)
728cb1e1211SMatthew G Knepley {
729cb1e1211SMatthew G Knepley   DM_Plex          *mesh = (DM_Plex *) dm->data;
7309a559087SMatthew G. Knepley   PetscFEM         *fem  = (PetscFEM *) user;
7310483ade4SMatthew G. Knepley   PetscFE          *fe   = fem->fe;
7320483ade4SMatthew G. Knepley   PetscQuadrature   quad;
7330483ade4SMatthew G. Knepley   PetscCellGeometry geom;
734cb1e1211SMatthew G Knepley   PetscSection      section;
735cb1e1211SMatthew G Knepley   JacActionCtx     *jctx;
736cb1e1211SMatthew G Knepley   PetscReal        *v0, *J, *invJ, *detJ;
737cb1e1211SMatthew G Knepley   PetscScalar      *elemVec, *u, *a;
7380483ade4SMatthew G. Knepley   PetscInt          dim, numFields, field, numCells, cStart, cEnd, c;
739cb1e1211SMatthew G Knepley   PetscInt          cellDof = 0;
740cb1e1211SMatthew G Knepley   PetscErrorCode    ierr;
741cb1e1211SMatthew G Knepley 
742cb1e1211SMatthew G Knepley   PetscFunctionBegin;
7430483ade4SMatthew G. Knepley   /* ierr = PetscLogEventBegin(DMPLEX_JacobianActionFEM,dm,0,0,0);CHKERRQ(ierr); */
744cb1e1211SMatthew G Knepley   ierr = MatShellGetContext(Jac, &jctx);CHKERRQ(ierr);
745cb1e1211SMatthew G Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
746cb1e1211SMatthew G Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
747cb1e1211SMatthew G Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
748cb1e1211SMatthew G Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
749cb1e1211SMatthew G Knepley   numCells = cEnd - cStart;
750cb1e1211SMatthew G Knepley   for (field = 0; field < numFields; ++field) {
7510483ade4SMatthew G. Knepley     PetscInt Nb, Nc;
7520483ade4SMatthew G. Knepley 
7530483ade4SMatthew G. Knepley     ierr = PetscFEGetDimension(fe[field], &Nb);CHKERRQ(ierr);
7540483ade4SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe[field], &Nc);CHKERRQ(ierr);
7550483ade4SMatthew G. Knepley     cellDof += Nb*Nc;
756cb1e1211SMatthew G Knepley   }
757cb1e1211SMatthew G Knepley   ierr = VecSet(F, 0.0);CHKERRQ(ierr);
758cb1e1211SMatthew 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);
759cb1e1211SMatthew G Knepley   for (c = cStart; c < cEnd; ++c) {
760a1e44745SMatthew G. Knepley     PetscScalar *x = NULL;
761cb1e1211SMatthew G Knepley     PetscInt     i;
762cb1e1211SMatthew G Knepley 
763cb1e1211SMatthew G Knepley     ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr);
764cb1e1211SMatthew G Knepley     if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c);
765cb1e1211SMatthew G Knepley     ierr = DMPlexVecGetClosure(dm, NULL, jctx->u, c, NULL, &x);CHKERRQ(ierr);
766cb1e1211SMatthew G Knepley     for (i = 0; i < cellDof; ++i) u[c*cellDof+i] = x[i];
767cb1e1211SMatthew G Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, jctx->u, c, NULL, &x);CHKERRQ(ierr);
768cb1e1211SMatthew G Knepley     ierr = DMPlexVecGetClosure(dm, NULL, X, c, NULL, &x);CHKERRQ(ierr);
769cb1e1211SMatthew G Knepley     for (i = 0; i < cellDof; ++i) a[c*cellDof+i] = x[i];
770cb1e1211SMatthew G Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, X, c, NULL, &x);CHKERRQ(ierr);
771cb1e1211SMatthew G Knepley   }
772cb1e1211SMatthew G Knepley   for (field = 0; field < numFields; ++field) {
7730483ade4SMatthew G. Knepley     PetscInt Nb;
774cb1e1211SMatthew G Knepley     /* Conforming batches */
775cb1e1211SMatthew G Knepley     PetscInt numBlocks  = 1;
7760483ade4SMatthew G. Knepley     PetscInt numBatches = 1;
7770483ade4SMatthew G. Knepley     PetscInt numChunks, Ne, blockSize, batchSize;
778cb1e1211SMatthew G Knepley     /* Remainder */
7790483ade4SMatthew G. Knepley     PetscInt Nr, offset;
780cb1e1211SMatthew G Knepley 
7810483ade4SMatthew G. Knepley     ierr = PetscFEGetQuadrature(fe[field], &quad);CHKERRQ(ierr);
7820483ade4SMatthew G. Knepley     ierr = PetscFEGetDimension(fe[field], &Nb);CHKERRQ(ierr);
7830483ade4SMatthew G. Knepley     blockSize = Nb*quad.numQuadPoints;
7840483ade4SMatthew G. Knepley     batchSize = numBlocks * blockSize;
7850483ade4SMatthew G. Knepley     numChunks = numCells / (numBatches*batchSize);
7860483ade4SMatthew G. Knepley     Ne        = numChunks*numBatches*batchSize;
7870483ade4SMatthew G. Knepley     Nr        = numCells % (numBatches*batchSize);
7880483ade4SMatthew G. Knepley     offset    = numCells - Nr;
7890483ade4SMatthew G. Knepley     geom.v0   = v0;
7900483ade4SMatthew G. Knepley     geom.J    = J;
7910483ade4SMatthew G. Knepley     geom.invJ = invJ;
7920483ade4SMatthew G. Knepley     geom.detJ = detJ;
7930483ade4SMatthew G. Knepley     ierr = PetscFEIntegrateJacobianAction(fe[field], Ne, numFields, fe, field, geom, u, a, fem->g0Funcs, fem->g1Funcs, fem->g2Funcs, fem->g3Funcs, elemVec);CHKERRQ(ierr);
7940483ade4SMatthew G. Knepley     geom.v0   = &v0[offset*dim];
7950483ade4SMatthew G. Knepley     geom.J    = &J[offset*dim*dim];
7960483ade4SMatthew G. Knepley     geom.invJ = &invJ[offset*dim*dim];
7970483ade4SMatthew G. Knepley     geom.detJ = &detJ[offset];
7980483ade4SMatthew G. Knepley     ierr = PetscFEIntegrateJacobianAction(fe[field], Nr, numFields, fe, field, geom, &u[offset*cellDof], &a[offset*cellDof],
799cb1e1211SMatthew G Knepley                                           fem->g0Funcs, fem->g1Funcs, fem->g2Funcs, fem->g3Funcs, &elemVec[offset*cellDof]);CHKERRQ(ierr);
800cb1e1211SMatthew G Knepley   }
801cb1e1211SMatthew G Knepley   for (c = cStart; c < cEnd; ++c) {
802cb1e1211SMatthew G Knepley     if (mesh->printFEM > 1) {ierr = DMPrintCellVector(c, "Jacobian Action", cellDof, &elemVec[c*cellDof]);CHKERRQ(ierr);}
803cb1e1211SMatthew G Knepley     ierr = DMPlexVecSetClosure(dm, NULL, F, c, &elemVec[c*cellDof], ADD_VALUES);CHKERRQ(ierr);
804cb1e1211SMatthew G Knepley   }
805cb1e1211SMatthew G Knepley   ierr = PetscFree7(u,a,v0,J,invJ,detJ,elemVec);CHKERRQ(ierr);
806cb1e1211SMatthew G Knepley   if (mesh->printFEM) {
807cb1e1211SMatthew G Knepley     PetscMPIInt rank, numProcs;
808cb1e1211SMatthew G Knepley     PetscInt    p;
809cb1e1211SMatthew G Knepley 
810cb1e1211SMatthew G Knepley     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr);
811cb1e1211SMatthew G Knepley     ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm), &numProcs);CHKERRQ(ierr);
81286a74ee0SMatthew G. Knepley     ierr = PetscPrintf(PetscObjectComm((PetscObject)dm), "Jacobian Action:\n");CHKERRQ(ierr);
813cb1e1211SMatthew G Knepley     for (p = 0; p < numProcs; ++p) {
814cb1e1211SMatthew G Knepley       if (p == rank) {ierr = VecView(F, PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);}
815cb1e1211SMatthew G Knepley       ierr = PetscBarrier((PetscObject) dm);CHKERRQ(ierr);
816cb1e1211SMatthew G Knepley     }
817cb1e1211SMatthew G Knepley   }
8180483ade4SMatthew G. Knepley   /* ierr = PetscLogEventEnd(DMPLEX_JacobianActionFEM,dm,0,0,0);CHKERRQ(ierr); */
819cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
820cb1e1211SMatthew G Knepley }
821cb1e1211SMatthew G Knepley 
822cb1e1211SMatthew G Knepley #undef __FUNCT__
823cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexComputeJacobianFEM"
824cb1e1211SMatthew G Knepley /*@
825cb1e1211SMatthew G Knepley   DMPlexComputeJacobianFEM - Form the local portion of the Jacobian matrix J at the local solution X using pointwise functions specified by the user.
826cb1e1211SMatthew G Knepley 
827cb1e1211SMatthew G Knepley   Input Parameters:
828cb1e1211SMatthew G Knepley + dm - The mesh
829cb1e1211SMatthew G Knepley . X  - Local input vector
830cb1e1211SMatthew G Knepley - user - The user context
831cb1e1211SMatthew G Knepley 
832cb1e1211SMatthew G Knepley   Output Parameter:
833cb1e1211SMatthew G Knepley . Jac  - Jacobian matrix
834cb1e1211SMatthew G Knepley 
835cb1e1211SMatthew G Knepley   Note:
836cb1e1211SMatthew G Knepley   The second member of the user context must be an FEMContext.
837cb1e1211SMatthew G Knepley 
838cb1e1211SMatthew G Knepley   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
839cb1e1211SMatthew G Knepley   like a GPU, or vectorize on a multicore machine.
840cb1e1211SMatthew G Knepley 
8410059ad2aSSatish Balay   Level: developer
8420059ad2aSSatish Balay 
843cb1e1211SMatthew G Knepley .seealso: FormFunctionLocal()
844878cb397SSatish Balay @*/
845cb1e1211SMatthew G Knepley PetscErrorCode DMPlexComputeJacobianFEM(DM dm, Vec X, Mat Jac, Mat JacP, MatStructure *str,void *user)
846cb1e1211SMatthew G Knepley {
847cb1e1211SMatthew G Knepley   DM_Plex          *mesh = (DM_Plex *) dm->data;
8489a559087SMatthew G. Knepley   PetscFEM         *fem  = (PetscFEM *) user;
849a319912fSMatthew G. Knepley   PetscFE          *fe   = fem->fe;
850a319912fSMatthew G. Knepley   const char       *name = "Jacobian";
851a319912fSMatthew G. Knepley   PetscQuadrature   quad;
852a319912fSMatthew G. Knepley   PetscCellGeometry geom;
853a319912fSMatthew G. Knepley   PetscSection      section, globalSection;
854cb1e1211SMatthew G Knepley   PetscReal        *v0, *J, *invJ, *detJ;
855cb1e1211SMatthew G Knepley   PetscScalar      *elemMat, *u;
856a319912fSMatthew G. Knepley   PetscInt          dim, numFields, f, fieldI, fieldJ, numCells, cStart, cEnd, c;
857cb1e1211SMatthew G Knepley   PetscInt          cellDof = 0, numComponents = 0;
858cb1e1211SMatthew G Knepley   PetscBool         isShell;
859cb1e1211SMatthew G Knepley   PetscErrorCode    ierr;
860cb1e1211SMatthew G Knepley 
861cb1e1211SMatthew G Knepley   PetscFunctionBegin;
862a319912fSMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_JacobianFEM,dm,0,0,0);CHKERRQ(ierr);
863cb1e1211SMatthew G Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
864cb1e1211SMatthew G Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
865a319912fSMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dm, &globalSection);CHKERRQ(ierr);
866cb1e1211SMatthew G Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
867cb1e1211SMatthew G Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
868cb1e1211SMatthew G Knepley   numCells = cEnd - cStart;
869a319912fSMatthew G. Knepley   for (f = 0; f < numFields; ++f) {
870a319912fSMatthew G. Knepley     PetscInt Nb, Nc;
871a319912fSMatthew G. Knepley 
872a319912fSMatthew G. Knepley     ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr);
873a319912fSMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe[f], &Nc);CHKERRQ(ierr);
874a319912fSMatthew G. Knepley     cellDof       += Nb*Nc;
875a319912fSMatthew G. Knepley     numComponents += Nc;
876cb1e1211SMatthew G Knepley   }
877cb1e1211SMatthew G Knepley   ierr = DMPlexProjectFunctionLocal(dm, numComponents, fem->bcFuncs, INSERT_BC_VALUES, X);CHKERRQ(ierr);
878cb1e1211SMatthew G Knepley   ierr = MatZeroEntries(JacP);CHKERRQ(ierr);
879cb1e1211SMatthew 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);
880cb1e1211SMatthew G Knepley   for (c = cStart; c < cEnd; ++c) {
881a1e44745SMatthew G. Knepley     PetscScalar *x = NULL;
882cb1e1211SMatthew G Knepley     PetscInt     i;
883cb1e1211SMatthew G Knepley 
884cb1e1211SMatthew G Knepley     ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr);
885cb1e1211SMatthew G Knepley     if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c);
886a319912fSMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr);
887cb1e1211SMatthew G Knepley     for (i = 0; i < cellDof; ++i) u[c*cellDof+i] = x[i];
888a319912fSMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr);
889cb1e1211SMatthew G Knepley   }
890cb1e1211SMatthew G Knepley   ierr = PetscMemzero(elemMat, numCells*cellDof*cellDof * sizeof(PetscScalar));CHKERRQ(ierr);
891cb1e1211SMatthew G Knepley   for (fieldI = 0; fieldI < numFields; ++fieldI) {
892a319912fSMatthew G. Knepley     PetscInt Nb;
893a319912fSMatthew G. Knepley     ierr = PetscFEGetQuadrature(fe[fieldI], &quad);CHKERRQ(ierr);
894a319912fSMatthew G. Knepley     ierr = PetscFEGetDimension(fe[fieldI], &Nb);CHKERRQ(ierr);
895cb1e1211SMatthew G Knepley     for (fieldJ = 0; fieldJ < numFields; ++fieldJ) {
896a0845e3aSMatthew G. Knepley       void   (*g0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g0Funcs[fieldI*numFields+fieldJ];
897a0845e3aSMatthew G. Knepley       void   (*g1)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g1Funcs[fieldI*numFields+fieldJ];
898a0845e3aSMatthew G. Knepley       void   (*g2)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g2Funcs[fieldI*numFields+fieldJ];
899a0845e3aSMatthew G. Knepley       void   (*g3)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g3Funcs[fieldI*numFields+fieldJ];
900cb1e1211SMatthew G Knepley       /* Conforming batches */
901cb1e1211SMatthew G Knepley       PetscInt numBlocks  = 1;
902a319912fSMatthew G. Knepley       PetscInt numBatches = 1;
903a319912fSMatthew G. Knepley       PetscInt numChunks, Ne, blockSize, batchSize;
904cb1e1211SMatthew G Knepley       /* Remainder */
905a319912fSMatthew G. Knepley       PetscInt Nr, offset;
906cb1e1211SMatthew G Knepley 
907a319912fSMatthew G. Knepley       blockSize = Nb*quad.numQuadPoints;
908a319912fSMatthew G. Knepley       batchSize = numBlocks * blockSize;
909a319912fSMatthew G. Knepley       numChunks = numCells / (numBatches*batchSize);
910a319912fSMatthew G. Knepley       Ne        = numChunks*numBatches*batchSize;
911a319912fSMatthew G. Knepley       Nr        = numCells % (numBatches*batchSize);
912a319912fSMatthew G. Knepley       offset    = numCells - Nr;
913a319912fSMatthew G. Knepley       geom.v0   = v0;
914a319912fSMatthew G. Knepley       geom.J    = J;
915a319912fSMatthew G. Knepley       geom.invJ = invJ;
916a319912fSMatthew G. Knepley       geom.detJ = detJ;
9170483ade4SMatthew G. Knepley       ierr = PetscFEIntegrateJacobian(fe[fieldI], Ne, numFields, fe, fieldI, fieldJ, geom, u, g0, g1, g2, g3, elemMat);CHKERRQ(ierr);
918a319912fSMatthew G. Knepley       geom.v0   = &v0[offset*dim];
919a319912fSMatthew G. Knepley       geom.J    = &J[offset*dim*dim];
920a319912fSMatthew G. Knepley       geom.invJ = &invJ[offset*dim*dim];
921a319912fSMatthew G. Knepley       geom.detJ = &detJ[offset];
9220483ade4SMatthew G. Knepley       ierr = PetscFEIntegrateJacobian(fe[fieldI], Nr, numFields, fe, fieldI, fieldJ, geom, &u[offset*cellDof], g0, g1, g2, g3, &elemMat[offset*cellDof*cellDof]);CHKERRQ(ierr);
923cb1e1211SMatthew G Knepley     }
924cb1e1211SMatthew G Knepley   }
925cb1e1211SMatthew G Knepley   for (c = cStart; c < cEnd; ++c) {
926a319912fSMatthew G. Knepley     if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(c, name, cellDof, cellDof, &elemMat[c*cellDof*cellDof]);CHKERRQ(ierr);}
927a319912fSMatthew G. Knepley     ierr = DMPlexMatSetClosure(dm, section, globalSection, JacP, c, &elemMat[c*cellDof*cellDof], ADD_VALUES);CHKERRQ(ierr);
928cb1e1211SMatthew G Knepley   }
929cb1e1211SMatthew G Knepley   ierr = PetscFree6(u,v0,J,invJ,detJ,elemMat);CHKERRQ(ierr);
930cb1e1211SMatthew G Knepley   ierr = MatAssemblyBegin(JacP, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
931cb1e1211SMatthew G Knepley   ierr = MatAssemblyEnd(JacP, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
932cb1e1211SMatthew G Knepley   if (mesh->printFEM) {
933a319912fSMatthew G. Knepley     ierr = PetscPrintf(PETSC_COMM_WORLD, "%s:\n", name);CHKERRQ(ierr);
934cb1e1211SMatthew G Knepley     ierr = MatChop(JacP, 1.0e-10);CHKERRQ(ierr);
935cb1e1211SMatthew G Knepley     ierr = MatView(JacP, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
936cb1e1211SMatthew G Knepley   }
937a319912fSMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_JacobianFEM,dm,0,0,0);CHKERRQ(ierr);
938cb1e1211SMatthew G Knepley   ierr = PetscObjectTypeCompare((PetscObject) Jac, MATSHELL, &isShell);CHKERRQ(ierr);
939cb1e1211SMatthew G Knepley   if (isShell) {
940cb1e1211SMatthew G Knepley     JacActionCtx *jctx;
941cb1e1211SMatthew G Knepley 
942cb1e1211SMatthew G Knepley     ierr = MatShellGetContext(Jac, &jctx);CHKERRQ(ierr);
943cb1e1211SMatthew G Knepley     ierr = VecCopy(X, jctx->u);CHKERRQ(ierr);
944cb1e1211SMatthew G Knepley   }
945cb1e1211SMatthew G Knepley   *str = SAME_NONZERO_PATTERN;
946cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
947cb1e1211SMatthew G Knepley }
948