xref: /petsc/src/dm/impls/plex/plexfvm.c (revision 5520554388890bd89a1c1cf7870aedf4e71d512f)
1 #include <petsc/private/dmpleximpl.h> /*I      "petscdmplex.h"   I*/
2 #include <petscsf.h>
3 
4 #include <petsc/private/petscfeimpl.h>
5 #include <petsc/private/petscfvimpl.h>
6 
7 static PetscErrorCode DMPlexApplyLimiter_Internal(DM dm, DM dmCell, PetscLimiter lim, PetscInt dim, PetscInt dof, PetscInt cell, PetscInt field, PetscInt face, PetscInt fStart, PetscInt fEnd, PetscReal *cellPhi, const PetscScalar *x, const PetscScalar *cellgeom, const PetscFVCellGeom *cg, const PetscScalar *cx, const PetscScalar *cgrad) {
8   const PetscInt *children;
9   PetscInt        numChildren;
10 
11   PetscFunctionBegin;
12   PetscCall(DMPlexGetTreeChildren(dm, face, &numChildren, &children));
13   if (numChildren) {
14     PetscInt c;
15 
16     for (c = 0; c < numChildren; c++) {
17       PetscInt childFace = children[c];
18 
19       if (childFace >= fStart && childFace < fEnd) PetscCall(DMPlexApplyLimiter_Internal(dm, dmCell, lim, dim, dof, cell, field, childFace, fStart, fEnd, cellPhi, x, cellgeom, cg, cx, cgrad));
20     }
21   } else {
22     PetscScalar     *ncx;
23     PetscFVCellGeom *ncg;
24     const PetscInt  *fcells;
25     PetscInt         ncell, d;
26     PetscReal        v[3];
27 
28     PetscCall(DMPlexGetSupport(dm, face, &fcells));
29     ncell = cell == fcells[0] ? fcells[1] : fcells[0];
30     if (field >= 0) {
31       PetscCall(DMPlexPointLocalFieldRead(dm, ncell, field, x, &ncx));
32     } else {
33       PetscCall(DMPlexPointLocalRead(dm, ncell, x, &ncx));
34     }
35     PetscCall(DMPlexPointLocalRead(dmCell, ncell, cellgeom, &ncg));
36     DMPlex_WaxpyD_Internal(dim, -1, cg->centroid, ncg->centroid, v);
37     for (d = 0; d < dof; ++d) {
38       /* We use the symmetric slope limited form of Berger, Aftosmis, and Murman 2005 */
39       PetscReal denom = DMPlex_DotD_Internal(dim, &cgrad[d * dim], v);
40       PetscReal phi, flim = 0.5 * PetscRealPart(ncx[d] - cx[d]) / denom;
41 
42       PetscCall(PetscLimiterLimit(lim, flim, &phi));
43       cellPhi[d] = PetscMin(cellPhi[d], phi);
44     }
45   }
46   PetscFunctionReturn(0);
47 }
48 
49 PetscErrorCode DMPlexReconstructGradients_Internal(DM dm, PetscFV fvm, PetscInt fStart, PetscInt fEnd, Vec faceGeometry, Vec cellGeometry, Vec locX, Vec grad) {
50   DM                 dmFace, dmCell, dmGrad;
51   DMLabel            ghostLabel;
52   PetscDS            prob;
53   PetscLimiter       lim;
54   const PetscScalar *facegeom, *cellgeom, *x;
55   PetscScalar       *gr;
56   PetscReal         *cellPhi;
57   PetscInt           dim, face, cell, field, dof, cStart, cEnd, nFields;
58 
59   PetscFunctionBegin;
60   PetscCall(DMGetDimension(dm, &dim));
61   PetscCall(DMGetDS(dm, &prob));
62   PetscCall(PetscDSGetNumFields(prob, &nFields));
63   PetscCall(PetscDSGetFieldIndex(prob, (PetscObject)fvm, &field));
64   PetscCall(PetscDSGetFieldSize(prob, field, &dof));
65   PetscCall(DMGetLabel(dm, "ghost", &ghostLabel));
66   PetscCall(PetscFVGetLimiter(fvm, &lim));
67   PetscCall(VecGetDM(faceGeometry, &dmFace));
68   PetscCall(VecGetArrayRead(faceGeometry, &facegeom));
69   PetscCall(VecGetDM(cellGeometry, &dmCell));
70   PetscCall(VecGetArrayRead(cellGeometry, &cellgeom));
71   PetscCall(VecGetArrayRead(locX, &x));
72   PetscCall(VecGetDM(grad, &dmGrad));
73   PetscCall(VecZeroEntries(grad));
74   PetscCall(VecGetArray(grad, &gr));
75   /* Reconstruct gradients */
76   for (face = fStart; face < fEnd; ++face) {
77     const PetscInt  *cells;
78     PetscFVFaceGeom *fg;
79     PetscScalar     *cx[2];
80     PetscScalar     *cgrad[2];
81     PetscBool        boundary;
82     PetscInt         ghost, c, pd, d, numChildren, numCells;
83 
84     PetscCall(DMLabelGetValue(ghostLabel, face, &ghost));
85     PetscCall(DMIsBoundaryPoint(dm, face, &boundary));
86     PetscCall(DMPlexGetTreeChildren(dm, face, &numChildren, NULL));
87     if (ghost >= 0 || boundary || numChildren) continue;
88     PetscCall(DMPlexGetSupportSize(dm, face, &numCells));
89     PetscCheck(numCells == 2, PETSC_COMM_SELF, PETSC_ERR_PLIB, "facet %" PetscInt_FMT " has %" PetscInt_FMT " support points: expected 2", face, numCells);
90     PetscCall(DMPlexGetSupport(dm, face, &cells));
91     PetscCall(DMPlexPointLocalRead(dmFace, face, facegeom, &fg));
92     for (c = 0; c < 2; ++c) {
93       if (nFields > 1) {
94         PetscCall(DMPlexPointLocalFieldRead(dm, cells[c], field, x, &cx[c]));
95       } else {
96         PetscCall(DMPlexPointLocalRead(dm, cells[c], x, &cx[c]));
97       }
98       PetscCall(DMPlexPointGlobalRef(dmGrad, cells[c], gr, &cgrad[c]));
99     }
100     for (pd = 0; pd < dof; ++pd) {
101       PetscScalar delta = cx[1][pd] - cx[0][pd];
102 
103       for (d = 0; d < dim; ++d) {
104         if (cgrad[0]) cgrad[0][pd * dim + d] += fg->grad[0][d] * delta;
105         if (cgrad[1]) cgrad[1][pd * dim + d] -= fg->grad[1][d] * delta;
106       }
107     }
108   }
109   /* Limit interior gradients (using cell-based loop because it generalizes better to vector limiters) */
110   PetscCall(DMPlexGetSimplexOrBoxCells(dm, 0, &cStart, &cEnd));
111   PetscCall(DMGetWorkArray(dm, dof, MPIU_REAL, &cellPhi));
112   for (cell = (dmGrad && lim) ? cStart : cEnd; cell < cEnd; ++cell) {
113     const PetscInt  *faces;
114     PetscScalar     *cx;
115     PetscFVCellGeom *cg;
116     PetscScalar     *cgrad;
117     PetscInt         coneSize, f, pd, d;
118 
119     PetscCall(DMPlexGetConeSize(dm, cell, &coneSize));
120     PetscCall(DMPlexGetCone(dm, cell, &faces));
121     if (nFields > 1) {
122       PetscCall(DMPlexPointLocalFieldRead(dm, cell, field, x, &cx));
123     } else {
124       PetscCall(DMPlexPointLocalRead(dm, cell, x, &cx));
125     }
126     PetscCall(DMPlexPointLocalRead(dmCell, cell, cellgeom, &cg));
127     PetscCall(DMPlexPointGlobalRef(dmGrad, cell, gr, &cgrad));
128     if (!cgrad) continue; /* Unowned overlap cell, we do not compute */
129     /* Limiter will be minimum value over all neighbors */
130     for (d = 0; d < dof; ++d) cellPhi[d] = PETSC_MAX_REAL;
131     for (f = 0; f < coneSize; ++f) PetscCall(DMPlexApplyLimiter_Internal(dm, dmCell, lim, dim, dof, cell, nFields > 1 ? field : -1, faces[f], fStart, fEnd, cellPhi, x, cellgeom, cg, cx, cgrad));
132     /* Apply limiter to gradient */
133     for (pd = 0; pd < dof; ++pd) /* Scalar limiter applied to each component separately */
134       for (d = 0; d < dim; ++d) cgrad[pd * dim + d] *= cellPhi[pd];
135   }
136   PetscCall(DMRestoreWorkArray(dm, dof, MPIU_REAL, &cellPhi));
137   PetscCall(VecRestoreArrayRead(faceGeometry, &facegeom));
138   PetscCall(VecRestoreArrayRead(cellGeometry, &cellgeom));
139   PetscCall(VecRestoreArrayRead(locX, &x));
140   PetscCall(VecRestoreArray(grad, &gr));
141   PetscFunctionReturn(0);
142 }
143 
144 /*@
145   DMPlexReconstructGradientsFVM - reconstruct the gradient of a vector using a finite volume method.
146 
147   Input Parameters:
148 + dm - the mesh
149 - locX - the local representation of the vector
150 
151   Output Parameter:
152 . grad - the global representation of the gradient
153 
154   Level: developer
155 
156 .seealso: `DMPlexGetGradientDM()`
157 @*/
158 PetscErrorCode DMPlexReconstructGradientsFVM(DM dm, Vec locX, Vec grad) {
159   PetscDS          prob;
160   PetscInt         Nf, f, fStart, fEnd;
161   PetscBool        useFVM = PETSC_FALSE;
162   PetscFV          fvm    = NULL;
163   Vec              faceGeometryFVM, cellGeometryFVM;
164   PetscFVCellGeom *cgeomFVM = NULL;
165   PetscFVFaceGeom *fgeomFVM = NULL;
166   DM               dmGrad   = NULL;
167 
168   PetscFunctionBegin;
169   PetscCall(DMGetDS(dm, &prob));
170   PetscCall(PetscDSGetNumFields(prob, &Nf));
171   for (f = 0; f < Nf; ++f) {
172     PetscObject  obj;
173     PetscClassId id;
174 
175     PetscCall(PetscDSGetDiscretization(prob, f, &obj));
176     PetscCall(PetscObjectGetClassId(obj, &id));
177     if (id == PETSCFV_CLASSID) {
178       useFVM = PETSC_TRUE;
179       fvm    = (PetscFV)obj;
180     }
181   }
182   PetscCheck(useFVM, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "This dm does not have a finite volume discretization");
183   PetscCall(DMPlexGetDataFVM(dm, fvm, &cellGeometryFVM, &faceGeometryFVM, &dmGrad));
184   PetscCheck(dmGrad, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "This dm's finite volume discretization does not reconstruct gradients");
185   PetscCall(VecGetArrayRead(faceGeometryFVM, (const PetscScalar **)&fgeomFVM));
186   PetscCall(VecGetArrayRead(cellGeometryFVM, (const PetscScalar **)&cgeomFVM));
187   PetscCall(DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd));
188   PetscCall(DMPlexReconstructGradients_Internal(dm, fvm, fStart, fEnd, faceGeometryFVM, cellGeometryFVM, locX, grad));
189   PetscFunctionReturn(0);
190 }
191