xref: /petsc/src/dm/dt/fe/interface/fegeom.c (revision 414189874cb4f5068fee41cf4cec5928a602cb44)
120cf1dd8SToby Isaac #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/
220cf1dd8SToby Isaac 
320cf1dd8SToby Isaac PetscErrorCode PetscFEGeomCreate(PetscQuadrature quad, PetscInt numCells, PetscInt dimEmbed, PetscBool faceData, PetscFEGeom **geom)
420cf1dd8SToby Isaac {
520cf1dd8SToby Isaac   PetscFEGeom     *g;
67489efa5SMatthew G. Knepley   PetscInt        dim, Nq, N;
720cf1dd8SToby Isaac   const PetscReal *p;
820cf1dd8SToby Isaac   PetscErrorCode  ierr;
920cf1dd8SToby Isaac 
1020cf1dd8SToby Isaac   PetscFunctionBegin;
117489efa5SMatthew G. Knepley   ierr = PetscQuadratureGetData(quad,&dim,NULL,&Nq,&p,NULL);CHKERRQ(ierr);
1220cf1dd8SToby Isaac   ierr = PetscNew(&g);CHKERRQ(ierr);
1320cf1dd8SToby Isaac   g->xi        = p;
1420cf1dd8SToby Isaac   g->numCells  = numCells;
157489efa5SMatthew G. Knepley   g->numPoints = Nq;
1620cf1dd8SToby Isaac   g->dim       = dim;
1720cf1dd8SToby Isaac   g->dimEmbed  = dimEmbed;
187489efa5SMatthew G. Knepley   N = numCells * Nq;
1920cf1dd8SToby Isaac   ierr = PetscCalloc3(N * dimEmbed, &g->v, N * dimEmbed * dimEmbed, &g->J, N, &g->detJ);CHKERRQ(ierr);
2020cf1dd8SToby Isaac   if (faceData) {
21*41418987SMatthew G. Knepley     ierr = PetscCalloc2(numCells, &g->face, N * dimEmbed, &g->n);CHKERRQ(ierr);
22*41418987SMatthew G. Knepley     ierr = PetscCalloc6(N * dimEmbed * dimEmbed, &(g->suppJ[0]),    N * dimEmbed * dimEmbed, &(g->suppJ[1]),
23*41418987SMatthew G. Knepley                         N * dimEmbed * dimEmbed, &(g->suppInvJ[0]), N * dimEmbed * dimEmbed, &(g->suppInvJ[1]),
24*41418987SMatthew G. Knepley                         N,                       &(g->suppDetJ[0]), N,                       &(g->suppDetJ[1]));CHKERRQ(ierr);
2520cf1dd8SToby Isaac   }
264e2f966aSMatthew G. Knepley   ierr = PetscCalloc1(N * dimEmbed * dimEmbed, &g->invJ);CHKERRQ(ierr);
2720cf1dd8SToby Isaac   *geom = g;
2820cf1dd8SToby Isaac   PetscFunctionReturn(0);
2920cf1dd8SToby Isaac }
3020cf1dd8SToby Isaac 
3120cf1dd8SToby Isaac PetscErrorCode PetscFEGeomDestroy(PetscFEGeom **geom)
3220cf1dd8SToby Isaac {
3320cf1dd8SToby Isaac   PetscErrorCode ierr;
3420cf1dd8SToby Isaac 
3520cf1dd8SToby Isaac   PetscFunctionBegin;
3620cf1dd8SToby Isaac   if (!*geom) PetscFunctionReturn(0);
3720cf1dd8SToby Isaac   ierr = PetscFree3((*geom)->v,(*geom)->J,(*geom)->detJ);CHKERRQ(ierr);
3820cf1dd8SToby Isaac   ierr = PetscFree((*geom)->invJ);CHKERRQ(ierr);
39*41418987SMatthew G. Knepley   ierr = PetscFree2((*geom)->face,(*geom)->n);CHKERRQ(ierr);
40*41418987SMatthew G. Knepley   ierr = PetscFree6((*geom)->suppJ[0],(*geom)->suppJ[1],(*geom)->suppInvJ[0],(*geom)->suppInvJ[1],(*geom)->suppDetJ[0],(*geom)->suppDetJ[1]);CHKERRQ(ierr);
4120cf1dd8SToby Isaac   ierr = PetscFree(*geom);CHKERRQ(ierr);
4220cf1dd8SToby Isaac   PetscFunctionReturn(0);
4320cf1dd8SToby Isaac }
4420cf1dd8SToby Isaac 
4520cf1dd8SToby Isaac PetscErrorCode PetscFEGeomGetChunk(PetscFEGeom *geom, PetscInt cStart, PetscInt cEnd, PetscFEGeom **chunkGeom)
4620cf1dd8SToby Isaac {
4720cf1dd8SToby Isaac   PetscInt       Nq;
4820cf1dd8SToby Isaac   PetscInt       dE;
4920cf1dd8SToby Isaac   PetscErrorCode ierr;
5020cf1dd8SToby Isaac 
5120cf1dd8SToby Isaac   PetscFunctionBegin;
5220cf1dd8SToby Isaac   PetscValidPointer(geom,1);
5320cf1dd8SToby Isaac   PetscValidPointer(chunkGeom,2);
5420cf1dd8SToby Isaac   if (!(*chunkGeom)) {
5520cf1dd8SToby Isaac     ierr = PetscNew(chunkGeom);CHKERRQ(ierr);
5620cf1dd8SToby Isaac   }
5720cf1dd8SToby Isaac   Nq = geom->numPoints;
5820cf1dd8SToby Isaac   dE= geom->dimEmbed;
5920cf1dd8SToby Isaac   (*chunkGeom)->dim = geom->dim;
6020cf1dd8SToby Isaac   (*chunkGeom)->dimEmbed = geom->dimEmbed;
6120cf1dd8SToby Isaac   (*chunkGeom)->numPoints = geom->numPoints;
6220cf1dd8SToby Isaac   (*chunkGeom)->numCells = cEnd - cStart;
6320cf1dd8SToby Isaac   (*chunkGeom)->xi = geom->xi;
6420cf1dd8SToby Isaac   (*chunkGeom)->v = &geom->v[Nq*dE*cStart];
6520cf1dd8SToby Isaac   (*chunkGeom)->J = &geom->J[Nq*dE*dE*cStart];
6620cf1dd8SToby Isaac   (*chunkGeom)->invJ = (geom->invJ) ? &geom->invJ[Nq*dE*dE*cStart] : NULL;
6720cf1dd8SToby Isaac   (*chunkGeom)->detJ = &geom->detJ[Nq*cStart];
6820cf1dd8SToby Isaac   (*chunkGeom)->n = geom->n ? &geom->n[Nq*dE*cStart] : NULL;
6920cf1dd8SToby Isaac   (*chunkGeom)->face = geom->face ? &geom->face[cStart] : NULL;
70*41418987SMatthew G. Knepley   (*chunkGeom)->suppJ[0]    = geom->suppJ[0]    ? &geom->suppJ[0][Nq*dE*dE*cStart]    : NULL;
71*41418987SMatthew G. Knepley   (*chunkGeom)->suppJ[1]    = geom->suppJ[1]    ? &geom->suppJ[1][Nq*dE*dE*cStart]    : NULL;
7220cf1dd8SToby Isaac   (*chunkGeom)->suppInvJ[0] = geom->suppInvJ[0] ? &geom->suppInvJ[0][Nq*dE*dE*cStart] : NULL;
7320cf1dd8SToby Isaac   (*chunkGeom)->suppInvJ[1] = geom->suppInvJ[1] ? &geom->suppInvJ[1][Nq*dE*dE*cStart] : NULL;
74*41418987SMatthew G. Knepley   (*chunkGeom)->suppDetJ[0] = geom->suppDetJ[0] ? &geom->suppDetJ[0][Nq*cStart]       : NULL;
75*41418987SMatthew G. Knepley   (*chunkGeom)->suppDetJ[1] = geom->suppDetJ[1] ? &geom->suppDetJ[1][Nq*cStart]       : NULL;
7620cf1dd8SToby Isaac   (*chunkGeom)->isAffine = geom->isAffine;
7720cf1dd8SToby Isaac   PetscFunctionReturn(0);
7820cf1dd8SToby Isaac }
7920cf1dd8SToby Isaac 
8020cf1dd8SToby Isaac PetscErrorCode PetscFEGeomRestoreChunk(PetscFEGeom *geom, PetscInt cStart, PetscInt cEnd, PetscFEGeom **chunkGeom)
8120cf1dd8SToby Isaac {
8220cf1dd8SToby Isaac   PetscErrorCode ierr;
8320cf1dd8SToby Isaac 
8420cf1dd8SToby Isaac   PetscFunctionBegin;
8520cf1dd8SToby Isaac   ierr = PetscFree(*chunkGeom);CHKERRQ(ierr);
8620cf1dd8SToby Isaac   PetscFunctionReturn(0);
8720cf1dd8SToby Isaac }
8820cf1dd8SToby Isaac 
8920cf1dd8SToby Isaac PetscErrorCode PetscFEGeomComplete(PetscFEGeom *geom)
9020cf1dd8SToby Isaac {
9120cf1dd8SToby Isaac   PetscInt i, j, N, dE;
9220cf1dd8SToby Isaac 
9320cf1dd8SToby Isaac   PetscFunctionBeginHot;
9420cf1dd8SToby Isaac   N = geom->numPoints * geom->numCells;
9520cf1dd8SToby Isaac   dE = geom->dimEmbed;
9620cf1dd8SToby Isaac   switch (dE) {
9720cf1dd8SToby Isaac   case 3:
9820cf1dd8SToby Isaac     for (i = 0; i < N; i++) {
9920cf1dd8SToby Isaac       DMPlex_Det3D_Internal(&geom->detJ[i], &geom->J[dE*dE*i]);
10020cf1dd8SToby Isaac       if (geom->invJ) DMPlex_Invert3D_Internal(&geom->invJ[dE*dE*i], &geom->J[dE*dE*i], geom->detJ[i]);
10120cf1dd8SToby Isaac     }
10220cf1dd8SToby Isaac     break;
10320cf1dd8SToby Isaac   case 2:
10420cf1dd8SToby Isaac     for (i = 0; i < N; i++) {
10520cf1dd8SToby Isaac       DMPlex_Det2D_Internal(&geom->detJ[i], &geom->J[dE*dE*i]);
10620cf1dd8SToby Isaac       if (geom->invJ) DMPlex_Invert2D_Internal(&geom->invJ[dE*dE*i], &geom->J[dE*dE*i], geom->detJ[i]);
10720cf1dd8SToby Isaac     }
10820cf1dd8SToby Isaac     break;
10920cf1dd8SToby Isaac   case 1:
11020cf1dd8SToby Isaac     for (i = 0; i < N; i++) {
11120cf1dd8SToby Isaac       geom->detJ[i] = PetscAbsReal(geom->J[i]);
11220cf1dd8SToby Isaac       if (geom->invJ) geom->invJ[i] = 1. / geom->J[i];
11320cf1dd8SToby Isaac     }
11420cf1dd8SToby Isaac     break;
11520cf1dd8SToby Isaac   }
11620cf1dd8SToby Isaac   if (geom->n) {
11720cf1dd8SToby Isaac     for (i = 0; i < N; i++) {
11820cf1dd8SToby Isaac       for (j = 0; j < dE; j++) {
11920cf1dd8SToby Isaac         geom->n[dE*i + j] = geom->J[dE*dE*i + dE*j + dE-1] * ((dE == 2) ? -1. : 1.);
12020cf1dd8SToby Isaac       }
12120cf1dd8SToby Isaac     }
12220cf1dd8SToby Isaac   }
12320cf1dd8SToby Isaac   PetscFunctionReturn(0);
12420cf1dd8SToby Isaac }
125