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