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

Functions to compute the hyperbolic flux for 3D Navier-Stokes equations. More...

#include <stdlib.h>
#include <basic.h>
#include <arrayfunctions.h>
#include <mathfunctions.h>
#include <matmult_native.h>
#include <physicalmodels/navierstokes3d.h>
#include <hypar.h>

Go to the source code of this file.

Functions

int NavierStokes3DFlux (double *f, double *u, int dir, void *s, double t)
 
int NavierStokes3DStiffFlux (double *f, double *u, int dir, void *s, double t)
 
int NavierStokes3DNonStiffFlux (double *f, double *u, int dir, void *s, double t)
 

Detailed Description

Functions to compute the hyperbolic flux for 3D Navier-Stokes equations.

Author
Debojyoti Ghosh

Definition in file NavierStokes3DFlux.c.

Function Documentation

◆ NavierStokes3DFlux()

int NavierStokes3DFlux ( double *  f,
double *  u,
int  dir,
void *  s,
double  t 
)

Compute the hyperbolic flux function for the 3D Navier-Stokes equations:

\begin{eqnarray} dir = x, & {\bf f}\left({\bf u}\right) = \left[\begin{array}{c} \rho u \\ \rho u^2 + p \\ \rho u v \\ \rho u w \\ (e+p)u \end{array}\right], \\ dir = y, & {\bf f}\left({\bf u}\right) = \left[\begin{array}{c} \rho v \\ \rho u v \\ \rho v^2 + p \\ \rho v w \\ (e+p)v \end{array}\right], \\ dir = z, & {\bf f}\left({\bf u}\right) = \left[\begin{array}{c} \rho w \\ \rho u w \\ \rho v w \\ \rho w^2 + p \\ (e+p)w \end{array}\right] \end{eqnarray}

Note: the flux function needs to be computed at the ghost points as well.

Parameters
fArray to hold the computed flux vector (same layout as u)
uArray with the solution vector
dirSpatial dimension (x, y, or z) for which to compute the flux
sSolver object of type HyPar
tCurrent simulation time

Definition at line 23 of file NavierStokes3DFlux.c.

30 {
31  HyPar *solver = (HyPar*) s;
32  NavierStokes3D *param = (NavierStokes3D*) solver->physics;
33  int i;
34 
35  int *dim = solver->dim_local;
36  int ghosts = solver->ghosts;
37  int ndims = solver->ndims;
38  int index[ndims], bounds[ndims], offset[ndims];
39 
40  /* set bounds for array index to include ghost points */
41  _ArrayCopy1D_(dim,bounds,ndims);
42  for (i=0; i<ndims; i++) bounds[i] += 2*ghosts;
43 
44  /* set offset such that index is compatible with ghost point arrangement */
45  _ArraySetValue_(offset,ndims,-ghosts);
46 
47  int done = 0; _ArraySetValue_(index,ndims,0);
48  while (!done) {
49  int p; _ArrayIndex1DWO_(ndims,dim,index,offset,ghosts,p);
50  double rho, vx, vy, vz, e, P;
53  _ArrayIncrementIndex_(ndims,bounds,index,done);
54  }
55 
56  return(0);
57 }
#define _NavierStokes3DSetFlux_(f, stride, rho, vx, vy, vz, e, P, dir)
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)

◆ NavierStokes3DStiffFlux()

int NavierStokes3DStiffFlux ( double *  f,
double *  u,
int  dir,
void *  s,
double  t 
)

Compute the stiff flux, given the solution vector. The stiff flux is the component of the total flux that represents the acoustic modes (see _NavierStokes3DSetStiffFlux_). Here, the linearized approximation to the stiff flux is computed as:

\begin{equation} {\bf f}\left({\bf u}\right) = A_f{\bf u} \end{equation}

where \(A_f = A_f\left({\bf u}_{ref}\right)\) is the fast Jacobian (NavierStokes3D::fast_jac) evaluated for the solution at the beginning of each time step ( \({\bf u}_{ref}\) is NavierStokes3D::solution). This is done in NavierStokes3DPreStep().

Note: the flux function needs to be computed at the ghost points as well.

Parameters
fArray to hold the computed flux vector (same layout as u)
uArray with the solution vector
dirSpatial dimension (x,y, or z) for which to compute the flux
sSolver object of type HyPar
tCurrent simulation time

Definition at line 70 of file NavierStokes3DFlux.c.

