HyPar  1.0
Finite-Difference Hyperbolic-Parabolic PDE Solver on Cartesian Grids
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
NavierStokes3DModifiedSolution.c File Reference

Compute the modified solution for the well-balanced treatment of gravitational source terms. More...

#include <basic.h>
#include <arrayfunctions.h>
#include <physicalmodels/navierstokes3d.h>
#include <mpivars.h>
#include <hypar.h>

Go to the source code of this file.

Functions

int NavierStokes3DModifiedSolution (double *uC, double *u, int d, void *s, void *m, double waqt)
 

Detailed Description

Compute the modified solution for the well-balanced treatment of gravitational source terms.

Author
Debojyoti Ghosh

Definition in file NavierStokes3DModifiedSolution.c.

Function Documentation

int NavierStokes3DModifiedSolution ( double *  uC,
double *  u,
int  d,
void *  s,
void *  m,
double  waqt 
)

This function computes the modified solution for the well-balanced treatment of the gravitational source terms. The modified solution vector is given by

\begin{equation} {\bf u}^* = \left[\begin{array}{c} \rho \varrho^{-1}\left(x,y\right) \\ \rho u \varrho^{-1}\left(x,y\right) \\ \rho v \varrho^{-1}\left(x,y\right) \\ \rho w \varrho^{-1}\left(x,y\right) \\ e^* \end{array}\right] \end{equation}

where

\begin{equation} e^* = \frac{p \varphi^{-1}\left(x,y\right)}{\gamma-1} + \frac{1}{2}\rho \varrho^{-1}\left(x,y\right) \left(u^2+v^2+w^2\right) \end{equation}

\(\varrho\) and \(\varphi\) are computed in NavierStokes3DGravityField(). For flows without gravity, \({\bf u}^* = {\bf u}\).

References:

  • Ghosh, D., Constantinescu, E.M., Well-Balanced Formulation of Gravitational Source Terms for Conservative Finite-Difference Atmospheric Flow Solvers, AIAA Paper 2015-2889, 7th AIAA Atmospheric and Space Environments Conference, June 22-26, 2015, Dallas, TX, http://dx.doi.org/10.2514/6.2015-2889
  • Ghosh, D., Constantinescu, E.M., A Well-Balanced, Conservative Finite-Difference Algorithm for Atmospheric Flows, AIAA Journal, 54 (4), 2016, pp. 1370-1385, http://dx.doi.org/10.2514/1.J054580.
Parameters
uCArray to hold the computed modified solution
uSolution vector array
dspatial dimension (not used)
sSolver object of type HyPar
mMPI object of time MPIVariables
waqtCurrent simulation time

Definition at line 31 of file NavierStokes3DModifiedSolution.c.

39 {
40  HyPar *solver = (HyPar*) s;
41  NavierStokes3D *param = (NavierStokes3D*) solver->physics;
42 
43  int ghosts = solver->ghosts;
44  int *dim = solver->dim_local;
45  int ndims = solver->ndims;
46  int index[ndims], bounds[ndims], offset[ndims];
47 
48  /* set bounds for array index to include ghost points */
49  _ArrayCopy1D_(dim,bounds,ndims);
50  int i; for (i=0; i<ndims; i++) bounds[i] += 2*ghosts;
51 
52  /* set offset such that index is compatible with ghost point arrangement */
53  _ArraySetValue_(offset,ndims,-ghosts);
54 
55  int done = 0; _ArraySetValue_(index,ndims,0);
56  double inv_gamma_m1 = 1.0 / (param->gamma-1.0);
57  while (!done) {
58  int p; _ArrayIndex1DWO_(ndims,dim,index,offset,ghosts,p);
59  double rho, uvel, vvel, wvel, E, P;
60  _NavierStokes3DGetFlowVar_((u+_MODEL_NVARS_*p),_NavierStokes3D_stride_,rho,uvel,vvel,wvel,E,P,param->gamma);
61  uC[_MODEL_NVARS_*p+0] = u[_MODEL_NVARS_*p+0] * param->grav_field_f[p];
62  uC[_MODEL_NVARS_*p+1] = u[_MODEL_NVARS_*p+1] * param->grav_field_f[p];
63  uC[_MODEL_NVARS_*p+2] = u[_MODEL_NVARS_*p+2] * param->grav_field_f[p];
64  uC[_MODEL_NVARS_*p+3] = u[_MODEL_NVARS_*p+3] * param->grav_field_f[p];
65  uC[_MODEL_NVARS_*p+4] = (P*inv_gamma_m1)*(1.0/param->grav_field_g[p]) + (0.5*rho*(uvel*uvel+vvel*vvel+wvel*wvel))*param->grav_field_f[p];
66  _ArrayIncrementIndex_(ndims,bounds,index,done);
67  }
68 
69  return(0);
70 }
#define _NavierStokes3DGetFlowVar_(u, stride, rho, vx, vy, vz, e, P, gamma)
#define _ArraySetValue_(x, size, value)
#define _ArrayIncrementIndex_(N, imax, i, done)
#define _MODEL_NVARS_
Definition: euler1d.h:58
void * physics
Definition: hypar.h:266
int * dim_local
Definition: hypar.h:37
int ghosts
Definition: hypar.h:52
double * grav_field_f
#define _ArrayIndex1DWO_(N, imax, i, offset, ghost, index)
#define _ArrayCopy1D_(x, y, size)
double * grav_field_g
Structure containing variables and parameters specific to the 3D Navier Stokes equations. This structure contains the physical parameters, variables, and function pointers specific to the 3D Navier-Stokes equations.
int ndims
Definition: hypar.h:26
Structure containing all solver-specific variables and functions.
Definition: hypar.h:23
static const int _NavierStokes3D_stride_