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