HyPar  1.0
Finite-Difference Hyperbolic-Parabolic PDE Solver on Cartesian Grids
NavierStokes3DFunctions.c File Reference

Miscellaneous functions for the 3D Navier Stokes equations. More...

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

Go to the source code of this file.

Functions

int NavierStokes3DRoeAverage (double *uavg, double *uL, double *uR, void *p)
 
int NavierStokes3DComputePressure (double *P, const double *const u, void *s)
 
int NavierStokes3DComputeTemperature (double *T, const double *const u, void *s)
 

Detailed Description

Miscellaneous functions for the 3D Navier Stokes equations.

Author
Debojyoti Ghosh

Definition in file NavierStokes3DFunctions.c.

Function Documentation

◆ NavierStokes3DRoeAverage()

int NavierStokes3DRoeAverage ( double *  uavg,
double *  uL,
double *  uR,
void *  p 
)

Compute the Roe-averaged state for the 3D Navier Stokes equations. This function just calls the macro _NavierStokes3DRoeAverage_ and is not used by any functions within the 3D Navier Stokes module. However, it's necessary to define it and provide it to the the solver object (HyPar) so that it can then send it to interpolation functions for a characteristic-based reconstruction.

Parameters
uavgThe computed Roe-averaged state
uLLeft state (conserved variables)
uRRight state (conserved variables)
pObject of type NavierStokes3D with physics-related variables

Definition at line 18 of file NavierStokes3DFunctions.c.

24 {
25  NavierStokes3D *param = (NavierStokes3D*) p;
27  return(0);
28 }
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.
#define _NavierStokes3DRoeAverage_(uavg, stride, uL, uR, gamma)
static const int _NavierStokes3D_stride_

◆ NavierStokes3DComputePressure()

int NavierStokes3DComputePressure ( double *  P,
const double *const  u,
void *  s 
)

Compute the pressure from the conserved solution on a grid

Parameters
PArray to hold the computed pressure (same layout as u)
uArray with the solution vector
sSolver object of type HyPar

Definition at line 31 of file NavierStokes3DFunctions.c.

35 {
36  HyPar *solver = (HyPar*) s;
37  NavierStokes3D *param = (NavierStokes3D*) solver->physics;
38  int i;
39 
40  int *dim = solver->dim_local;
41  int ghosts = solver->ghosts;
42  int ndims = solver->ndims;
43  int index[ndims], bounds[ndims], offset[ndims];
44 
45  /* set bounds for array index to include ghost points */
46  _ArrayCopy1D_(dim,bounds,ndims);
47  for (i=0; i<ndims; i++) bounds[i] += 2*ghosts;
48 
49  /* set offset such that index is compatible with ghost point arrangement */
50  _ArraySetValue_(offset,ndims,-ghosts);
51 
52  int done = 0; _ArraySetValue_(index,ndims,0);
53  while (!done) {
54  int idx; _ArrayIndex1DWO_(ndims,dim,index,offset,ghosts,idx);
55  double rho, vx, vy, vz, e, pressure;
56  _NavierStokes3DGetFlowVar_((u+_MODEL_NVARS_*idx),_NavierStokes3D_stride_,rho,vx,vy,vz,e,pressure,param->gamma);
57  P[idx] = pressure;
58  _ArrayIncrementIndex_(ndims,bounds,index,done);
59  }
60 
61  return(0);
62 }
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
#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)
#define _NavierStokes3DGetFlowVar_(u, stride, rho, vx, vy, vz, e, P, gamma)

◆ NavierStokes3DComputeTemperature()

int NavierStokes3DComputeTemperature ( double *  T,
const double *const  u,
void *  s 
)

Compute the temperature from the conserved solution on a grid

Parameters
TArray to hold the computed pressure (same layout as u)
uArray with the solution vector
sSolver object of type HyPar

Definition at line 65 of file NavierStokes3DFunctions.c.

69 {
70  HyPar *solver = (HyPar*) s;
71  NavierStokes3D *param = (NavierStokes3D*) solver->physics;
72  int i;
73 
74  int *dim = solver->dim_local;
75  int ghosts = solver->ghosts;
76  int ndims = solver->ndims;
77  int index[ndims], bounds[ndims], offset[ndims];
78 
79  /* set bounds for array index to include ghost points */
80  _ArrayCopy1D_(dim,bounds,ndims);
81  for (i=0; i<ndims; i++) bounds[i] += 2*ghosts;
82 
83  /* set offset such that index is compatible with ghost point arrangement */
84  _ArraySetValue_(offset,ndims,-ghosts);
85 
86  int done = 0; _ArraySetValue_(index,ndims,0);
87  while (!done) {
88  int idx; _ArrayIndex1DWO_(ndims,dim,index,offset,ghosts,idx);
89  double rho, vx, vy, vz, e, pressure;
90  _NavierStokes3DGetFlowVar_((u+_MODEL_NVARS_*idx),_NavierStokes3D_stride_,rho,vx,vy,vz,e,pressure,param->gamma);
91  T[idx] = pressure/rho;
92  _ArrayIncrementIndex_(ndims,bounds,index,done);
93  }
94 
95  return(0);
96 }
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
#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)
#define _NavierStokes3DGetFlowVar_(u, stride, rho, vx, vy, vz, e, P, gamma)