Skip to content

Commit

Permalink
SWE: Update DMPlexCreateSphereMesh() call after change in PETSc's API…
Browse files Browse the repository at this point in the history
… and add compatibility macro
  • Loading branch information
valeriabarra committed Aug 4, 2020
1 parent 854251d commit a653fb1
Show file tree
Hide file tree
Showing 3 changed files with 16 additions and 11 deletions.
20 changes: 10 additions & 10 deletions examples/fluids/shallow-water/qfunctions/geostrophic.h
Original file line number Diff line number Diff line change
Expand Up @@ -83,9 +83,9 @@ CEED_QFUNCTION(ICsSW)(void *ctx, CeedInt Q,
const CeedScalar x[] = {X[0][i], X[1][i], X[2][i]};
CeedScalar q[5];

Exact_SW(2, 0., x, 5, q, ctx);
Exact_SW(2, 0., x, 3, q, ctx);

for (CeedInt j=0; j<5; j++)
for (CeedInt j=0; j<3; j++)
q0[j][i] = q[j];
} // End of Quadrature Point Loop

Expand Down Expand Up @@ -153,15 +153,15 @@ CEED_QFUNCTION(SWExplicit)(void *ctx, CeedInt Q, const CeedScalar *const *in,
// The Physics
// Explicit spatial terms of G_1(t,q):
// Explicit terms multiplying v
// - (omega + f) * khat curl u - grad(|u|^2/2)
// - (omega + f) * khat curl u - grad(|u|^2/2) // TODO: needs fix with weak form
v[0][i] = - wdetJ*(u[0]*du[0][0] + u[1]*du[0][1] + f*u[1]);
// No explicit terms multiplying dv
dv[0][0][i] = 0;
dv[1][0][i] = 0;

// Explicit spatial terms of G_2(t,q):
// Explicit terms multiplying v
// - (omega + f) * khat curl u - grad(|u|^2/2)
// - (omega + f) * khat curl u - grad(|u|^2/2) // TODO: needs fix with weak form
v[1][i] = - wdetJ*(u[0]*du[1][0] + u[1]*du[1][1] - f*u[0]);
// No explicit terms multiplying dv
dv[0][1][i] = 0;
Expand Down Expand Up @@ -210,9 +210,9 @@ CEED_QFUNCTION(SWImplicit)(void *ctx, CeedInt Q, const CeedScalar *const *in,
(*dv)[3][CEED_Q_VLA] = (CeedScalar(*)[3][CEED_Q_VLA])out[1];
// *INDENT-ON*
// Context
const PhysicsContext context = (PhysicsContext)ctx;
const CeedScalar g = context->g;
const CeedScalar H0 = context->H0;
const PhysicsContext context = (PhysicsContext)ctx;
const CeedScalar g = context->g;
const CeedScalar H0 = context->H0;

CeedPragmaSIMD
// Quadrature Point Loop
Expand Down Expand Up @@ -299,9 +299,9 @@ CEED_QFUNCTION(SWJacobian)(void *ctx, CeedInt Q, const CeedScalar *const *in,
CeedScalar (*deltadvdX)[3][CEED_Q_VLA] = (CeedScalar(*)[3][CEED_Q_VLA])out[0];
// *INDENT-ON*
// Context
const PhysicsContext context = (PhysicsContext)ctx;
const CeedScalar g = context->g;
const CeedScalar H0 = context->H0;
const PhysicsContext context = (PhysicsContext)ctx;
const CeedScalar g = context->g;
const CeedScalar H0 = context->H0;

CeedPragmaSIMD
// Quadrature Point Loop
Expand Down
3 changes: 2 additions & 1 deletion examples/fluids/shallow-water/shallowwater.c
Original file line number Diff line number Diff line change
Expand Up @@ -223,7 +223,8 @@ int main(int argc, char **argv) {
} else {
// Create the mesh as a 0-refined sphere. This will create a cubic surface, not a box.
PetscBool simplex = PETSC_FALSE;
ierr = DMPlexCreateSphereMesh(PETSC_COMM_WORLD, topodim, simplex, &dm);
ierr = DMPlexCreateSphereMesh(PETSC_COMM_WORLD, topodim, simplex,
phys_ctx.R, &dm);
CHKERRQ(ierr);
// Set the object name
ierr = PetscObjectSetName((PetscObject)dm, "Sphere"); CHKERRQ(ierr);
Expand Down
4 changes: 4 additions & 0 deletions examples/fluids/shallow-water/src/setup.c
Original file line number Diff line number Diff line change
Expand Up @@ -34,6 +34,10 @@
# define DMPlexRestoreClosureIndices(a,b,c,d,e,f,g,h,i) DMPlexRestoreClosureIndices(a,b,c,d,f,g,i)
#endif

#if PETSC_VERSION_LT(3,14,0)
# define DMPlexCreateSphereMesh(a,b,c,d,e) DMPlexCreateSphereMesh(a,b,c,e)
#endif

problemData problemOptions[] = {
[SWE_ADVECTION] = {
.topodim = 2,
Expand Down

0 comments on commit a653fb1

Please sign in to comment.