77 {
78  HyPar *solver = (HyPar*) s;
79  NavierStokes3D *param = (NavierStokes3D*) solver->physics;
80  int *dim = solver->dim_local;
81  int ghosts = solver->ghosts;
82  static const int JacSize = _MODEL_NVARS_*_MODEL_NVARS_;
83  static int index[_MODEL_NDIMS_], bounds[_MODEL_NDIMS_], offset[_MODEL_NDIMS_];
84 
85  /* set bounds for array index to include ghost points */
86  _ArrayAddCopy1D_(dim,(2*ghosts),bounds,_MODEL_NDIMS_);
87  /* set offset such that index is compatible with ghost point arrangement */
88  _ArraySetValue_(offset,_MODEL_NDIMS_,-ghosts);
89 
90  int done = 0; _ArraySetValue_(index,_MODEL_NDIMS_,0);
91  while (!done) {
92  int p; _ArrayIndex1DWO_(_MODEL_NDIMS_,dim,index,offset,ghosts,p);
93  double *Af = param->fast_jac+(_MODEL_NDIMS_*p+dir)*JacSize;
95  _ArrayIncrementIndex_(_MODEL_NDIMS_,bounds,index,done);
96  }
97 
98  return(0);
99 }
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 _ArrayAddCopy1D_(x, a, y, size)
Structure containing all solver-specific variables and functions.
Definition: hypar.h:23
#define _MODEL_NDIMS_
Definition: euler1d.h:56
#define _ArrayIndex1DWO_(N, imax, i, offset, ghost, index)
#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 MatVecMult5(N, y, A, x)

◆ NavierStokes3DNonStiffFlux()

int NavierStokes3DNonStiffFlux ( double *  f,
double *  u,
int  dir,
void *  s,
double  t 
)

Compute the non-stiff flux, given the solution vector. The non-stiff flux is the component of the total flux that represents the entropy modes (see _NavierStokes3DSetNonStiffFlux_). Here, the linearized approximation to the non-stiff flux is computed as:

\begin{equation} {\bf f}\left({\bf u}\right) - A_f{\bf u} \end{equation}

where \({\bf f}\left({\bf u}\right)\) is the total flux computed in NavierStokes3DFlux(), and \(A_f{\bf u}\) is the linearized stiff flux computed in NavierStokes3DStiffFlux().

Note: the flux function needs to be computed at the ghost points as well.

Parameters
fArray to hold the computed flux vector (same layout as u)
uArray with the solution vector
dirSpatial dimension (x,y, or z) for which to compute the flux
sSolver object of type HyPar
tCurrent simulation time

Definition at line 111 of file NavierStokes3DFlux.c.

118 {
119  HyPar *solver = (HyPar*) s;
120  NavierStokes3D *param = (NavierStokes3D*) solver->physics;
121  int *dim = solver->dim_local;
122  int ghosts = solver->ghosts;
123  static const int JacSize = _MODEL_NVARS_*_MODEL_NVARS_;
124  static int index[_MODEL_NDIMS_], bounds[_MODEL_NDIMS_], offset[_MODEL_NDIMS_];
125  static double ftot[_MODEL_NVARS_], fstiff[_MODEL_NVARS_];
126 
127  /* set bounds for array index to include ghost points */
128  _ArrayAddCopy1D_(dim,(2*ghosts),bounds,_MODEL_NDIMS_);
129  /* set offset such that index is compatible with ghost point arrangement */
130  _ArraySetValue_(offset,_MODEL_NDIMS_,-ghosts);
131 
132  int done = 0; _ArraySetValue_(index,_MODEL_NDIMS_,0);
133  while (!done) {
134  int p; _ArrayIndex1DWO_(_MODEL_NDIMS_,dim,index,offset,ghosts,p);
135  /* compute total flux */
136  double rho, vx, vy, vz, e, P;
138  _NavierStokes3DSetFlux_(ftot,_NavierStokes3D_stride_,rho,vx,vy,vz,e,P,dir);
139  /* compute stiff stuff */
140  double *Af = param->fast_jac+(_MODEL_NDIMS_*p+dir)*JacSize;
141  MatVecMult5(_MODEL_NVARS_,fstiff,Af,(u+_MODEL_NVARS_*p));
142  /* subtract stiff flux from total flux */
143  _ArraySubtract1D_((f+_MODEL_NVARS_*p),ftot,fstiff,_MODEL_NVARS_);
144  /* Done */
145  _ArrayIncrementIndex_(_MODEL_NDIMS_,bounds,index,done);
146  }
147 
148  return(0);
149 }
#define _NavierStokes3DSetFlux_(f, stride, rho, vx, vy, vz, e, P, dir)
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 _ArraySubtract1D_(x, a, b, size)
#define _ArrayAddCopy1D_(x, a, y, size)
Structure containing all solver-specific variables and functions.
Definition: hypar.h:23
#define _MODEL_NDIMS_
Definition: euler1d.h:56
#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 MatVecMult5(N, y, A, x)
#define _NavierStokes3DGetFlowVar_(u, stride, rho, vx, vy, vz, e, P, gamma)