11faf85eaSMatthew G. Knepley #include <petsc-private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2254c1ad2SMatthew G. Knepley #include <petsc-private/tsimpl.h> /*I "petscts.h" I*/ 3924a1b8fSMatthew G. Knepley #include <petscds.h> 46dbbd306SMatthew G. Knepley #include <petscfv.h> 56dbbd306SMatthew G. Knepley 6254c1ad2SMatthew G. Knepley #undef __FUNCT__ 7a0ac79e7SMatthew G. Knepley #define __FUNCT__ "DMPlexTSGetGeometry" 8c510411aSMatthew G. Knepley /*@ 9c510411aSMatthew G. Knepley DMPlexTSGetGeometry - Return precomputed geometric data 10c510411aSMatthew G. Knepley 11c510411aSMatthew G. Knepley Input Parameter: 12c510411aSMatthew G. Knepley . dm - The DM 13c510411aSMatthew G. Knepley 14c510411aSMatthew G. Knepley Output Parameters: 15c510411aSMatthew G. Knepley + facegeom - The values precomputed from face geometry 16c510411aSMatthew G. Knepley . cellgeom - The values precomputed from cell geometry 17c510411aSMatthew G. Knepley - minRadius - The minimum radius over the mesh of an inscribed sphere in a cell 18c510411aSMatthew G. Knepley 19c510411aSMatthew G. Knepley Level: developer 20c510411aSMatthew G. Knepley 21c510411aSMatthew G. Knepley .seealso: DMPlexTSSetRHSFunctionLocal() 22c510411aSMatthew G. Knepley @*/ 23a0ac79e7SMatthew G. Knepley PetscErrorCode DMPlexTSGetGeometry(DM dm, Vec *facegeom, Vec *cellgeom, PetscReal *minRadius) 24a0ac79e7SMatthew G. Knepley { 25a0ac79e7SMatthew G. Knepley DMTS dmts; 261faf85eaSMatthew G. Knepley PetscObject obj; 27a0ac79e7SMatthew G. Knepley PetscErrorCode ierr; 28a0ac79e7SMatthew G. Knepley 29a0ac79e7SMatthew G. Knepley PetscFunctionBegin; 30924a1b8fSMatthew G. Knepley PetscValidHeaderSpecific(dm,DM_CLASSID,1); 31a0ac79e7SMatthew G. Knepley ierr = DMGetDMTS(dm, &dmts);CHKERRQ(ierr); 321faf85eaSMatthew G. Knepley ierr = PetscObjectQuery((PetscObject) dmts, "DMPlexTS_facegeom", &obj);CHKERRQ(ierr); 331faf85eaSMatthew G. Knepley if (!obj) { 341faf85eaSMatthew G. Knepley Vec cellgeom, facegeom; 351faf85eaSMatthew G. Knepley 361faf85eaSMatthew G. Knepley ierr = DMPlexComputeGeometryFVM(dm, &cellgeom, &facegeom);CHKERRQ(ierr); 371faf85eaSMatthew G. Knepley ierr = PetscObjectCompose((PetscObject) dmts, "DMPlexTS_facegeom", (PetscObject) facegeom);CHKERRQ(ierr); 381faf85eaSMatthew G. Knepley ierr = PetscObjectCompose((PetscObject) dmts, "DMPlexTS_cellgeom", (PetscObject) cellgeom);CHKERRQ(ierr); 391faf85eaSMatthew G. Knepley ierr = VecDestroy(&facegeom);CHKERRQ(ierr); 401faf85eaSMatthew G. Knepley ierr = VecDestroy(&cellgeom);CHKERRQ(ierr); 411faf85eaSMatthew G. Knepley } 42924a1b8fSMatthew G. Knepley if (facegeom) {PetscValidPointer(facegeom, 2); ierr = PetscObjectQuery((PetscObject) dmts, "DMPlexTS_facegeom", (PetscObject *) facegeom);CHKERRQ(ierr);} 43924a1b8fSMatthew G. Knepley if (cellgeom) {PetscValidPointer(cellgeom, 3); ierr = PetscObjectQuery((PetscObject) dmts, "DMPlexTS_cellgeom", (PetscObject *) cellgeom);CHKERRQ(ierr);} 44113c68e6SMatthew G. Knepley if (minRadius) {ierr = DMPlexGetMinRadius(dm, minRadius);CHKERRQ(ierr);} 45924a1b8fSMatthew G. Knepley PetscFunctionReturn(0); 46924a1b8fSMatthew G. Knepley } 47924a1b8fSMatthew G. Knepley 48924a1b8fSMatthew G. Knepley #undef __FUNCT__ 49924a1b8fSMatthew G. Knepley #define __FUNCT__ "DMPlexTSGetGradientDM" 50924a1b8fSMatthew G. Knepley /*@ 51924a1b8fSMatthew G. Knepley DMPlexTSGetGradientDM - Return gradient data layout 52924a1b8fSMatthew G. Knepley 53924a1b8fSMatthew G. Knepley Input Parameter: 54924a1b8fSMatthew G. Knepley . dm - The DM 55924a1b8fSMatthew G. Knepley 56924a1b8fSMatthew G. Knepley Output Parameter: 57924a1b8fSMatthew G. Knepley . dmGrad - The layout for gradient values 58924a1b8fSMatthew G. Knepley 59924a1b8fSMatthew G. Knepley Level: developer 60924a1b8fSMatthew G. Knepley 61924a1b8fSMatthew G. Knepley .seealso: DMPlexTSGetGeometry(), DMPlexTSSetRHSFunctionLocal() 62924a1b8fSMatthew G. Knepley @*/ 63924a1b8fSMatthew G. Knepley PetscErrorCode DMPlexTSGetGradientDM(DM dm, DM *dmGrad) 64924a1b8fSMatthew G. Knepley { 65924a1b8fSMatthew G. Knepley DMTS dmts; 66924a1b8fSMatthew G. Knepley PetscErrorCode ierr; 67924a1b8fSMatthew G. Knepley 68924a1b8fSMatthew G. Knepley PetscFunctionBegin; 69924a1b8fSMatthew G. Knepley PetscValidHeaderSpecific(dm,DM_CLASSID,1); 70924a1b8fSMatthew G. Knepley PetscValidPointer(dmGrad,2); 71924a1b8fSMatthew G. Knepley ierr = DMGetDMTS(dm, &dmts);CHKERRQ(ierr); 72924a1b8fSMatthew G. Knepley ierr = PetscObjectQuery((PetscObject) dmts, "DMPlexTS_dmgrad", (PetscObject *) dmGrad);CHKERRQ(ierr); 73a0ac79e7SMatthew G. Knepley PetscFunctionReturn(0); 74a0ac79e7SMatthew G. Knepley } 75a0ac79e7SMatthew G. Knepley 76a0ac79e7SMatthew G. Knepley #undef __FUNCT__ 77c5148223SMatthew G. Knepley #define __FUNCT__ "BuildGradientReconstruction" 78c5148223SMatthew G. Knepley static PetscErrorCode BuildGradientReconstruction(DM dm, PetscFV fvm, DM dmFace, PetscScalar *fgeom, DM dmCell, PetscScalar *cgeom) 79254c1ad2SMatthew G. Knepley { 808cd7fcbbSMatthew G. Knepley DMLabel ghostLabel; 81c5148223SMatthew G. Knepley PetscScalar *dx, *grad, **gref; 82c5148223SMatthew G. Knepley PetscInt dim, cStart, cEnd, c, cEndInterior, maxNumFaces; 83254c1ad2SMatthew G. Knepley PetscErrorCode ierr; 84254c1ad2SMatthew G. Knepley 85254c1ad2SMatthew G. Knepley PetscFunctionBegin; 86c73cfb54SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 87254c1ad2SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 88c5148223SMatthew G. Knepley ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 89254c1ad2SMatthew G. Knepley ierr = DMPlexGetMaxSizes(dm, &maxNumFaces, NULL);CHKERRQ(ierr); 90c5148223SMatthew G. Knepley ierr = PetscFVLeastSquaresSetMaxFaces(fvm, maxNumFaces);CHKERRQ(ierr); 91254c1ad2SMatthew G. Knepley ierr = DMPlexGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr); 92c5148223SMatthew G. Knepley ierr = PetscMalloc3(maxNumFaces*dim, &dx, maxNumFaces*dim, &grad, maxNumFaces, &gref);CHKERRQ(ierr); 93254c1ad2SMatthew G. Knepley for (c = cStart; c < cEndInterior; c++) { 94254c1ad2SMatthew G. Knepley const PetscInt *faces; 95c5148223SMatthew G. Knepley PetscInt numFaces, usedFaces, f, d; 961faf85eaSMatthew G. Knepley const PetscFVCellGeom *cg; 978cd7fcbbSMatthew G. Knepley PetscBool boundary; 988cd7fcbbSMatthew G. Knepley PetscInt ghost; 998cd7fcbbSMatthew G. Knepley 100254c1ad2SMatthew G. Knepley ierr = DMPlexPointLocalRead(dmCell, c, cgeom, &cg);CHKERRQ(ierr); 101c5148223SMatthew G. Knepley ierr = DMPlexGetConeSize(dm, c, &numFaces);CHKERRQ(ierr); 102c5148223SMatthew G. Knepley ierr = DMPlexGetCone(dm, c, &faces);CHKERRQ(ierr); 103c5148223SMatthew G. Knepley if (numFaces < dim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Cell %D has only %D faces, not enough for gradient reconstruction", c, numFaces); 104c5148223SMatthew G. Knepley for (f = 0, usedFaces = 0; f < numFaces; ++f) { 1051faf85eaSMatthew G. Knepley const PetscFVCellGeom *cg1; 1061faf85eaSMatthew G. Knepley PetscFVFaceGeom *fg; 107254c1ad2SMatthew G. Knepley const PetscInt *fcells; 108254c1ad2SMatthew G. Knepley PetscInt ncell, side; 109254c1ad2SMatthew G. Knepley 110254c1ad2SMatthew G. Knepley ierr = DMLabelGetValue(ghostLabel, faces[f], &ghost);CHKERRQ(ierr); 1118cd7fcbbSMatthew G. Knepley ierr = DMPlexIsBoundaryPoint(dm, faces[f], &boundary);CHKERRQ(ierr); 1128cd7fcbbSMatthew G. Knepley if ((ghost >= 0) || boundary) continue; 113254c1ad2SMatthew G. Knepley ierr = DMPlexGetSupport(dm, faces[f], &fcells);CHKERRQ(ierr); 114254c1ad2SMatthew G. Knepley side = (c != fcells[0]); /* c is on left=0 or right=1 of face */ 115254c1ad2SMatthew G. Knepley ncell = fcells[!side]; /* the neighbor */ 116254c1ad2SMatthew G. Knepley ierr = DMPlexPointLocalRef(dmFace, faces[f], fgeom, &fg);CHKERRQ(ierr); 117254c1ad2SMatthew G. Knepley ierr = DMPlexPointLocalRead(dmCell, ncell, cgeom, &cg1);CHKERRQ(ierr); 118c5148223SMatthew G. Knepley for (d = 0; d < dim; ++d) dx[usedFaces*dim+d] = cg1->centroid[d] - cg->centroid[d]; 119254c1ad2SMatthew G. Knepley gref[usedFaces++] = fg->grad[side]; /* Gradient reconstruction term will go here */ 120254c1ad2SMatthew G. Knepley } 121254c1ad2SMatthew G. Knepley if (!usedFaces) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_USER, "Mesh contains isolated cell (no neighbors). Is it intentional?"); 122c5148223SMatthew G. Knepley ierr = PetscFVComputeGradient(fvm, usedFaces, dx, grad);CHKERRQ(ierr); 123c5148223SMatthew G. Knepley for (f = 0, usedFaces = 0; f < numFaces; ++f) { 124254c1ad2SMatthew G. Knepley ierr = DMLabelGetValue(ghostLabel, faces[f], &ghost);CHKERRQ(ierr); 1258cd7fcbbSMatthew G. Knepley ierr = DMPlexIsBoundaryPoint(dm, faces[f], &boundary);CHKERRQ(ierr); 1268cd7fcbbSMatthew G. Knepley if ((ghost >= 0) || boundary) continue; 127c5148223SMatthew G. Knepley for (d = 0; d < dim; ++d) gref[usedFaces][d] = grad[usedFaces*dim+d]; 128c5148223SMatthew G. Knepley ++usedFaces; 129254c1ad2SMatthew G. Knepley } 130254c1ad2SMatthew G. Knepley } 131c5148223SMatthew G. Knepley ierr = PetscFree3(dx, grad, gref);CHKERRQ(ierr); 132254c1ad2SMatthew G. Knepley PetscFunctionReturn(0); 133254c1ad2SMatthew G. Knepley } 134254c1ad2SMatthew G. Knepley 135254c1ad2SMatthew G. Knepley #undef __FUNCT__ 136254c1ad2SMatthew G. Knepley #define __FUNCT__ "DMPlexTSSetupGradient" 137924a1b8fSMatthew G. Knepley static PetscErrorCode DMPlexTSSetupGradient(DM dm, PetscFV fvm, DMTS dmts) 138254c1ad2SMatthew G. Knepley { 139924a1b8fSMatthew G. Knepley DM dmFace, dmCell, dmGrad; 140924a1b8fSMatthew G. Knepley Vec facegeom, cellgeom; 141254c1ad2SMatthew G. Knepley PetscScalar *fgeom, *cgeom; 142254c1ad2SMatthew G. Knepley PetscSection sectionGrad; 143254c1ad2SMatthew G. Knepley PetscInt dim, pdim, cStart, cEnd, cEndInterior, c; 144254c1ad2SMatthew G. Knepley PetscErrorCode ierr; 145254c1ad2SMatthew G. Knepley 146254c1ad2SMatthew G. Knepley PetscFunctionBegin; 147c73cfb54SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 148254c1ad2SMatthew G. Knepley ierr = PetscFVGetNumComponents(fvm, &pdim);CHKERRQ(ierr); 149254c1ad2SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 150254c1ad2SMatthew G. Knepley ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 151254c1ad2SMatthew G. Knepley /* Construct the interpolant corresponding to each face from the leat-square solution over the cell neighborhood */ 152924a1b8fSMatthew G. Knepley ierr = DMPlexTSGetGeometry(dm, &facegeom, &cellgeom, NULL);CHKERRQ(ierr); 153924a1b8fSMatthew G. Knepley ierr = VecGetDM(facegeom, &dmFace);CHKERRQ(ierr); 154924a1b8fSMatthew G. Knepley ierr = VecGetDM(cellgeom, &dmCell);CHKERRQ(ierr); 155924a1b8fSMatthew G. Knepley ierr = VecGetArray(facegeom, &fgeom);CHKERRQ(ierr); 156924a1b8fSMatthew G. Knepley ierr = VecGetArray(cellgeom, &cgeom);CHKERRQ(ierr); 157c5148223SMatthew G. Knepley ierr = BuildGradientReconstruction(dm, fvm, dmFace, fgeom, dmCell, cgeom);CHKERRQ(ierr); 158924a1b8fSMatthew G. Knepley ierr = VecRestoreArray(facegeom, &fgeom);CHKERRQ(ierr); 159924a1b8fSMatthew G. Knepley ierr = VecRestoreArray(cellgeom, &cgeom);CHKERRQ(ierr); 160254c1ad2SMatthew G. Knepley /* Create storage for gradients */ 161924a1b8fSMatthew G. Knepley ierr = DMClone(dm, &dmGrad);CHKERRQ(ierr); 162254c1ad2SMatthew G. Knepley ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), §ionGrad);CHKERRQ(ierr); 163254c1ad2SMatthew G. Knepley ierr = PetscSectionSetChart(sectionGrad, cStart, cEnd);CHKERRQ(ierr); 164254c1ad2SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) {ierr = PetscSectionSetDof(sectionGrad, c, pdim*dim);CHKERRQ(ierr);} 165254c1ad2SMatthew G. Knepley ierr = PetscSectionSetUp(sectionGrad);CHKERRQ(ierr); 166924a1b8fSMatthew G. Knepley ierr = DMSetDefaultSection(dmGrad, sectionGrad);CHKERRQ(ierr); 167254c1ad2SMatthew G. Knepley ierr = PetscSectionDestroy(§ionGrad);CHKERRQ(ierr); 168924a1b8fSMatthew G. Knepley ierr = PetscObjectCompose((PetscObject) dmts, "DMPlexTS_dmgrad", (PetscObject) dmGrad);CHKERRQ(ierr); 169924a1b8fSMatthew G. Knepley ierr = DMDestroy(&dmGrad);CHKERRQ(ierr); 170254c1ad2SMatthew G. Knepley PetscFunctionReturn(0); 171254c1ad2SMatthew G. Knepley } 172254c1ad2SMatthew G. Knepley 173254c1ad2SMatthew G. Knepley #undef __FUNCT__ 174254c1ad2SMatthew G. Knepley #define __FUNCT__ "DMPlexInsertBoundaryValuesFVM_Static" 175924a1b8fSMatthew G. Knepley static PetscErrorCode DMPlexInsertBoundaryValuesFVM_Static(DM dm, PetscFV fvm, PetscReal time, Vec locX, Vec Grad) 176254c1ad2SMatthew G. Knepley { 177924a1b8fSMatthew G. Knepley Vec faceGeometry, cellGeometry; 178924a1b8fSMatthew G. Knepley DM dmFace, dmCell, dmGrad; 179254c1ad2SMatthew G. Knepley const PetscScalar *facegeom, *cellgeom, *grad; 180254c1ad2SMatthew G. Knepley PetscScalar *x, *fx; 181254c1ad2SMatthew G. Knepley PetscInt numBd, b, dim, pdim, fStart, fEnd; 182254c1ad2SMatthew G. Knepley PetscErrorCode ierr; 183254c1ad2SMatthew G. Knepley 184254c1ad2SMatthew G. Knepley PetscFunctionBegin; 185c73cfb54SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 186924a1b8fSMatthew G. Knepley ierr = DMPlexTSGetGeometry(dm, &faceGeometry, &cellGeometry, NULL);CHKERRQ(ierr); 187924a1b8fSMatthew G. Knepley ierr = DMPlexTSGetGradientDM(dm, &dmGrad);CHKERRQ(ierr); 188254c1ad2SMatthew G. Knepley ierr = PetscFVGetNumComponents(fvm, &pdim);CHKERRQ(ierr); 189254c1ad2SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 190254c1ad2SMatthew G. Knepley ierr = DMPlexGetNumBoundary(dm, &numBd);CHKERRQ(ierr); 191254c1ad2SMatthew G. Knepley if (Grad) { 192254c1ad2SMatthew G. Knepley ierr = VecGetDM(cellGeometry, &dmCell);CHKERRQ(ierr); 193254c1ad2SMatthew G. Knepley ierr = VecGetArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr); 194254c1ad2SMatthew G. Knepley ierr = DMGetWorkArray(dm, pdim, PETSC_SCALAR, &fx);CHKERRQ(ierr); 195254c1ad2SMatthew G. Knepley ierr = VecGetArrayRead(Grad, &grad);CHKERRQ(ierr); 196254c1ad2SMatthew G. Knepley } 197254c1ad2SMatthew G. Knepley ierr = VecGetDM(faceGeometry, &dmFace);CHKERRQ(ierr); 198254c1ad2SMatthew G. Knepley ierr = VecGetArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr); 199254c1ad2SMatthew G. Knepley ierr = VecGetArray(locX, &x);CHKERRQ(ierr); 200254c1ad2SMatthew G. Knepley for (b = 0; b < numBd; ++b) { 2011475bf87SMatthew G. Knepley PetscErrorCode (*func)(PetscReal,const PetscReal*,const PetscReal*,const PetscScalar*,PetscScalar*,void*); 2028cd7fcbbSMatthew G. Knepley DMLabel label; 2038cd7fcbbSMatthew G. Knepley const char *labelname; 204254c1ad2SMatthew G. Knepley const PetscInt *ids; 205254c1ad2SMatthew G. Knepley PetscInt numids, i; 206254c1ad2SMatthew G. Knepley void *ctx; 207254c1ad2SMatthew G. Knepley 2088cd7fcbbSMatthew G. Knepley ierr = DMPlexGetBoundary(dm, b, NULL, NULL, &labelname, NULL, (void (**)()) &func, &numids, &ids, &ctx);CHKERRQ(ierr); 2098cd7fcbbSMatthew G. Knepley ierr = DMPlexGetLabel(dm, labelname, &label);CHKERRQ(ierr); 210254c1ad2SMatthew G. Knepley for (i = 0; i < numids; ++i) { 211254c1ad2SMatthew G. Knepley IS faceIS; 212254c1ad2SMatthew G. Knepley const PetscInt *faces; 213254c1ad2SMatthew G. Knepley PetscInt numFaces, f; 214254c1ad2SMatthew G. Knepley 215254c1ad2SMatthew G. Knepley ierr = DMLabelGetStratumIS(label, ids[i], &faceIS);CHKERRQ(ierr); 216254c1ad2SMatthew G. Knepley if (!faceIS) continue; /* No points with that id on this process */ 217254c1ad2SMatthew G. Knepley ierr = ISGetLocalSize(faceIS, &numFaces);CHKERRQ(ierr); 218254c1ad2SMatthew G. Knepley ierr = ISGetIndices(faceIS, &faces);CHKERRQ(ierr); 219254c1ad2SMatthew G. Knepley for (f = 0; f < numFaces; ++f) { 220254c1ad2SMatthew G. Knepley const PetscInt face = faces[f], *cells; 2211faf85eaSMatthew G. Knepley const PetscFVFaceGeom *fg; 222254c1ad2SMatthew G. Knepley 223254c1ad2SMatthew G. Knepley if ((face < fStart) || (face >= fEnd)) continue; /* Refinement adds non-faces to labels */ 224254c1ad2SMatthew G. Knepley ierr = DMPlexPointLocalRead(dmFace, face, facegeom, &fg);CHKERRQ(ierr); 225254c1ad2SMatthew G. Knepley ierr = DMPlexGetSupport(dm, face, &cells);CHKERRQ(ierr); 226254c1ad2SMatthew G. Knepley if (Grad) { 2271faf85eaSMatthew G. Knepley const PetscFVCellGeom *cg; 228254c1ad2SMatthew G. Knepley const PetscScalar *cx, *cgrad; 2291475bf87SMatthew G. Knepley PetscScalar *xG; 2301475bf87SMatthew G. Knepley PetscReal dx[3]; 231254c1ad2SMatthew G. Knepley PetscInt d; 232254c1ad2SMatthew G. Knepley 233254c1ad2SMatthew G. Knepley ierr = DMPlexPointLocalRead(dmCell, cells[0], cellgeom, &cg);CHKERRQ(ierr); 234254c1ad2SMatthew G. Knepley ierr = DMPlexPointLocalRead(dm, cells[0], x, &cx);CHKERRQ(ierr); 235924a1b8fSMatthew G. Knepley ierr = DMPlexPointLocalRead(dmGrad, cells[0], grad, &cgrad);CHKERRQ(ierr); 236254c1ad2SMatthew G. Knepley ierr = DMPlexPointLocalRef(dm, cells[1], x, &xG);CHKERRQ(ierr); 2371faf85eaSMatthew G. Knepley DMPlex_WaxpyD_Internal(dim, -1, cg->centroid, fg->centroid, dx); 2381faf85eaSMatthew G. Knepley for (d = 0; d < pdim; ++d) fx[d] = cx[d] + DMPlex_DotD_Internal(dim, &cgrad[d*dim], dx); 239254c1ad2SMatthew G. Knepley ierr = (*func)(time, fg->centroid, fg->normal, fx, xG, ctx);CHKERRQ(ierr); 240254c1ad2SMatthew G. Knepley } else { 241254c1ad2SMatthew G. Knepley const PetscScalar *xI; 242254c1ad2SMatthew G. Knepley PetscScalar *xG; 243254c1ad2SMatthew G. Knepley 244254c1ad2SMatthew G. Knepley ierr = DMPlexPointLocalRead(dm, cells[0], x, &xI);CHKERRQ(ierr); 245254c1ad2SMatthew G. Knepley ierr = DMPlexPointLocalRef(dm, cells[1], x, &xG);CHKERRQ(ierr); 246254c1ad2SMatthew G. Knepley ierr = (*func)(time, fg->centroid, fg->normal, xI, xG, ctx);CHKERRQ(ierr); 247254c1ad2SMatthew G. Knepley } 248254c1ad2SMatthew G. Knepley } 249254c1ad2SMatthew G. Knepley ierr = ISRestoreIndices(faceIS, &faces);CHKERRQ(ierr); 250254c1ad2SMatthew G. Knepley ierr = ISDestroy(&faceIS);CHKERRQ(ierr); 251254c1ad2SMatthew G. Knepley } 252254c1ad2SMatthew G. Knepley } 253254c1ad2SMatthew G. Knepley ierr = VecRestoreArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr); 254254c1ad2SMatthew G. Knepley ierr = VecRestoreArray(locX, &x);CHKERRQ(ierr); 255254c1ad2SMatthew G. Knepley if (Grad) { 256254c1ad2SMatthew G. Knepley ierr = DMRestoreWorkArray(dm, pdim, PETSC_SCALAR, &fx);CHKERRQ(ierr); 257254c1ad2SMatthew G. Knepley ierr = VecRestoreArrayRead(Grad, &grad);CHKERRQ(ierr); 258254c1ad2SMatthew G. Knepley } 259254c1ad2SMatthew G. Knepley PetscFunctionReturn(0); 260254c1ad2SMatthew G. Knepley } 261254c1ad2SMatthew G. Knepley 262254c1ad2SMatthew G. Knepley #undef __FUNCT__ 263924a1b8fSMatthew G. Knepley #define __FUNCT__ "DMPlexTSComputeRHSFunctionFVM" 264924a1b8fSMatthew G. Knepley PetscErrorCode DMPlexTSComputeRHSFunctionFVM(DM dm, PetscReal time, Vec locX, Vec F, void *ctx) 265254c1ad2SMatthew G. Knepley { 266924a1b8fSMatthew G. Knepley PetscDS prob; 267924a1b8fSMatthew G. Knepley void (*riemann)(const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscScalar[], void *); 268924a1b8fSMatthew G. Knepley void *rctx; 269254c1ad2SMatthew G. Knepley PetscFV fvm; 270254c1ad2SMatthew G. Knepley PetscLimiter lim; 271924a1b8fSMatthew G. Knepley Vec faceGeometry, cellGeometry; 272924a1b8fSMatthew G. Knepley Vec Grad = NULL, locGrad; 273924a1b8fSMatthew G. Knepley DM dmFace, dmCell, dmGrad; 2748cd7fcbbSMatthew G. Knepley DMLabel ghostLabel; 275254c1ad2SMatthew G. Knepley PetscCellGeometry fgeom, cgeom; 276254c1ad2SMatthew G. Knepley const PetscScalar *facegeom, *cellgeom, *x, *lgrad; 277254c1ad2SMatthew G. Knepley PetscScalar *grad, *f, *uL, *uR, *fluxL, *fluxR; 278254c1ad2SMatthew G. Knepley PetscReal *centroid, *normal, *vol, *cellPhi; 279254c1ad2SMatthew G. Knepley PetscBool computeGradients; 280254c1ad2SMatthew G. Knepley PetscInt Nf, dim, pdim, fStart, fEnd, numFaces = 0, face, iface, cell, cStart, cEnd, cEndInterior; 281254c1ad2SMatthew G. Knepley PetscErrorCode ierr; 282254c1ad2SMatthew G. Knepley 283254c1ad2SMatthew G. Knepley PetscFunctionBegin; 284924a1b8fSMatthew G. Knepley PetscValidHeaderSpecific(dm,DM_CLASSID,1); 285924a1b8fSMatthew G. Knepley PetscValidHeaderSpecific(locX,VEC_CLASSID,3); 286254c1ad2SMatthew G. Knepley PetscValidHeaderSpecific(F,VEC_CLASSID,5); 287924a1b8fSMatthew G. Knepley ierr = DMPlexTSGetGeometry(dm, &faceGeometry, &cellGeometry, NULL);CHKERRQ(ierr); 288924a1b8fSMatthew G. Knepley ierr = DMPlexTSGetGradientDM(dm, &dmGrad);CHKERRQ(ierr); 289924a1b8fSMatthew G. Knepley ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 290924a1b8fSMatthew G. Knepley ierr = PetscDSGetRiemannSolver(prob, 0, &riemann);CHKERRQ(ierr); 291924a1b8fSMatthew G. Knepley ierr = PetscDSGetContext(prob, 0, &rctx);CHKERRQ(ierr); 292c73cfb54SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 2936dbbd306SMatthew G. Knepley ierr = DMGetNumFields(dm, &Nf);CHKERRQ(ierr); 2946dbbd306SMatthew G. Knepley ierr = DMGetField(dm, 0, (PetscObject *) &fvm);CHKERRQ(ierr); 295254c1ad2SMatthew G. Knepley ierr = PetscFVGetLimiter(fvm, &lim);CHKERRQ(ierr); 2966dbbd306SMatthew G. Knepley ierr = PetscFVGetNumComponents(fvm, &pdim);CHKERRQ(ierr); 297254c1ad2SMatthew G. Knepley ierr = PetscFVGetComputeGradients(fvm, &computeGradients);CHKERRQ(ierr); 298254c1ad2SMatthew G. Knepley if (computeGradients) { 299924a1b8fSMatthew G. Knepley ierr = DMGetGlobalVector(dmGrad, &Grad);CHKERRQ(ierr); 300254c1ad2SMatthew G. Knepley ierr = VecZeroEntries(Grad);CHKERRQ(ierr); 301254c1ad2SMatthew G. Knepley ierr = VecGetArray(Grad, &grad);CHKERRQ(ierr); 302254c1ad2SMatthew G. Knepley } 3036dbbd306SMatthew G. Knepley ierr = DMPlexGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr); 3046dbbd306SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 3056dbbd306SMatthew G. Knepley ierr = VecGetDM(faceGeometry, &dmFace);CHKERRQ(ierr); 3066dbbd306SMatthew G. Knepley ierr = VecGetDM(cellGeometry, &dmCell);CHKERRQ(ierr); 3076dbbd306SMatthew G. Knepley ierr = VecGetArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr); 3086dbbd306SMatthew G. Knepley ierr = VecGetArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr); 3096dbbd306SMatthew G. Knepley ierr = VecGetArrayRead(locX, &x);CHKERRQ(ierr); 310254c1ad2SMatthew G. Knepley /* Count faces and reconstruct gradients */ 3116dbbd306SMatthew G. Knepley for (face = fStart; face < fEnd; ++face) { 312254c1ad2SMatthew G. Knepley const PetscInt *cells; 3131faf85eaSMatthew G. Knepley const PetscFVFaceGeom *fg; 314254c1ad2SMatthew G. Knepley const PetscScalar *cx[2]; 315254c1ad2SMatthew G. Knepley PetscScalar *cgrad[2]; 3168cd7fcbbSMatthew G. Knepley PetscBool boundary; 3178cd7fcbbSMatthew G. Knepley PetscInt ghost, c, pd, d; 3186dbbd306SMatthew G. Knepley 3196dbbd306SMatthew G. Knepley ierr = DMLabelGetValue(ghostLabel, face, &ghost);CHKERRQ(ierr); 3206dbbd306SMatthew G. Knepley if (ghost >= 0) continue; 3216dbbd306SMatthew G. Knepley ++numFaces; 322254c1ad2SMatthew G. Knepley if (!computeGradients) continue; 3238cd7fcbbSMatthew G. Knepley ierr = DMPlexIsBoundaryPoint(dm, face, &boundary);CHKERRQ(ierr); 3248cd7fcbbSMatthew G. Knepley if (boundary) continue; 325254c1ad2SMatthew G. Knepley ierr = DMPlexGetSupport(dm, face, &cells);CHKERRQ(ierr); 326254c1ad2SMatthew G. Knepley ierr = DMPlexPointLocalRead(dmFace, face, facegeom, &fg);CHKERRQ(ierr); 327254c1ad2SMatthew G. Knepley for (c = 0; c < 2; ++c) { 328254c1ad2SMatthew G. Knepley ierr = DMPlexPointLocalRead(dm, cells[c], x, &cx[c]);CHKERRQ(ierr); 329924a1b8fSMatthew G. Knepley ierr = DMPlexPointGlobalRef(dmGrad, cells[c], grad, &cgrad[c]);CHKERRQ(ierr); 330254c1ad2SMatthew G. Knepley } 331254c1ad2SMatthew G. Knepley for (pd = 0; pd < pdim; ++pd) { 332254c1ad2SMatthew G. Knepley PetscScalar delta = cx[1][pd] - cx[0][pd]; 333254c1ad2SMatthew G. Knepley 334254c1ad2SMatthew G. Knepley for (d = 0; d < dim; ++d) { 335254c1ad2SMatthew G. Knepley if (cgrad[0]) cgrad[0][pd*dim+d] += fg->grad[0][d] * delta; 336254c1ad2SMatthew G. Knepley if (cgrad[1]) cgrad[1][pd*dim+d] -= fg->grad[1][d] * delta; 337254c1ad2SMatthew G. Knepley } 338254c1ad2SMatthew G. Knepley } 339254c1ad2SMatthew G. Knepley } 340254c1ad2SMatthew G. Knepley /* Limit interior gradients (using cell-based loop because it generalizes better to vector limiters) */ 341254c1ad2SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 342254c1ad2SMatthew G. Knepley ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 343254c1ad2SMatthew G. Knepley ierr = DMGetWorkArray(dm, pdim, PETSC_REAL, &cellPhi);CHKERRQ(ierr); 344254c1ad2SMatthew G. Knepley for (cell = computeGradients && lim ? cStart : cEnd; cell < cEndInterior; ++cell) { 345254c1ad2SMatthew G. Knepley const PetscInt *faces; 346254c1ad2SMatthew G. Knepley const PetscScalar *cx; 3471faf85eaSMatthew G. Knepley const PetscFVCellGeom *cg; 348254c1ad2SMatthew G. Knepley PetscScalar *cgrad; 349254c1ad2SMatthew G. Knepley PetscInt coneSize, f, pd, d; 350254c1ad2SMatthew G. Knepley 351254c1ad2SMatthew G. Knepley ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr); 352254c1ad2SMatthew G. Knepley ierr = DMPlexGetCone(dm, cell, &faces);CHKERRQ(ierr); 353254c1ad2SMatthew G. Knepley ierr = DMPlexPointLocalRead(dm, cell, x, &cx);CHKERRQ(ierr); 354254c1ad2SMatthew G. Knepley ierr = DMPlexPointLocalRead(dmCell, cell, cellgeom, &cg);CHKERRQ(ierr); 355924a1b8fSMatthew G. Knepley ierr = DMPlexPointGlobalRef(dmGrad, cell, grad, &cgrad);CHKERRQ(ierr); 35671ea8452SMatthew G. Knepley if (!cgrad) continue; /* Unowned overlap cell, we do not compute */ 357254c1ad2SMatthew G. Knepley /* Limiter will be minimum value over all neighbors */ 358254c1ad2SMatthew G. Knepley for (d = 0; d < pdim; ++d) cellPhi[d] = PETSC_MAX_REAL; 359254c1ad2SMatthew G. Knepley for (f = 0; f < coneSize; ++f) { 360254c1ad2SMatthew G. Knepley const PetscScalar *ncx; 3611faf85eaSMatthew G. Knepley const PetscFVCellGeom *ncg; 362254c1ad2SMatthew G. Knepley const PetscInt *fcells; 3638cd7fcbbSMatthew G. Knepley PetscInt face = faces[f], ncell, ghost; 3641475bf87SMatthew G. Knepley PetscReal v[3]; 3658cd7fcbbSMatthew G. Knepley PetscBool boundary; 366254c1ad2SMatthew G. Knepley 367254c1ad2SMatthew G. Knepley ierr = DMLabelGetValue(ghostLabel, face, &ghost);CHKERRQ(ierr); 3688cd7fcbbSMatthew G. Knepley ierr = DMPlexIsBoundaryPoint(dm, face, &boundary);CHKERRQ(ierr); 3698cd7fcbbSMatthew G. Knepley if ((ghost >= 0) || boundary) continue; 370254c1ad2SMatthew G. Knepley ierr = DMPlexGetSupport(dm, face, &fcells);CHKERRQ(ierr); 371254c1ad2SMatthew G. Knepley ncell = cell == fcells[0] ? fcells[1] : fcells[0]; 372254c1ad2SMatthew G. Knepley ierr = DMPlexPointLocalRead(dm, ncell, x, &ncx);CHKERRQ(ierr); 373254c1ad2SMatthew G. Knepley ierr = DMPlexPointLocalRead(dmCell, ncell, cellgeom, &ncg);CHKERRQ(ierr); 3741faf85eaSMatthew G. Knepley DMPlex_WaxpyD_Internal(dim, -1, cg->centroid, ncg->centroid, v); 375254c1ad2SMatthew G. Knepley for (d = 0; d < pdim; ++d) { 376254c1ad2SMatthew G. Knepley /* We use the symmetric slope limited form of Berger, Aftosmis, and Murman 2005 */ 3771faf85eaSMatthew G. Knepley PetscReal phi, flim = 0.5 * PetscRealPart(ncx[d] - cx[d]) / DMPlex_DotD_Internal(dim, &cgrad[d*dim], v); 378254c1ad2SMatthew G. Knepley 379254c1ad2SMatthew G. Knepley ierr = PetscLimiterLimit(lim, flim, &phi);CHKERRQ(ierr); 380254c1ad2SMatthew G. Knepley cellPhi[d] = PetscMin(cellPhi[d], phi); 381254c1ad2SMatthew G. Knepley } 382254c1ad2SMatthew G. Knepley } 383254c1ad2SMatthew G. Knepley /* Apply limiter to gradient */ 384254c1ad2SMatthew G. Knepley for (pd = 0; pd < pdim; ++pd) 385254c1ad2SMatthew G. Knepley /* Scalar limiter applied to each component separately */ 386254c1ad2SMatthew G. Knepley for (d = 0; d < dim; ++d) cgrad[pd*dim+d] *= cellPhi[pd]; 387254c1ad2SMatthew G. Knepley } 388254c1ad2SMatthew G. Knepley ierr = DMRestoreWorkArray(dm, pdim, PETSC_REAL, &cellPhi);CHKERRQ(ierr); 389924a1b8fSMatthew G. Knepley ierr = DMPlexInsertBoundaryValuesFVM_Static(dm, fvm, time, locX, Grad);CHKERRQ(ierr); 390254c1ad2SMatthew G. Knepley if (computeGradients) { 391254c1ad2SMatthew G. Knepley ierr = VecRestoreArray(Grad, &grad);CHKERRQ(ierr); 392924a1b8fSMatthew G. Knepley ierr = DMGetLocalVector(dmGrad, &locGrad);CHKERRQ(ierr); 393924a1b8fSMatthew G. Knepley ierr = DMGlobalToLocalBegin(dmGrad, Grad, INSERT_VALUES, locGrad);CHKERRQ(ierr); 394924a1b8fSMatthew G. Knepley ierr = DMGlobalToLocalEnd(dmGrad, Grad, INSERT_VALUES, locGrad);CHKERRQ(ierr); 395924a1b8fSMatthew G. Knepley ierr = DMRestoreGlobalVector(dmGrad, &Grad);CHKERRQ(ierr); 396254c1ad2SMatthew G. Knepley ierr = VecGetArrayRead(locGrad, &lgrad);CHKERRQ(ierr); 3976dbbd306SMatthew G. Knepley } 3986dbbd306SMatthew G. Knepley ierr = PetscMalloc7(numFaces*dim,¢roid,numFaces*dim,&normal,numFaces*2,&vol,numFaces*pdim,&uL,numFaces*pdim,&uR,numFaces*pdim,&fluxL,numFaces*pdim,&fluxR);CHKERRQ(ierr); 399254c1ad2SMatthew G. Knepley /* Read out values */ 4006dbbd306SMatthew G. Knepley for (face = fStart, iface = 0; face < fEnd; ++face) { 4016dbbd306SMatthew G. Knepley const PetscInt *cells; 4021faf85eaSMatthew G. Knepley const PetscFVFaceGeom *fg; 4031faf85eaSMatthew G. Knepley const PetscFVCellGeom *cgL, *cgR; 404254c1ad2SMatthew G. Knepley const PetscScalar *xL, *xR, *gL, *gR; 4056dbbd306SMatthew G. Knepley PetscInt ghost, d; 4066dbbd306SMatthew G. Knepley 4076dbbd306SMatthew G. Knepley ierr = DMLabelGetValue(ghostLabel, face, &ghost);CHKERRQ(ierr); 4086dbbd306SMatthew G. Knepley if (ghost >= 0) continue; 4096dbbd306SMatthew G. Knepley ierr = DMPlexPointLocalRead(dmFace, face, facegeom, &fg);CHKERRQ(ierr); 4106dbbd306SMatthew G. Knepley ierr = DMPlexGetSupport(dm, face, &cells);CHKERRQ(ierr); 4116dbbd306SMatthew G. Knepley ierr = DMPlexPointLocalRead(dmCell, cells[0], cellgeom, &cgL);CHKERRQ(ierr); 4126dbbd306SMatthew G. Knepley ierr = DMPlexPointLocalRead(dmCell, cells[1], cellgeom, &cgR);CHKERRQ(ierr); 4136dbbd306SMatthew G. Knepley ierr = DMPlexPointLocalRead(dm, cells[0], x, &xL);CHKERRQ(ierr); 4146dbbd306SMatthew G. Knepley ierr = DMPlexPointLocalRead(dm, cells[1], x, &xR);CHKERRQ(ierr); 415254c1ad2SMatthew G. Knepley if (computeGradients) { 4161475bf87SMatthew G. Knepley PetscReal dxL[3], dxR[3]; 417254c1ad2SMatthew G. Knepley 418924a1b8fSMatthew G. Knepley ierr = DMPlexPointLocalRead(dmGrad, cells[0], lgrad, &gL);CHKERRQ(ierr); 419924a1b8fSMatthew G. Knepley ierr = DMPlexPointLocalRead(dmGrad, cells[1], lgrad, &gR);CHKERRQ(ierr); 4201faf85eaSMatthew G. Knepley DMPlex_WaxpyD_Internal(dim, -1, cgL->centroid, fg->centroid, dxL); 4211faf85eaSMatthew G. Knepley DMPlex_WaxpyD_Internal(dim, -1, cgR->centroid, fg->centroid, dxR); 422254c1ad2SMatthew G. Knepley for (d = 0; d < pdim; ++d) { 4231faf85eaSMatthew G. Knepley uL[iface*pdim+d] = xL[d] + DMPlex_DotD_Internal(dim, &gL[d*dim], dxL); 4241faf85eaSMatthew G. Knepley uR[iface*pdim+d] = xR[d] + DMPlex_DotD_Internal(dim, &gR[d*dim], dxR); 425254c1ad2SMatthew G. Knepley } 426254c1ad2SMatthew G. Knepley } else { 4276dbbd306SMatthew G. Knepley for (d = 0; d < pdim; ++d) { 4286dbbd306SMatthew G. Knepley uL[iface*pdim+d] = xL[d]; 4296dbbd306SMatthew G. Knepley uR[iface*pdim+d] = xR[d]; 4306dbbd306SMatthew G. Knepley } 431254c1ad2SMatthew G. Knepley } 4326dbbd306SMatthew G. Knepley for (d = 0; d < dim; ++d) { 4336dbbd306SMatthew G. Knepley centroid[iface*dim+d] = fg->centroid[d]; 4346dbbd306SMatthew G. Knepley normal[iface*dim+d] = fg->normal[d]; 4356dbbd306SMatthew G. Knepley } 4366dbbd306SMatthew G. Knepley vol[iface*2+0] = cgL->volume; 4376dbbd306SMatthew G. Knepley vol[iface*2+1] = cgR->volume; 4386dbbd306SMatthew G. Knepley ++iface; 4396dbbd306SMatthew G. Knepley } 440254c1ad2SMatthew G. Knepley if (computeGradients) { 441254c1ad2SMatthew G. Knepley ierr = VecRestoreArrayRead(locGrad,&lgrad);CHKERRQ(ierr); 442924a1b8fSMatthew G. Knepley ierr = DMRestoreLocalVector(dmGrad, &locGrad);CHKERRQ(ierr); 443254c1ad2SMatthew G. Knepley } 4446dbbd306SMatthew G. Knepley ierr = VecRestoreArrayRead(locX, &x);CHKERRQ(ierr); 4456dbbd306SMatthew G. Knepley ierr = VecRestoreArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr); 4466dbbd306SMatthew G. Knepley ierr = VecRestoreArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr); 4476dbbd306SMatthew G. Knepley fgeom.v0 = centroid; 4486dbbd306SMatthew G. Knepley fgeom.n = normal; 4496dbbd306SMatthew G. Knepley cgeom.vol = vol; 450254c1ad2SMatthew G. Knepley /* Riemann solve */ 451924a1b8fSMatthew G. Knepley ierr = PetscFVIntegrateRHSFunction(fvm, numFaces, Nf, &fvm, 0, fgeom, cgeom, uL, uR, riemann, fluxL, fluxR, rctx);CHKERRQ(ierr); 452254c1ad2SMatthew G. Knepley /* Insert fluxes */ 4536dbbd306SMatthew G. Knepley ierr = VecGetArray(F, &f);CHKERRQ(ierr); 4546dbbd306SMatthew G. Knepley for (face = fStart, iface = 0; face < fEnd; ++face) { 4556dbbd306SMatthew G. Knepley const PetscInt *cells; 4566dbbd306SMatthew G. Knepley PetscScalar *fL, *fR; 4576dbbd306SMatthew G. Knepley PetscInt ghost, d; 4586dbbd306SMatthew G. Knepley 4596dbbd306SMatthew G. Knepley ierr = DMLabelGetValue(ghostLabel, face, &ghost);CHKERRQ(ierr); 4606dbbd306SMatthew G. Knepley if (ghost >= 0) continue; 4616dbbd306SMatthew G. Knepley ierr = DMPlexGetSupport(dm, face, &cells);CHKERRQ(ierr); 4626dbbd306SMatthew G. Knepley ierr = DMPlexPointGlobalRef(dm, cells[0], f, &fL);CHKERRQ(ierr); 4636dbbd306SMatthew G. Knepley ierr = DMPlexPointGlobalRef(dm, cells[1], f, &fR);CHKERRQ(ierr); 4646dbbd306SMatthew G. Knepley for (d = 0; d < pdim; ++d) { 4656dbbd306SMatthew G. Knepley if (fL) fL[d] -= fluxL[iface*pdim+d]; 4666dbbd306SMatthew G. Knepley if (fR) fR[d] += fluxR[iface*pdim+d]; 4676dbbd306SMatthew G. Knepley } 4686dbbd306SMatthew G. Knepley ++iface; 4696dbbd306SMatthew G. Knepley } 4706dbbd306SMatthew G. Knepley ierr = VecRestoreArray(F, &f);CHKERRQ(ierr); 4716dbbd306SMatthew G. Knepley ierr = PetscFree7(centroid,normal,vol,uL,uR,fluxL,fluxR);CHKERRQ(ierr); 472254c1ad2SMatthew G. Knepley PetscFunctionReturn(0); 473254c1ad2SMatthew G. Knepley } 474254c1ad2SMatthew G. Knepley 475254c1ad2SMatthew G. Knepley #undef __FUNCT__ 476254c1ad2SMatthew G. Knepley #define __FUNCT__ "DMPlexTSSetRHSFunctionLocal" 477254c1ad2SMatthew G. Knepley /*@C 478254c1ad2SMatthew G. Knepley DMPlexTSSetRHSFunctionLocal - set a local residual evaluation function 479254c1ad2SMatthew G. Knepley 480254c1ad2SMatthew G. Knepley Logically Collective 481254c1ad2SMatthew G. Knepley 482254c1ad2SMatthew G. Knepley Input Arguments: 483254c1ad2SMatthew G. Knepley + dm - DM to associate callback with 484924a1b8fSMatthew G. Knepley . func - the RHS evaluation function 485924a1b8fSMatthew G. Knepley - ctx - an optional user context 486254c1ad2SMatthew G. Knepley 487254c1ad2SMatthew G. Knepley Level: beginner 488254c1ad2SMatthew G. Knepley 489254c1ad2SMatthew G. Knepley .seealso: DMTSSetRHSFunctionLocal() 490254c1ad2SMatthew G. Knepley @*/ 491924a1b8fSMatthew G. Knepley PetscErrorCode DMPlexTSSetRHSFunctionLocal(DM dm, PetscErrorCode (*func)(DM dm, PetscReal time, Vec X, Vec F, void *ctx), void *ctx) 492254c1ad2SMatthew G. Knepley { 493254c1ad2SMatthew G. Knepley DMTS dmts; 494254c1ad2SMatthew G. Knepley PetscFV fvm; 495254c1ad2SMatthew G. Knepley PetscInt Nf; 496254c1ad2SMatthew G. Knepley PetscBool computeGradients; 497254c1ad2SMatthew G. Knepley PetscErrorCode ierr; 498254c1ad2SMatthew G. Knepley 499254c1ad2SMatthew G. Knepley PetscFunctionBegin; 500254c1ad2SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 501254c1ad2SMatthew G. Knepley ierr = DMGetDMTSWrite(dm, &dmts);CHKERRQ(ierr); 502254c1ad2SMatthew G. Knepley ierr = DMGetNumFields(dm, &Nf);CHKERRQ(ierr); 503254c1ad2SMatthew G. Knepley ierr = DMGetField(dm, 0, (PetscObject *) &fvm);CHKERRQ(ierr); 504254c1ad2SMatthew G. Knepley ierr = PetscFVGetComputeGradients(fvm, &computeGradients);CHKERRQ(ierr); 505924a1b8fSMatthew G. Knepley if (computeGradients) {ierr = DMPlexTSSetupGradient(dm, fvm, dmts);CHKERRQ(ierr);} 506924a1b8fSMatthew G. Knepley ierr = DMTSSetRHSFunctionLocal(dm, func, ctx);CHKERRQ(ierr); 5076dbbd306SMatthew G. Knepley PetscFunctionReturn(0); 5086dbbd306SMatthew G. Knepley } 509*24cdb843SMatthew G. Knepley 510*24cdb843SMatthew G. Knepley #undef __FUNCT__ 511*24cdb843SMatthew G. Knepley #define __FUNCT__ "DMPlexTSComputeIFunctionFEM" 512*24cdb843SMatthew G. Knepley /*@ 513*24cdb843SMatthew G. Knepley DMPlexTSComputeIFunctionFEM - Form the local residual F from the local input X using pointwise functions specified by the user 514*24cdb843SMatthew G. Knepley 515*24cdb843SMatthew G. Knepley Input Parameters: 516*24cdb843SMatthew G. Knepley + dm - The mesh 517*24cdb843SMatthew G. Knepley . t - The time 518*24cdb843SMatthew G. Knepley . X - Local solution 519*24cdb843SMatthew G. Knepley . X_t - Local solution time derivative, or NULL 520*24cdb843SMatthew G. Knepley - user - The user context 521*24cdb843SMatthew G. Knepley 522*24cdb843SMatthew G. Knepley Output Parameter: 523*24cdb843SMatthew G. Knepley . F - Local output vector 524*24cdb843SMatthew G. Knepley 525*24cdb843SMatthew G. Knepley Level: developer 526*24cdb843SMatthew G. Knepley 527*24cdb843SMatthew G. Knepley .seealso: DMPlexComputeJacobianActionFEM() 528*24cdb843SMatthew G. Knepley @*/ 529*24cdb843SMatthew G. Knepley PetscErrorCode DMPlexTSComputeIFunctionFEM(DM dm, PetscReal time, Vec X, Vec X_t, Vec F, void *user) 530*24cdb843SMatthew G. Knepley { 531*24cdb843SMatthew G. Knepley PetscErrorCode ierr; 532*24cdb843SMatthew G. Knepley 533*24cdb843SMatthew G. Knepley PetscFunctionBegin; 534*24cdb843SMatthew G. Knepley ierr = DMPlexComputeResidualFEM_Internal(dm, X, X_t, F, user);CHKERRQ(ierr); 535*24cdb843SMatthew G. Knepley PetscFunctionReturn(0); 536*24cdb843SMatthew G. Knepley } 537