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

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

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

Go to the source code of this file.

Functions

int NavierStokes2DFlux (double *f, double *u, int dir, void *s, double t)
 
int NavierStokes2DStiffFlux (double *f, double *u, int dir, void *s, double t)
 
int NavierStokes2DNonStiffFlux (double *f, double *u, int dir, void *s, double t)
 

Detailed Description

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

Author
Debojyoti Ghosh

Definition in file NavierStokes2DFlux.c.

Function Documentation

◆ NavierStokes2DFlux()

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

Compute the hyperbolic flux function for the 2D 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 \\ (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 \\ (e+p)v \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 or y) for which to compute the flux
sSolver object of type HyPar
tCurrent simulation time

Definition at line 21 of file NavierStokes2DFlux.c.

28 {
29  HyPar *solver = (HyPar*) s;
30  NavierStokes2D *param = (NavierStokes2D*) solver->physics;
31  int *dim = solver->dim_local;
32  int ghosts = solver->ghosts;
33  static int index[_MODEL_NDIMS_], bounds[_MODEL_NDIMS_], offset[_MODEL_NDIMS_];
34 
35  /* set bounds for array index to include ghost points */
36  _ArrayAddCopy1D_(dim,(2*ghosts),bounds,_MODEL_NDIMS_);
37  /* set offset such that index is compatible with ghost point arrangement */
38  _ArraySetValue_(offset,_MODEL_NDIMS_,-ghosts);
39 
40  int done = 0; _ArraySetValue_(index,_MODEL_NDIMS_,0);
41  while (!done) {
42  int p; _ArrayIndex1DWO_(_MODEL_NDIMS_,dim,index,offset,ghosts,p);
43  double rho, vx, vy, e, P;
44  _NavierStokes2DGetFlowVar_((u+_MODEL_NVARS_*p),rho,vx,vy,e,P,param->gamma);
45  _NavierStokes2DSetFlux_((f+_MODEL_NVARS_*p),rho,vx,vy,e,P,dir);
46  _ArrayIncrementIndex_(_MODEL_NDIMS_,bounds,index,done);
47  }
48 
49  return(0);
50 }
#define _NavierStokes2DGetFlowVar_(u, rho, vx, vy, e, P, gamma)
#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)
Structure containing variables and parameters specific to the 2D Navier Stokes equations. This structure contains the physical parameters, variables, and function pointers specific to the 2D Navier-Stokes equations.
#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 _NavierStokes2DSetFlux_(f, rho, vx, vy, e, P, dir)
#define _MODEL_NVARS_
Definition: euler1d.h:58

◆ NavierStokes2DStiffFlux()

int NavierStokes2DStiffFlux ( 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 _NavierStokes2DSetStiffFlux_). 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 (NavierStokes2D::fast_jac) evaluated for the solution at the beginning of each time step ( \({\bf u}_{ref}\) is NavierStokes2D::solution). This is done in NavierStokes2DPreStep().

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 or y) for which to compute the flux
sSolver object of type HyPar
tCurrent simulation time

Definition at line 63 of file NavierStokes2DFlux.c.

70 {
71  HyPar *solver = (HyPar*) s;
72  NavierStokes2D *param = (NavierStokes2D*) solver->physics;
73  int *dim = solver->dim_local;
74  int ghosts = solver->ghosts;
75  static const int JacSize = _MODEL_NVARS_*_MODEL_NVARS_;
76  static int index[_MODEL_NDIMS_], bounds[_MODEL_NDIMS_], offset[_MODEL_NDIMS_];
77 
78  /* set bounds for array index to include ghost points */
79  _ArrayAddCopy1D_(dim,(2*ghosts),bounds,_MODEL_NDIMS_);
80  /* set offset such that index is compatible with ghost point arrangement */
81  _ArraySetValue_(offset,_MODEL_NDIMS_,-ghosts);
82 
83  int done = 0; _ArraySetValue_(index,_MODEL_NDIMS_,0);
84  while (!done) {
85  int p; _ArrayIndex1DWO_(_MODEL_NDIMS_,dim,index,offset,ghosts,p);
86  double *Af = param->fast_jac+(2*p+dir)*JacSize;
88  _ArrayIncrementIndex_(_MODEL_NDIMS_,bounds,index,done);
89  }
90 
91  return(0);
92 }
#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 MatVecMult4(N, y, A, x)
#define _ArrayIndex1DWO_(N, imax, i, offset, ghost, index)
Structure containing variables and parameters specific to the 2D Navier Stokes equations. This structure contains the physical parameters, variables, and function pointers specific to the 2D Navier-Stokes equations.
#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

◆ NavierStokes2DNonStiffFlux()

int NavierStokes2DNonStiffFlux ( 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 _NavierStokes2DSetNonStiffFlux_). 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 NavierStokes2DFlux(), and \(A_f{\bf u}\) is the linearized stiff flux computed in NavierStokes2DStiffFlux().

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 or y) for which to compute the flux
sSolver object of type HyPar
tCurrent simulation time

Definition at line 104 of file NavierStokes2DFlux.c.

111 {
112  HyPar *solver = (HyPar*) s;
113  NavierStokes2D *param = (NavierStokes2D*) solver->physics;
114  int *dim = solver->dim_local;
115  int ghosts = solver->ghosts;
116  static const int JacSize = _MODEL_NVARS_*_MODEL_NVARS_;
117  static int index[_MODEL_NDIMS_], bounds[_MODEL_NDIMS_], offset[_MODEL_NDIMS_];
118  static double ftot[_MODEL_NVARS_], fstiff[_MODEL_NVARS_];
119 
120  /* set bounds for array index to include ghost points */
121  _ArrayAddCopy1D_(dim,(2*ghosts),bounds,_MODEL_NDIMS_);
122  /* set offset such that index is compatible with ghost point arrangement */
123  _ArraySetValue_(offset,_MODEL_NDIMS_,-ghosts);
124 
125  int done = 0; _ArraySetValue_(index,_MODEL_NDIMS_,0);
126  while (!done) {
127  int p; _ArrayIndex1DWO_(_MODEL_NDIMS_,dim,index,offset,ghosts,p);
128  /* compute total flux */
129  double rho, vx, vy, e, P;
130  _NavierStokes2DGetFlowVar_((u+_MODEL_NVARS_*p),rho,vx,vy,e,P,param->gamma);
131  _NavierStokes2DSetFlux_(ftot,rho,vx,vy,e,P,dir);
132  /* compute stiff stuff */
133  double *Af = param->fast_jac+(2*p+dir)*JacSize;
134  MatVecMult4(_MODEL_NVARS_,fstiff,Af,(u+_MODEL_NVARS_*p));
135  /* subtract stiff flux from total flux */
136  _ArraySubtract1D_((f+_MODEL_NVARS_*p),ftot,fstiff,_MODEL_NVARS_);
137  /* Done */
138  _ArrayIncrementIndex_(_MODEL_NDIMS_,bounds,index,done);
139  }
140 
141  return(0);
142 }
#define _NavierStokes2DGetFlowVar_(u, rho, vx, vy, e, P, gamma)
#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 MatVecMult4(N, y, A, x)
#define _ArrayIndex1DWO_(N, imax, i, offset, ghost, index)
Structure containing variables and parameters specific to the 2D Navier Stokes equations. This structure contains the physical parameters, variables, and function pointers specific to the 2D Navier-Stokes equations.
#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 _NavierStokes2DSetFlux_(f, rho, vx, vy, e, P, dir)
#define _MODEL_NVARS_
Definition: euler1d.h:58