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