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