120cf1dd8SToby Isaac #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/ 220cf1dd8SToby Isaac 3*2b99622eSMatthew G. Knepley /*@C 4*2b99622eSMatthew G. Knepley PetscFEGeomCreate - Create a PetscFEGeom object to manage geometry for a group of cells 5*2b99622eSMatthew G. Knepley 6*2b99622eSMatthew G. Knepley Input Parameters: 7*2b99622eSMatthew G. Knepley + quad - A PetscQuadrature determining the tabulation 8*2b99622eSMatthew G. Knepley . numCells - The number of cells in the group 9*2b99622eSMatthew G. Knepley . dimEmbed - The coordinate dimension 10*2b99622eSMatthew G. Knepley - faceData - Flag to construct geometry data for the faces 11*2b99622eSMatthew G. Knepley 12*2b99622eSMatthew G. Knepley Output Parameter: 13*2b99622eSMatthew G. Knepley . geom - The PetscFEGeom object 14*2b99622eSMatthew G. Knepley 15*2b99622eSMatthew G. Knepley Level: beginner 16*2b99622eSMatthew G. Knepley 17*2b99622eSMatthew G. Knepley .seealso: PetscFEGeomDestroy(), PetscFEGeomComplete() 18*2b99622eSMatthew G. Knepley @*/ 1920cf1dd8SToby Isaac PetscErrorCode PetscFEGeomCreate(PetscQuadrature quad, PetscInt numCells, PetscInt dimEmbed, PetscBool faceData, PetscFEGeom **geom) 2020cf1dd8SToby Isaac { 2120cf1dd8SToby Isaac PetscFEGeom *g; 227489efa5SMatthew G. Knepley PetscInt dim, Nq, N; 2320cf1dd8SToby Isaac const PetscReal *p; 2420cf1dd8SToby Isaac PetscErrorCode ierr; 2520cf1dd8SToby Isaac 2620cf1dd8SToby Isaac PetscFunctionBegin; 277489efa5SMatthew G. Knepley ierr = PetscQuadratureGetData(quad,&dim,NULL,&Nq,&p,NULL);CHKERRQ(ierr); 2820cf1dd8SToby Isaac ierr = PetscNew(&g);CHKERRQ(ierr); 2920cf1dd8SToby Isaac g->xi = p; 3020cf1dd8SToby Isaac g->numCells = numCells; 317489efa5SMatthew G. Knepley g->numPoints = Nq; 3220cf1dd8SToby Isaac g->dim = dim; 3320cf1dd8SToby Isaac g->dimEmbed = dimEmbed; 347489efa5SMatthew G. Knepley N = numCells * Nq; 3520cf1dd8SToby Isaac ierr = PetscCalloc3(N * dimEmbed, &g->v, N * dimEmbed * dimEmbed, &g->J, N, &g->detJ);CHKERRQ(ierr); 3620cf1dd8SToby Isaac if (faceData) { 3741418987SMatthew G. Knepley ierr = PetscCalloc2(numCells, &g->face, N * dimEmbed, &g->n);CHKERRQ(ierr); 3841418987SMatthew G. Knepley ierr = PetscCalloc6(N * dimEmbed * dimEmbed, &(g->suppJ[0]), N * dimEmbed * dimEmbed, &(g->suppJ[1]), 3941418987SMatthew G. Knepley N * dimEmbed * dimEmbed, &(g->suppInvJ[0]), N * dimEmbed * dimEmbed, &(g->suppInvJ[1]), 4041418987SMatthew G. Knepley N, &(g->suppDetJ[0]), N, &(g->suppDetJ[1]));CHKERRQ(ierr); 4120cf1dd8SToby Isaac } 424e2f966aSMatthew G. Knepley ierr = PetscCalloc1(N * dimEmbed * dimEmbed, &g->invJ);CHKERRQ(ierr); 4320cf1dd8SToby Isaac *geom = g; 4420cf1dd8SToby Isaac PetscFunctionReturn(0); 4520cf1dd8SToby Isaac } 4620cf1dd8SToby Isaac 47*2b99622eSMatthew G. Knepley /*@C 48*2b99622eSMatthew G. Knepley PetscFEGeomDestroy - Destroy a PetscFEGeom object 49*2b99622eSMatthew G. Knepley 50*2b99622eSMatthew G. Knepley Input Parameter: 51*2b99622eSMatthew G. Knepley . geom - PetscFEGeom object 52*2b99622eSMatthew G. Knepley 53*2b99622eSMatthew G. Knepley Level: beginner 54*2b99622eSMatthew G. Knepley 55*2b99622eSMatthew G. Knepley .seealso: PetscFEGeomCreate() 56*2b99622eSMatthew G. Knepley @*/ 5720cf1dd8SToby Isaac PetscErrorCode PetscFEGeomDestroy(PetscFEGeom **geom) 5820cf1dd8SToby Isaac { 5920cf1dd8SToby Isaac PetscErrorCode ierr; 6020cf1dd8SToby Isaac 6120cf1dd8SToby Isaac PetscFunctionBegin; 6220cf1dd8SToby Isaac if (!*geom) PetscFunctionReturn(0); 6320cf1dd8SToby Isaac ierr = PetscFree3((*geom)->v,(*geom)->J,(*geom)->detJ);CHKERRQ(ierr); 6420cf1dd8SToby Isaac ierr = PetscFree((*geom)->invJ);CHKERRQ(ierr); 6541418987SMatthew G. Knepley ierr = PetscFree2((*geom)->face,(*geom)->n);CHKERRQ(ierr); 6641418987SMatthew G. Knepley ierr = PetscFree6((*geom)->suppJ[0],(*geom)->suppJ[1],(*geom)->suppInvJ[0],(*geom)->suppInvJ[1],(*geom)->suppDetJ[0],(*geom)->suppDetJ[1]);CHKERRQ(ierr); 6720cf1dd8SToby Isaac ierr = PetscFree(*geom);CHKERRQ(ierr); 6820cf1dd8SToby Isaac PetscFunctionReturn(0); 6920cf1dd8SToby Isaac } 7020cf1dd8SToby Isaac 71*2b99622eSMatthew G. Knepley /*@C 72*2b99622eSMatthew G. Knepley PetscFEGeomGetChunk - Get a chunk of cells in the group as a PetscFEGeom 73*2b99622eSMatthew G. Knepley 74*2b99622eSMatthew G. Knepley Input Parameters: 75*2b99622eSMatthew G. Knepley + geom - PetscFEGeom object 76*2b99622eSMatthew G. Knepley . cStart - The first cell in the chunk 77*2b99622eSMatthew G. Knepley - cEnd - The first cell not in the chunk 78*2b99622eSMatthew G. Knepley 79*2b99622eSMatthew G. Knepley Output Parameter: 80*2b99622eSMatthew G. Knepley . chunkGeom - The chunk of cells 81*2b99622eSMatthew G. Knepley 82*2b99622eSMatthew G. Knepley Level: intermediate 83*2b99622eSMatthew G. Knepley 84*2b99622eSMatthew G. Knepley .seealso: PetscFEGeomRestoreChunk(), PetscFEGeomCreate() 85*2b99622eSMatthew G. Knepley @*/ 8620cf1dd8SToby Isaac PetscErrorCode PetscFEGeomGetChunk(PetscFEGeom *geom, PetscInt cStart, PetscInt cEnd, PetscFEGeom **chunkGeom) 8720cf1dd8SToby Isaac { 8820cf1dd8SToby Isaac PetscInt Nq; 8920cf1dd8SToby Isaac PetscInt dE; 9020cf1dd8SToby Isaac PetscErrorCode ierr; 9120cf1dd8SToby Isaac 9220cf1dd8SToby Isaac PetscFunctionBegin; 9320cf1dd8SToby Isaac PetscValidPointer(geom,1); 9420cf1dd8SToby Isaac PetscValidPointer(chunkGeom,2); 9520cf1dd8SToby Isaac if (!(*chunkGeom)) { 9620cf1dd8SToby Isaac ierr = PetscNew(chunkGeom);CHKERRQ(ierr); 9720cf1dd8SToby Isaac } 9820cf1dd8SToby Isaac Nq = geom->numPoints; 9920cf1dd8SToby Isaac dE= geom->dimEmbed; 10020cf1dd8SToby Isaac (*chunkGeom)->dim = geom->dim; 10120cf1dd8SToby Isaac (*chunkGeom)->dimEmbed = geom->dimEmbed; 10220cf1dd8SToby Isaac (*chunkGeom)->numPoints = geom->numPoints; 10320cf1dd8SToby Isaac (*chunkGeom)->numCells = cEnd - cStart; 10420cf1dd8SToby Isaac (*chunkGeom)->xi = geom->xi; 10520cf1dd8SToby Isaac (*chunkGeom)->v = &geom->v[Nq*dE*cStart]; 10620cf1dd8SToby Isaac (*chunkGeom)->J = &geom->J[Nq*dE*dE*cStart]; 10720cf1dd8SToby Isaac (*chunkGeom)->invJ = (geom->invJ) ? &geom->invJ[Nq*dE*dE*cStart] : NULL; 10820cf1dd8SToby Isaac (*chunkGeom)->detJ = &geom->detJ[Nq*cStart]; 10920cf1dd8SToby Isaac (*chunkGeom)->n = geom->n ? &geom->n[Nq*dE*cStart] : NULL; 11020cf1dd8SToby Isaac (*chunkGeom)->face = geom->face ? &geom->face[cStart] : NULL; 11141418987SMatthew G. Knepley (*chunkGeom)->suppJ[0] = geom->suppJ[0] ? &geom->suppJ[0][Nq*dE*dE*cStart] : NULL; 11241418987SMatthew G. Knepley (*chunkGeom)->suppJ[1] = geom->suppJ[1] ? &geom->suppJ[1][Nq*dE*dE*cStart] : NULL; 11320cf1dd8SToby Isaac (*chunkGeom)->suppInvJ[0] = geom->suppInvJ[0] ? &geom->suppInvJ[0][Nq*dE*dE*cStart] : NULL; 11420cf1dd8SToby Isaac (*chunkGeom)->suppInvJ[1] = geom->suppInvJ[1] ? &geom->suppInvJ[1][Nq*dE*dE*cStart] : NULL; 11541418987SMatthew G. Knepley (*chunkGeom)->suppDetJ[0] = geom->suppDetJ[0] ? &geom->suppDetJ[0][Nq*cStart] : NULL; 11641418987SMatthew G. Knepley (*chunkGeom)->suppDetJ[1] = geom->suppDetJ[1] ? &geom->suppDetJ[1][Nq*cStart] : NULL; 11720cf1dd8SToby Isaac (*chunkGeom)->isAffine = geom->isAffine; 11820cf1dd8SToby Isaac PetscFunctionReturn(0); 11920cf1dd8SToby Isaac } 12020cf1dd8SToby Isaac 121*2b99622eSMatthew G. Knepley /*@C 122*2b99622eSMatthew G. Knepley PetscFEGeomRestoreChunk - Restore the chunk 123*2b99622eSMatthew G. Knepley 124*2b99622eSMatthew G. Knepley Input Parameters: 125*2b99622eSMatthew G. Knepley + geom - PetscFEGeom object 126*2b99622eSMatthew G. Knepley . cStart - The first cell in the chunk 127*2b99622eSMatthew G. Knepley . cEnd - The first cell not in the chunk 128*2b99622eSMatthew G. Knepley - chunkGeom - The chunk of cells 129*2b99622eSMatthew G. Knepley 130*2b99622eSMatthew G. Knepley Level: intermediate 131*2b99622eSMatthew G. Knepley 132*2b99622eSMatthew G. Knepley .seealso: PetscFEGeomGetChunk(), PetscFEGeomCreate() 133*2b99622eSMatthew G. Knepley @*/ 13420cf1dd8SToby Isaac PetscErrorCode PetscFEGeomRestoreChunk(PetscFEGeom *geom, PetscInt cStart, PetscInt cEnd, PetscFEGeom **chunkGeom) 13520cf1dd8SToby Isaac { 13620cf1dd8SToby Isaac PetscErrorCode ierr; 13720cf1dd8SToby Isaac 13820cf1dd8SToby Isaac PetscFunctionBegin; 13920cf1dd8SToby Isaac ierr = PetscFree(*chunkGeom);CHKERRQ(ierr); 14020cf1dd8SToby Isaac PetscFunctionReturn(0); 14120cf1dd8SToby Isaac } 14220cf1dd8SToby Isaac 143*2b99622eSMatthew G. Knepley /*@ 144*2b99622eSMatthew G. Knepley PetscFEGeomComplete - Calculate derived quntites from base geometry specification 145*2b99622eSMatthew G. Knepley 146*2b99622eSMatthew G. Knepley Input Parameter: 147*2b99622eSMatthew G. Knepley . geom - PetscFEGeom object 148*2b99622eSMatthew G. Knepley 149*2b99622eSMatthew G. Knepley Level: intermediate 150*2b99622eSMatthew G. Knepley 151*2b99622eSMatthew G. Knepley .seealso: PetscFEGeomCreate() 152*2b99622eSMatthew G. Knepley @*/ 15320cf1dd8SToby Isaac PetscErrorCode PetscFEGeomComplete(PetscFEGeom *geom) 15420cf1dd8SToby Isaac { 15520cf1dd8SToby Isaac PetscInt i, j, N, dE; 15620cf1dd8SToby Isaac 15720cf1dd8SToby Isaac PetscFunctionBeginHot; 15820cf1dd8SToby Isaac N = geom->numPoints * geom->numCells; 15920cf1dd8SToby Isaac dE = geom->dimEmbed; 16020cf1dd8SToby Isaac switch (dE) { 16120cf1dd8SToby Isaac case 3: 16220cf1dd8SToby Isaac for (i = 0; i < N; i++) { 16320cf1dd8SToby Isaac DMPlex_Det3D_Internal(&geom->detJ[i], &geom->J[dE*dE*i]); 16420cf1dd8SToby Isaac if (geom->invJ) DMPlex_Invert3D_Internal(&geom->invJ[dE*dE*i], &geom->J[dE*dE*i], geom->detJ[i]); 16520cf1dd8SToby Isaac } 16620cf1dd8SToby Isaac break; 16720cf1dd8SToby Isaac case 2: 16820cf1dd8SToby Isaac for (i = 0; i < N; i++) { 16920cf1dd8SToby Isaac DMPlex_Det2D_Internal(&geom->detJ[i], &geom->J[dE*dE*i]); 17020cf1dd8SToby Isaac if (geom->invJ) DMPlex_Invert2D_Internal(&geom->invJ[dE*dE*i], &geom->J[dE*dE*i], geom->detJ[i]); 17120cf1dd8SToby Isaac } 17220cf1dd8SToby Isaac break; 17320cf1dd8SToby Isaac case 1: 17420cf1dd8SToby Isaac for (i = 0; i < N; i++) { 17520cf1dd8SToby Isaac geom->detJ[i] = PetscAbsReal(geom->J[i]); 17620cf1dd8SToby Isaac if (geom->invJ) geom->invJ[i] = 1. / geom->J[i]; 17720cf1dd8SToby Isaac } 17820cf1dd8SToby Isaac break; 17920cf1dd8SToby Isaac } 18020cf1dd8SToby Isaac if (geom->n) { 18120cf1dd8SToby Isaac for (i = 0; i < N; i++) { 18220cf1dd8SToby Isaac for (j = 0; j < dE; j++) { 18320cf1dd8SToby Isaac geom->n[dE*i + j] = geom->J[dE*dE*i + dE*j + dE-1] * ((dE == 2) ? -1. : 1.); 18420cf1dd8SToby Isaac } 18520cf1dd8SToby Isaac } 18620cf1dd8SToby Isaac } 18720cf1dd8SToby Isaac PetscFunctionReturn(0); 18820cf1dd8SToby Isaac } 189