xref: /petsc/src/ts/utils/dmplexts.c (revision 6dbbd3060fde0bb4b0dff80567db5704b5e5a680)
1*6dbbd306SMatthew G. Knepley #include <petscdmplex.h>          /*I "petscdmplex.h" I*/
2*6dbbd306SMatthew G. Knepley #include <petscfv.h>
3*6dbbd306SMatthew G. Knepley 
4*6dbbd306SMatthew G. Knepley struct _FVPhysics {
5*6dbbd306SMatthew G. Knepley   void (*riemann)(const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], PetscScalar flux[], void *ctx);
6*6dbbd306SMatthew G. Knepley };
7*6dbbd306SMatthew G. Knepley 
8*6dbbd306SMatthew G. Knepley #undef __FUNCT__
9*6dbbd306SMatthew G. Knepley #define __FUNCT__ "TSComputeRHSFunctionLocalDMPlex"
10*6dbbd306SMatthew G. Knepley PetscErrorCode TSComputeRHSFunctionLocalDMPlex(DM dm, PetscReal time, Vec locX, Vec F, void *ctx)
11*6dbbd306SMatthew G. Knepley {
12*6dbbd306SMatthew G. Knepley   void             (*riemann)(const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], PetscScalar flux[], void *ctx) = ((struct _FVPhysics *) ctx)->riemann;
13*6dbbd306SMatthew G. Knepley   PetscFV            fvm;
14*6dbbd306SMatthew G. Knepley   Vec                faceGeometry, cellGeometry;
15*6dbbd306SMatthew G. Knepley   DM                 dmFace, dmCell;
16*6dbbd306SMatthew G. Knepley   DMLabel            ghostLabel;
17*6dbbd306SMatthew G. Knepley   PetscCellGeometry  fgeom, cgeom;
18*6dbbd306SMatthew G. Knepley   const PetscScalar *facegeom, *cellgeom, *x;
19*6dbbd306SMatthew G. Knepley   PetscScalar       *f, *uL, *uR, *fluxL, *fluxR;
20*6dbbd306SMatthew G. Knepley   PetscReal         *centroid, *normal, *vol;
21*6dbbd306SMatthew G. Knepley   PetscInt           Nf, dim, pdim, fStart, fEnd, numFaces = 0, face, iface;
22*6dbbd306SMatthew G. Knepley   PetscErrorCode     ierr;
23*6dbbd306SMatthew G. Knepley 
24*6dbbd306SMatthew G. Knepley   PetscFunctionBeginUser;
25*6dbbd306SMatthew G. Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
26*6dbbd306SMatthew G. Knepley   ierr = DMGetNumFields(dm, &Nf);CHKERRQ(ierr);
27*6dbbd306SMatthew G. Knepley   ierr = DMGetField(dm, 0, (PetscObject *) &fvm);CHKERRQ(ierr);
28*6dbbd306SMatthew G. Knepley   ierr = PetscFVGetNumComponents(fvm, &pdim);CHKERRQ(ierr);
29*6dbbd306SMatthew G. Knepley   ierr = DMPlexInsertBoundaryValuesFVM(dm, time, locX);CHKERRQ(ierr);
30*6dbbd306SMatthew G. Knepley   ierr = DMPlexGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr);
31*6dbbd306SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr);
32*6dbbd306SMatthew G. Knepley   /* TODO Pull this geometry calculation into the library */
33*6dbbd306SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "FaceGeometry", (PetscObject *) &faceGeometry);CHKERRQ(ierr);
34*6dbbd306SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "CellGeometry", (PetscObject *) &cellGeometry);CHKERRQ(ierr);
35*6dbbd306SMatthew G. Knepley   ierr = VecGetDM(faceGeometry, &dmFace);CHKERRQ(ierr);
36*6dbbd306SMatthew G. Knepley   ierr = VecGetDM(cellGeometry, &dmCell);CHKERRQ(ierr);
37*6dbbd306SMatthew G. Knepley   ierr = VecGetArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr);
38*6dbbd306SMatthew G. Knepley   ierr = VecGetArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr);
39*6dbbd306SMatthew G. Knepley   ierr = VecGetArrayRead(locX, &x);CHKERRQ(ierr);
40*6dbbd306SMatthew G. Knepley   for (face = fStart; face < fEnd; ++face) {
41*6dbbd306SMatthew G. Knepley     PetscInt ghost;
42*6dbbd306SMatthew G. Knepley 
43*6dbbd306SMatthew G. Knepley     ierr = DMLabelGetValue(ghostLabel, face, &ghost);CHKERRQ(ierr);
44*6dbbd306SMatthew G. Knepley     if (ghost >= 0) continue;
45*6dbbd306SMatthew G. Knepley     ++numFaces;
46*6dbbd306SMatthew G. Knepley   }
47*6dbbd306SMatthew G. Knepley   ierr = PetscMalloc7(numFaces*dim,&centroid,numFaces*dim,&normal,numFaces*2,&vol,numFaces*pdim,&uL,numFaces*pdim,&uR,numFaces*pdim,&fluxL,numFaces*pdim,&fluxR);CHKERRQ(ierr);
48*6dbbd306SMatthew G. Knepley   for (face = fStart, iface = 0; face < fEnd; ++face) {
49*6dbbd306SMatthew G. Knepley     const PetscInt    *cells;
50*6dbbd306SMatthew G. Knepley     const FaceGeom    *fg;
51*6dbbd306SMatthew G. Knepley     const CellGeom    *cgL, *cgR;
52*6dbbd306SMatthew G. Knepley     const PetscScalar *xL, *xR;
53*6dbbd306SMatthew G. Knepley     PetscInt           ghost, d;
54*6dbbd306SMatthew G. Knepley 
55*6dbbd306SMatthew G. Knepley     ierr = DMLabelGetValue(ghostLabel, face, &ghost);CHKERRQ(ierr);
56*6dbbd306SMatthew G. Knepley     if (ghost >= 0) continue;
57*6dbbd306SMatthew G. Knepley     ierr = DMPlexPointLocalRead(dmFace, face, facegeom, &fg);CHKERRQ(ierr);
58*6dbbd306SMatthew G. Knepley     ierr = DMPlexGetSupport(dm, face, &cells);CHKERRQ(ierr);
59*6dbbd306SMatthew G. Knepley     ierr = DMPlexPointLocalRead(dmCell, cells[0], cellgeom, &cgL);CHKERRQ(ierr);
60*6dbbd306SMatthew G. Knepley     ierr = DMPlexPointLocalRead(dmCell, cells[1], cellgeom, &cgR);CHKERRQ(ierr);
61*6dbbd306SMatthew G. Knepley     ierr = DMPlexPointLocalRead(dm, cells[0], x, &xL);CHKERRQ(ierr);
62*6dbbd306SMatthew G. Knepley     ierr = DMPlexPointLocalRead(dm, cells[1], x, &xR);CHKERRQ(ierr);
63*6dbbd306SMatthew G. Knepley     for (d = 0; d < pdim; ++d) {
64*6dbbd306SMatthew G. Knepley       uL[iface*pdim+d] = xL[d];
65*6dbbd306SMatthew G. Knepley       uR[iface*pdim+d] = xR[d];
66*6dbbd306SMatthew G. Knepley     }
67*6dbbd306SMatthew G. Knepley     for (d = 0; d < dim; ++d) {
68*6dbbd306SMatthew G. Knepley       centroid[iface*dim+d] = fg->centroid[d];
69*6dbbd306SMatthew G. Knepley       normal[iface*dim+d]   = fg->normal[d];
70*6dbbd306SMatthew G. Knepley     }
71*6dbbd306SMatthew G. Knepley     vol[iface*2+0] = cgL->volume;
72*6dbbd306SMatthew G. Knepley     vol[iface*2+1] = cgR->volume;
73*6dbbd306SMatthew G. Knepley     ++iface;
74*6dbbd306SMatthew G. Knepley   }
75*6dbbd306SMatthew G. Knepley   ierr = VecRestoreArrayRead(locX, &x);CHKERRQ(ierr);
76*6dbbd306SMatthew G. Knepley   ierr = VecRestoreArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr);
77*6dbbd306SMatthew G. Knepley   ierr = VecRestoreArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr);
78*6dbbd306SMatthew G. Knepley   fgeom.v0  = centroid;
79*6dbbd306SMatthew G. Knepley   fgeom.n   = normal;
80*6dbbd306SMatthew G. Knepley   cgeom.vol = vol;
81*6dbbd306SMatthew G. Knepley   ierr = PetscFVIntegrateRHSFunction(fvm, numFaces, Nf, &fvm, 0, fgeom, cgeom, uL, uR, riemann, fluxL, fluxR, ctx);CHKERRQ(ierr);
82*6dbbd306SMatthew G. Knepley   ierr = VecGetArray(F, &f);CHKERRQ(ierr);
83*6dbbd306SMatthew G. Knepley   for (face = fStart, iface = 0; face < fEnd; ++face) {
84*6dbbd306SMatthew G. Knepley     const PetscInt *cells;
85*6dbbd306SMatthew G. Knepley     PetscScalar    *fL, *fR;
86*6dbbd306SMatthew G. Knepley     PetscInt        ghost, d;
87*6dbbd306SMatthew G. Knepley 
88*6dbbd306SMatthew G. Knepley     ierr = DMLabelGetValue(ghostLabel, face, &ghost);CHKERRQ(ierr);
89*6dbbd306SMatthew G. Knepley     if (ghost >= 0) continue;
90*6dbbd306SMatthew G. Knepley     ierr = DMPlexGetSupport(dm, face, &cells);CHKERRQ(ierr);
91*6dbbd306SMatthew G. Knepley     ierr = DMPlexPointGlobalRef(dm, cells[0], f, &fL);CHKERRQ(ierr);
92*6dbbd306SMatthew G. Knepley     ierr = DMPlexPointGlobalRef(dm, cells[1], f, &fR);CHKERRQ(ierr);
93*6dbbd306SMatthew G. Knepley     for (d = 0; d < pdim; ++d) {
94*6dbbd306SMatthew G. Knepley       if (fL) fL[d] -= fluxL[iface*pdim+d];
95*6dbbd306SMatthew G. Knepley       if (fR) fR[d] += fluxR[iface*pdim+d];
96*6dbbd306SMatthew G. Knepley     }
97*6dbbd306SMatthew G. Knepley     ++iface;
98*6dbbd306SMatthew G. Knepley   }
99*6dbbd306SMatthew G. Knepley   ierr = VecRestoreArray(F, &f);CHKERRQ(ierr);
100*6dbbd306SMatthew G. Knepley   ierr = PetscFree7(centroid,normal,vol,uL,uR,fluxL,fluxR);CHKERRQ(ierr);
101*6dbbd306SMatthew G. Knepley   PetscFunctionReturn(0);
102*6dbbd306SMatthew G. Knepley }
103