HyPar  1.0
Finite-Difference Hyperbolic-Parabolic PDE Solver on Cartesian Grids
NavierStokes3DModifiedSolution.c
Go to the documentation of this file.
1 
5 #include <basic.h>
6 #include <arrayfunctions.h>
8 #include <mpivars.h>
9 #include <hypar.h>
10 
32  double *uC,
33  double *u,
34  int d,
35  void *s,
36  void *m,
37  double waqt
38  )
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 }
MPI related function definitions.
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.
Some basic definitions and macros.
double * grav_field_g
3D Navier Stokes equations (compressible flows)
int ndims
Definition: hypar.h:26
Structure containing all solver-specific variables and functions.
Definition: hypar.h:23
Contains structure definition for hypar.
#define _ArrayIndex1DWO_(N, imax, i, offset, ghost, index)
static const int _NavierStokes3D_stride_
#define _ArraySetValue_(x, size, value)
int * dim_local
Definition: hypar.h:37
#define _ArrayIncrementIndex_(N, imax, i, done)
void * physics
Definition: hypar.h:266
int ghosts
Definition: hypar.h:52
#define _MODEL_NVARS_
Definition: euler1d.h:58
#define _ArrayCopy1D_(x, y, size)
Contains macros and function definitions for common array operations.
int NavierStokes3DModifiedSolution(double *uC, double *u, int d, void *s, void *m, double waqt)
double * grav_field_f
#define _NavierStokes3DGetFlowVar_(u, stride, rho, vx, vy, vz, e, P, gamma)