xref: /petsc/src/dm/dt/fe/interface/fegeom.c (revision 20cf1dd8cd656e27510885a71ea4be028bf48813)
1*20cf1dd8SToby Isaac #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/
2*20cf1dd8SToby Isaac 
3*20cf1dd8SToby Isaac PetscErrorCode PetscFEGeomCreate(PetscQuadrature quad, PetscInt numCells, PetscInt dimEmbed, PetscBool faceData, PetscFEGeom **geom)
4*20cf1dd8SToby Isaac {
5*20cf1dd8SToby Isaac   PetscFEGeom     *g;
6*20cf1dd8SToby Isaac   PetscInt        dim, numPoints, N;
7*20cf1dd8SToby Isaac   const PetscReal *p;
8*20cf1dd8SToby Isaac   PetscErrorCode  ierr;
9*20cf1dd8SToby Isaac 
10*20cf1dd8SToby Isaac   PetscFunctionBegin;
11*20cf1dd8SToby Isaac   ierr = PetscQuadratureGetData(quad,&dim,NULL,&numPoints,&p,NULL);CHKERRQ(ierr);
12*20cf1dd8SToby Isaac   ierr = PetscNew(&g);CHKERRQ(ierr);
13*20cf1dd8SToby Isaac   g->xi        = p;
14*20cf1dd8SToby Isaac   g->numCells  = numCells;
15*20cf1dd8SToby Isaac   g->numPoints = numPoints;
16*20cf1dd8SToby Isaac   g->dim       = dim;
17*20cf1dd8SToby Isaac   g->dimEmbed  = dimEmbed;
18*20cf1dd8SToby Isaac   N = numCells * numPoints;
19*20cf1dd8SToby Isaac   ierr = PetscCalloc3(N * dimEmbed, &g->v, N * dimEmbed * dimEmbed, &g->J, N, &g->detJ);CHKERRQ(ierr);
20*20cf1dd8SToby Isaac   if (faceData) {
21*20cf1dd8SToby Isaac     ierr = PetscCalloc4(numCells, &g->face, N * dimEmbed, &g->n, N * dimEmbed * dimEmbed, &(g->suppInvJ[0]), N * dimEmbed * dimEmbed, &(g->suppInvJ[1]));CHKERRQ(ierr);
22*20cf1dd8SToby Isaac   } else {
23*20cf1dd8SToby Isaac     ierr = PetscCalloc1(N * dimEmbed * dimEmbed, &g->invJ);CHKERRQ(ierr);
24*20cf1dd8SToby Isaac   }
25*20cf1dd8SToby Isaac   *geom = g;
26*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
27*20cf1dd8SToby Isaac }
28*20cf1dd8SToby Isaac 
29*20cf1dd8SToby Isaac PetscErrorCode PetscFEGeomDestroy(PetscFEGeom **geom)
30*20cf1dd8SToby Isaac {
31*20cf1dd8SToby Isaac   PetscErrorCode ierr;
32*20cf1dd8SToby Isaac 
33*20cf1dd8SToby Isaac   PetscFunctionBegin;
34*20cf1dd8SToby Isaac   if (!*geom) PetscFunctionReturn(0);
35*20cf1dd8SToby Isaac   ierr = PetscFree3((*geom)->v,(*geom)->J,(*geom)->detJ);CHKERRQ(ierr);
36*20cf1dd8SToby Isaac   ierr = PetscFree((*geom)->invJ);CHKERRQ(ierr);
37*20cf1dd8SToby Isaac   ierr = PetscFree4((*geom)->face,(*geom)->n,(*geom)->suppInvJ[0],(*geom)->suppInvJ[1]);CHKERRQ(ierr);
38*20cf1dd8SToby Isaac   ierr = PetscFree(*geom);CHKERRQ(ierr);
39*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
40*20cf1dd8SToby Isaac }
41*20cf1dd8SToby Isaac 
42*20cf1dd8SToby Isaac PetscErrorCode PetscFEGeomGetChunk(PetscFEGeom *geom, PetscInt cStart, PetscInt cEnd, PetscFEGeom **chunkGeom)
43*20cf1dd8SToby Isaac {
44*20cf1dd8SToby Isaac   PetscInt       Nq;
45*20cf1dd8SToby Isaac   PetscInt       dE;
46*20cf1dd8SToby Isaac   PetscErrorCode ierr;
47*20cf1dd8SToby Isaac 
48*20cf1dd8SToby Isaac   PetscFunctionBegin;
49*20cf1dd8SToby Isaac   PetscValidPointer(geom,1);
50*20cf1dd8SToby Isaac   PetscValidPointer(chunkGeom,2);
51*20cf1dd8SToby Isaac   if (!(*chunkGeom)) {
52*20cf1dd8SToby Isaac     ierr = PetscNew(chunkGeom);CHKERRQ(ierr);
53*20cf1dd8SToby Isaac   }
54*20cf1dd8SToby Isaac   Nq = geom->numPoints;
55*20cf1dd8SToby Isaac   dE= geom->dimEmbed;
56*20cf1dd8SToby Isaac   (*chunkGeom)->dim = geom->dim;
57*20cf1dd8SToby Isaac   (*chunkGeom)->dimEmbed = geom->dimEmbed;
58*20cf1dd8SToby Isaac   (*chunkGeom)->numPoints = geom->numPoints;
59*20cf1dd8SToby Isaac   (*chunkGeom)->numCells = cEnd - cStart;
60*20cf1dd8SToby Isaac   (*chunkGeom)->xi = geom->xi;
61*20cf1dd8SToby Isaac   (*chunkGeom)->v = &geom->v[Nq*dE*cStart];
62*20cf1dd8SToby Isaac   (*chunkGeom)->J = &geom->J[Nq*dE*dE*cStart];
63*20cf1dd8SToby Isaac   (*chunkGeom)->invJ = (geom->invJ) ? &geom->invJ[Nq*dE*dE*cStart] : NULL;
64*20cf1dd8SToby Isaac   (*chunkGeom)->detJ = &geom->detJ[Nq*cStart];
65*20cf1dd8SToby Isaac   (*chunkGeom)->n = geom->n ? &geom->n[Nq*dE*cStart] : NULL;
66*20cf1dd8SToby Isaac   (*chunkGeom)->face = geom->face ? &geom->face[cStart] : NULL;
67*20cf1dd8SToby Isaac   (*chunkGeom)->suppInvJ[0] = geom->suppInvJ[0] ? &geom->suppInvJ[0][Nq*dE*dE*cStart] : NULL;
68*20cf1dd8SToby Isaac   (*chunkGeom)->suppInvJ[1] = geom->suppInvJ[1] ? &geom->suppInvJ[1][Nq*dE*dE*cStart] : NULL;
69*20cf1dd8SToby Isaac   (*chunkGeom)->isAffine = geom->isAffine;
70*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
71*20cf1dd8SToby Isaac }
72*20cf1dd8SToby Isaac 
73*20cf1dd8SToby Isaac PetscErrorCode PetscFEGeomRestoreChunk(PetscFEGeom *geom, PetscInt cStart, PetscInt cEnd, PetscFEGeom **chunkGeom)
74*20cf1dd8SToby Isaac {
75*20cf1dd8SToby Isaac   PetscErrorCode ierr;
76*20cf1dd8SToby Isaac 
77*20cf1dd8SToby Isaac   PetscFunctionBegin;
78*20cf1dd8SToby Isaac   ierr = PetscFree(*chunkGeom);CHKERRQ(ierr);
79*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
80*20cf1dd8SToby Isaac }
81*20cf1dd8SToby Isaac 
82*20cf1dd8SToby Isaac PetscErrorCode PetscFEGeomComplete(PetscFEGeom *geom)
83*20cf1dd8SToby Isaac {
84*20cf1dd8SToby Isaac   PetscInt i, j, N, dE;
85*20cf1dd8SToby Isaac 
86*20cf1dd8SToby Isaac   PetscFunctionBeginHot;
87*20cf1dd8SToby Isaac   N = geom->numPoints * geom->numCells;
88*20cf1dd8SToby Isaac   dE = geom->dimEmbed;
89*20cf1dd8SToby Isaac   switch (dE) {
90*20cf1dd8SToby Isaac   case 3:
91*20cf1dd8SToby Isaac     for (i = 0; i < N; i++) {
92*20cf1dd8SToby Isaac       DMPlex_Det3D_Internal(&geom->detJ[i], &geom->J[dE*dE*i]);
93*20cf1dd8SToby Isaac       if (geom->invJ) DMPlex_Invert3D_Internal(&geom->invJ[dE*dE*i], &geom->J[dE*dE*i], geom->detJ[i]);
94*20cf1dd8SToby Isaac     }
95*20cf1dd8SToby Isaac     break;
96*20cf1dd8SToby Isaac   case 2:
97*20cf1dd8SToby Isaac     for (i = 0; i < N; i++) {
98*20cf1dd8SToby Isaac       DMPlex_Det2D_Internal(&geom->detJ[i], &geom->J[dE*dE*i]);
99*20cf1dd8SToby Isaac       if (geom->invJ) DMPlex_Invert2D_Internal(&geom->invJ[dE*dE*i], &geom->J[dE*dE*i], geom->detJ[i]);
100*20cf1dd8SToby Isaac     }
101*20cf1dd8SToby Isaac     break;
102*20cf1dd8SToby Isaac   case 1:
103*20cf1dd8SToby Isaac     for (i = 0; i < N; i++) {
104*20cf1dd8SToby Isaac       geom->detJ[i] = PetscAbsReal(geom->J[i]);
105*20cf1dd8SToby Isaac       if (geom->invJ) geom->invJ[i] = 1. / geom->J[i];
106*20cf1dd8SToby Isaac     }
107*20cf1dd8SToby Isaac     break;
108*20cf1dd8SToby Isaac   }
109*20cf1dd8SToby Isaac   if (geom->n) {
110*20cf1dd8SToby Isaac     for (i = 0; i < N; i++) {
111*20cf1dd8SToby Isaac       for (j = 0; j < dE; j++) {
112*20cf1dd8SToby Isaac         geom->n[dE*i + j] = geom->J[dE*dE*i + dE*j + dE-1] * ((dE == 2) ? -1. : 1.);
113*20cf1dd8SToby Isaac       }
114*20cf1dd8SToby Isaac     }
115*20cf1dd8SToby Isaac   }
116*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
117*20cf1dd8SToby Isaac }
118*20cf1dd8SToby Isaac 
119