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