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) { 2120cf1dd8SToby Isaac ierr = PetscCalloc4(numCells, &g->face, N * dimEmbed, &g->n, N * dimEmbed * dimEmbed, &(g->suppInvJ[0]), N * dimEmbed * dimEmbed, &(g->suppInvJ[1]));CHKERRQ(ierr); 2220cf1dd8SToby Isaac } 23*4e2f966aSMatthew G. Knepley ierr = PetscCalloc1(N * dimEmbed * dimEmbed, &g->invJ);CHKERRQ(ierr); 2420cf1dd8SToby Isaac *geom = g; 2520cf1dd8SToby Isaac PetscFunctionReturn(0); 2620cf1dd8SToby Isaac } 2720cf1dd8SToby Isaac 2820cf1dd8SToby Isaac PetscErrorCode PetscFEGeomDestroy(PetscFEGeom **geom) 2920cf1dd8SToby Isaac { 3020cf1dd8SToby Isaac PetscErrorCode ierr; 3120cf1dd8SToby Isaac 3220cf1dd8SToby Isaac PetscFunctionBegin; 3320cf1dd8SToby Isaac if (!*geom) PetscFunctionReturn(0); 3420cf1dd8SToby Isaac ierr = PetscFree3((*geom)->v,(*geom)->J,(*geom)->detJ);CHKERRQ(ierr); 3520cf1dd8SToby Isaac ierr = PetscFree((*geom)->invJ);CHKERRQ(ierr); 3620cf1dd8SToby Isaac ierr = PetscFree4((*geom)->face,(*geom)->n,(*geom)->suppInvJ[0],(*geom)->suppInvJ[1]);CHKERRQ(ierr); 3720cf1dd8SToby Isaac ierr = PetscFree(*geom);CHKERRQ(ierr); 3820cf1dd8SToby Isaac PetscFunctionReturn(0); 3920cf1dd8SToby Isaac } 4020cf1dd8SToby Isaac 4120cf1dd8SToby Isaac PetscErrorCode PetscFEGeomGetChunk(PetscFEGeom *geom, PetscInt cStart, PetscInt cEnd, PetscFEGeom **chunkGeom) 4220cf1dd8SToby Isaac { 4320cf1dd8SToby Isaac PetscInt Nq; 4420cf1dd8SToby Isaac PetscInt dE; 4520cf1dd8SToby Isaac PetscErrorCode ierr; 4620cf1dd8SToby Isaac 4720cf1dd8SToby Isaac PetscFunctionBegin; 4820cf1dd8SToby Isaac PetscValidPointer(geom,1); 4920cf1dd8SToby Isaac PetscValidPointer(chunkGeom,2); 5020cf1dd8SToby Isaac if (!(*chunkGeom)) { 5120cf1dd8SToby Isaac ierr = PetscNew(chunkGeom);CHKERRQ(ierr); 5220cf1dd8SToby Isaac } 5320cf1dd8SToby Isaac Nq = geom->numPoints; 5420cf1dd8SToby Isaac dE= geom->dimEmbed; 5520cf1dd8SToby Isaac (*chunkGeom)->dim = geom->dim; 5620cf1dd8SToby Isaac (*chunkGeom)->dimEmbed = geom->dimEmbed; 5720cf1dd8SToby Isaac (*chunkGeom)->numPoints = geom->numPoints; 5820cf1dd8SToby Isaac (*chunkGeom)->numCells = cEnd - cStart; 5920cf1dd8SToby Isaac (*chunkGeom)->xi = geom->xi; 6020cf1dd8SToby Isaac (*chunkGeom)->v = &geom->v[Nq*dE*cStart]; 6120cf1dd8SToby Isaac (*chunkGeom)->J = &geom->J[Nq*dE*dE*cStart]; 6220cf1dd8SToby Isaac (*chunkGeom)->invJ = (geom->invJ) ? &geom->invJ[Nq*dE*dE*cStart] : NULL; 6320cf1dd8SToby Isaac (*chunkGeom)->detJ = &geom->detJ[Nq*cStart]; 6420cf1dd8SToby Isaac (*chunkGeom)->n = geom->n ? &geom->n[Nq*dE*cStart] : NULL; 6520cf1dd8SToby Isaac (*chunkGeom)->face = geom->face ? &geom->face[cStart] : NULL; 6620cf1dd8SToby Isaac (*chunkGeom)->suppInvJ[0] = geom->suppInvJ[0] ? &geom->suppInvJ[0][Nq*dE*dE*cStart] : NULL; 6720cf1dd8SToby Isaac (*chunkGeom)->suppInvJ[1] = geom->suppInvJ[1] ? &geom->suppInvJ[1][Nq*dE*dE*cStart] : NULL; 6820cf1dd8SToby Isaac (*chunkGeom)->isAffine = geom->isAffine; 6920cf1dd8SToby Isaac PetscFunctionReturn(0); 7020cf1dd8SToby Isaac } 7120cf1dd8SToby Isaac 7220cf1dd8SToby Isaac PetscErrorCode PetscFEGeomRestoreChunk(PetscFEGeom *geom, PetscInt cStart, PetscInt cEnd, PetscFEGeom **chunkGeom) 7320cf1dd8SToby Isaac { 7420cf1dd8SToby Isaac PetscErrorCode ierr; 7520cf1dd8SToby Isaac 7620cf1dd8SToby Isaac PetscFunctionBegin; 7720cf1dd8SToby Isaac ierr = PetscFree(*chunkGeom);CHKERRQ(ierr); 7820cf1dd8SToby Isaac PetscFunctionReturn(0); 7920cf1dd8SToby Isaac } 8020cf1dd8SToby Isaac 8120cf1dd8SToby Isaac PetscErrorCode PetscFEGeomComplete(PetscFEGeom *geom) 8220cf1dd8SToby Isaac { 8320cf1dd8SToby Isaac PetscInt i, j, N, dE; 8420cf1dd8SToby Isaac 8520cf1dd8SToby Isaac PetscFunctionBeginHot; 8620cf1dd8SToby Isaac N = geom->numPoints * geom->numCells; 8720cf1dd8SToby Isaac dE = geom->dimEmbed; 8820cf1dd8SToby Isaac switch (dE) { 8920cf1dd8SToby Isaac case 3: 9020cf1dd8SToby Isaac for (i = 0; i < N; i++) { 9120cf1dd8SToby Isaac DMPlex_Det3D_Internal(&geom->detJ[i], &geom->J[dE*dE*i]); 9220cf1dd8SToby Isaac if (geom->invJ) DMPlex_Invert3D_Internal(&geom->invJ[dE*dE*i], &geom->J[dE*dE*i], geom->detJ[i]); 9320cf1dd8SToby Isaac } 9420cf1dd8SToby Isaac break; 9520cf1dd8SToby Isaac case 2: 9620cf1dd8SToby Isaac for (i = 0; i < N; i++) { 9720cf1dd8SToby Isaac DMPlex_Det2D_Internal(&geom->detJ[i], &geom->J[dE*dE*i]); 9820cf1dd8SToby Isaac if (geom->invJ) DMPlex_Invert2D_Internal(&geom->invJ[dE*dE*i], &geom->J[dE*dE*i], geom->detJ[i]); 9920cf1dd8SToby Isaac } 10020cf1dd8SToby Isaac break; 10120cf1dd8SToby Isaac case 1: 10220cf1dd8SToby Isaac for (i = 0; i < N; i++) { 10320cf1dd8SToby Isaac geom->detJ[i] = PetscAbsReal(geom->J[i]); 10420cf1dd8SToby Isaac if (geom->invJ) geom->invJ[i] = 1. / geom->J[i]; 10520cf1dd8SToby Isaac } 10620cf1dd8SToby Isaac break; 10720cf1dd8SToby Isaac } 10820cf1dd8SToby Isaac if (geom->n) { 10920cf1dd8SToby Isaac for (i = 0; i < N; i++) { 11020cf1dd8SToby Isaac for (j = 0; j < dE; j++) { 11120cf1dd8SToby Isaac geom->n[dE*i + j] = geom->J[dE*dE*i + dE*j + dE-1] * ((dE == 2) ? -1. : 1.); 11220cf1dd8SToby Isaac } 11320cf1dd8SToby Isaac } 11420cf1dd8SToby Isaac } 11520cf1dd8SToby Isaac PetscFunctionReturn(0); 11620cf1dd8SToby Isaac } 11720cf1dd8SToby Isaac 118