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

2D Navier Stokes equations (compressible flows) More...

#include <basic.h>

Go to the source code of this file.

Data Structures

struct  NavierStokes2D
 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. More...
 

Macros

#define _NAVIER_STOKES_2D_   "navierstokes2d"
 
#define _MODEL_NDIMS_   2
 
#define _MODEL_NVARS_   4
 
#define _ROE_   "roe"
 
#define _RF_   "rf-char"
 
#define _LLF_   "llf-char"
 
#define _SWFS_   "steger-warming"
 
#define _RUSANOV_   "rusanov"
 
#define _XDIR_   0
 
#define _YDIR_   1
 
#define _NavierStokes2DGetFlowVar_(u, rho, vx, vy, e, P, gamma)
 
#define _NavierStokes2DSetFlux_(f, rho, vx, vy, e, P, dir)
 
#define _NavierStokes2DSetStiffFlux_(f, rho, vx, vy, e, P, dir, gamma)
 
#define _NavierStokes2DSetNonStiffFlux_(f, rho, vx, vy, e, P, dir, gamma)
 
#define _NavierStokes2DRoeAverage_(uavg, uL, uR, gamma)
 
#define _NavierStokes2DEigenvalues_(u, D, gamma, dir)
 
#define _NavierStokes2DLeftEigenvectors_(u, L, ga, dir)
 
#define _NavierStokes2DRightEigenvectors_(u, R, ga, dir)
 

Functions

int NavierStokes2DInitialize (void *, void *)
 
int NavierStokes2DCleanup (void *)
 

Detailed Description

2D Navier Stokes equations (compressible flows)

Author
Debojyoti Ghosh

2D Navier-Stokes equations for viscous and inviscid compressible flows (with gravitational terms)

\begin{equation} \frac {\partial} {\partial t} \left[\begin{array}{c} \rho \\ \rho u \\ \rho v \\ e \end{array}\right] + \frac {\partial} {\partial x} \left[\begin{array}{c} \rho u \\ \rho u^2 + p \\ \rho u v \\ (e+p) u\end{array}\right] + \frac {\partial} {\partial y} \left[\begin{array}{c} \rho v \\ \rho u v \\ \rho v^2 + p \\ (e+p) v \end{array}\right] = \frac {\partial} {\partial x} \left[\begin{array}{c} 0 \\ \tau_{xx} \\ \tau_{yx} \\ u \tau_{xx} + v \tau_{yx} - q_x \end{array}\right] + \frac {\partial} {\partial y} \left[\begin{array}{c} 0 \\ \tau_{xy} \\ \tau_{yy} \\ u \tau_{xy} + v \tau_{yy} - q_y \end{array}\right] + \left[\begin{array}{c} 0 \\ -\rho {\bf g}\cdot{\bf \hat{i}} \\ -\rho {\bf g}\cdot{\bf \hat{j}} \\ -\rho u {\bf g}\cdot{\bf \hat{i}} - \rho v {\bf g}\cdot{\bf \hat{j}} \end{array}\right] \end{equation}

where \({\bf g}\) is the gravitational force vector per unit mass, \({\bf \hat{i}},{\bf \hat{j}}\) are the unit vectors along the x and y, the viscous terms are given by

\begin{align} \tau_{ij} &= \frac{\mu}{Re_\infty} \left[ \left( \frac{\partial u_i}{\partial x_j} + \frac{\partial u_j}{\partial x_i}\right) - \frac{2}{3}\frac{\partial u_k}{\partial x_k} \delta_{ij} \right], \\ q_i &= - \frac{\mu}{\left(\gamma-1\right)Re_\infty Pr} \frac{\partial T}{\partial x_i} \end{align}

with \(\mu\) being the viscosity coefficient (computed using Sutherland's law), and the equation of state is

\begin{equation} e = \frac {p} {\gamma-1} + \frac{1}{2} \rho \left(u^2 + v^2\right) \end{equation}

References for the governing equations (as well as non-dimensional form):-

  • Tannehill, Anderson and Pletcher, Computational Fluid Mechanics and Heat Transfer, Chapter 5, Section 5.1.7 (However, the non-dimensional velocity and the Reynolds number is based on speed of sound, instead of the freestream velocity).

Reference for the well-balanced treatment of gravitational source term:

  • 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

Reference for the partitioning of the flux into its stiff (acoustic) and non-stiff (convective) components:

  • Ghosh, D., Constantinescu, E. M., Semi-Implicit Time Integration of Atmospheric Flows with Characteristic-Based Flux Partitioning, SIAM Journal on Scientific Computing, 38 (3), 2016, A1848-A1875, http://dx.doi.org/10.1137/15M1044369.

Definition in file navierstokes2d.h.


Data Structure Documentation

struct NavierStokes2D

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.

Definition at line 374 of file navierstokes2d.h.

Data Fields
double gamma

Ratio of heat capacities

char upw_choice[_MAX_STRING_SIZE_]

choice of upwinding

double grav_x

acceleration due to gravity in x

double grav_y

acceleration due to gravity in y

double rho0

reference density at zero altitude for flows with gravity

double p0

reference pressure at zero altitude for flows with gravity

double Re

Reynolds number

double Pr

Prandtl number

double Minf

Freestream Mach number

double C1

Sutherlands law constants

double C2

Sutherlands law constants

double R

universal Gas constant

double * grav_field_f

density variation function ( \(\varrho\)) for hydrostatic equilibrium for flows with gravity

double * grav_field_g

pressure variation function ( \(\varrho\)) for hydrostatic equilibrium for flows with gravity

double * fast_jac

"Fast" Jacobian of the flux function (comprising the acoustic modes)

double * solution

array to store the solution at the beginning of each time step

int HB

< Choice of hydrostatic balance for flows with gravity (1 - isothermal equilibrium, 2 - constant potential temperature 3 - stratified atmosphere with a Brunt-Vaisala frequency)

double N_bv

the Brunt-Vaisala frequency for NavierStokes2D::HB = 3

double * gpu_grav_field_f
double * gpu_grav_field_g
double * gpu_fast_jac
double * gpu_solution

Macro Definition Documentation

#define _NAVIER_STOKES_2D_   "navierstokes2d"

2D Navier Stokes equations

Definition at line 47 of file navierstokes2d.h.

#define _MODEL_NDIMS_   2

Number of spatial dimensions

Definition at line 53 of file navierstokes2d.h.

#define _MODEL_NVARS_   4

Number of variables per grid point

Definition at line 55 of file navierstokes2d.h.

#define _ROE_   "roe"

Roe's upwinding scheme

Definition at line 59 of file navierstokes2d.h.

#define _RF_   "rf-char"

Characteristic-based Roe-fixed scheme

Definition at line 61 of file navierstokes2d.h.

#define _LLF_   "llf-char"

Characteristic-based local Lax-Friedrich scheme

Definition at line 63 of file navierstokes2d.h.

#define _SWFS_   "steger-warming"

Steger-Warming flux splitting scheme

Definition at line 65 of file navierstokes2d.h.

#define _RUSANOV_   "rusanov"

Rusanov's upwinding scheme

Definition at line 67 of file navierstokes2d.h.

#define _XDIR_   0

dimension corresponding to the x spatial dimension

Definition at line 71 of file navierstokes2d.h.

#define _YDIR_   1

dimension corresponding to the y spatial dimension

Definition at line 73 of file navierstokes2d.h.

#define _NavierStokes2DGetFlowVar_ (   u,
  rho,
  vx,
  vy,
  e,
  P,
  gamma 
)
Value:
{ \
double vsq; \
rho = u[0]; \
vx = (rho==0) ? 0 : u[1] / rho; \
vy = (rho==0) ? 0 : u[2] / rho; \
e = u[3]; \
vsq = (vx*vx) + (vy*vy); \
P = (e - 0.5*rho*vsq) * (gamma-1.0); \
}

Get the flow variables from the conserved solution vector.

\begin{equation} {\bf u} = \left[\begin{array}{c} \rho \\ \rho u \\ \rho v \\ e \end{array}\right] \end{equation}

Definition at line 81 of file navierstokes2d.h.

#define _NavierStokes2DSetFlux_ (   f,
  rho,
  vx,
  vy,
  e,
  P,
  dir 
)
Value:
{ \
if (dir == _XDIR_) { \
f[0] = rho * vx; \
f[1] = rho * vx * vx + P; \
f[2] = rho * vx * vy; \
f[3] = (e + P) * vx; \
} else if (dir == _YDIR_) { \
f[0] = rho * vy; \
f[1] = rho * vy * vx; \
f[2] = rho * vy * vy + P; \
f[3] = (e + P) * vy; \
} \
}
#define _XDIR_
#define _YDIR_

Compute the flux vector, given the flow variables

\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}

Definition at line 99 of file navierstokes2d.h.

#define _NavierStokes2DSetStiffFlux_ (   f,
  rho,
  vx,
  vy,
  e,
  P,
  dir,
  gamma 
)
Value:
{ \
double gamma_inv = 1.0/gamma; \
if (dir == _XDIR_) { \
f[0] = gamma_inv * rho * vx; \
f[1] = gamma_inv * rho * vx * vx + P; \
f[2] = gamma_inv * rho * vx * vy; \
f[3] = (e + P) * vx - 0.5 * gamma_inv * (gamma-1.0) * rho * (vx*vx+vy*vy) * vx; \
} else if (dir == _YDIR_) { \
f[0] = gamma_inv * rho * vy; \
f[1] = gamma_inv * rho * vy * vx; \
f[2] = gamma_inv * rho * vy * vy + P; \
f[3] = (e + P) * vy - 0.5 * gamma_inv * (gamma-1.0) * rho * (vx*vx+vy*vy) * vy; \
} \
}
#define _XDIR_
#define _YDIR_

Compute the stiff flux vector (comprising the acoustic modes only), given the flow variables

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

Reference:

  • Ghosh, D., Constantinescu, E. M., Semi-Implicit Time Integration of Atmospheric Flows with Characteristic-Based Flux Partitioning, SIAM Journal on Scientific Computing, 38 (3), 2016, A1848-A1875, http://dx.doi.org/10.1137/15M1044369.

Definition at line 125 of file navierstokes2d.h.

#define _NavierStokes2DSetNonStiffFlux_ (   f,
  rho,
  vx,
  vy,
  e,
  P,
  dir,
  gamma 
)
Value:
{ \
double gamma_inv = 1.0/gamma; \
if (dir == _XDIR_) { \
f[0] = (gamma-1.0) * gamma_inv * rho * vx; \
f[1] = (gamma-1.0) * gamma_inv * rho * vx * vx; \
f[2] = (gamma-1.0) * gamma_inv * rho * vx * vy; \
f[3] = 0.5 * gamma_inv * (gamma-1.0) * rho * (vx*vx+vy*vy) * vx; \
} else if (dir == _YDIR_) { \
f[0] = (gamma-1.0) * gamma_inv * rho * vy; \
f[1] = (gamma-1.0) * gamma_inv * rho * vy * vx; \
f[2] = (gamma-1.0) * gamma_inv * rho * vy * vy; \
f[3] = 0.5 * gamma_inv * (gamma-1.0) * rho * (vx*vx+vy*vy) * vy; \
} \
}
#define _XDIR_
#define _YDIR_

Compute the non-stiff flux vector (comprising the entropy modes only), given the flow variables

\begin{eqnarray} dir = x, & {\bf f}\left({\bf u}\right) = \left[\begin{array}{c} \frac{\gamma-1}{\gamma}\rho u \\ \frac{\gamma-1}{\gamma}\rho u^2 \\ \frac{\gamma-11}{\gamma}\rho u v \\ \frac{1}{2} \frac{\gamma-1}{\gamma}\rho\left(u^2+v^2\right)u \end{array}\right], \\ dir = y, & {\bf f}\left({\bf u}\right) = \left[\begin{array}{c} \frac{\gamma-1}{\gamma}\rho v \\ \frac{\gamma-1}{\gamma}\rho u v \\ \frac{\gamma-1}{\gamma}\rho v^2 \\ \frac{1}{2} \frac{\gamma-1}{\gamma}\rho\left(u^2+v^2\right)v \end{array}\right] \end{eqnarray}

Reference:

  • Ghosh, D., Constantinescu, E. M., Semi-Implicit Time Integration of Atmospheric Flows with Characteristic-Based Flux Partitioning, SIAM Journal on Scientific Computing, 38 (3), 2016, A1848-A1875, http://dx.doi.org/10.1137/15M1044369.

Definition at line 152 of file navierstokes2d.h.

#define _NavierStokes2DRoeAverage_ (   uavg,
  uL,
  uR,
  gamma 
)

Compute the Roe-average of two solutions.

Definition at line 171 of file navierstokes2d.h.

#define _NavierStokes2DEigenvalues_ (   u,
  D,
  gamma,
  dir 
)
Value:
{ \
double rho,vx,vy,e,P,c,vn,vsq; \
rho = u[0]; \
vx = u[1] / rho; \
vy = u[2] / rho; \
e = u[3]; \
vsq = (vx*vx) + (vy*vy); \
P = (e - 0.5*rho*vsq) * (gamma-1.0); \
c = sqrt(gamma*P/rho); \
if (dir == _XDIR_) vn = vx; \
else if (dir == _YDIR_) vn = vy; \
else vn = 0.0; \
if (dir == _XDIR_) {\
D[0*_MODEL_NVARS_+0] = vn-c; D[0*_MODEL_NVARS_+1] = 0; D[0*_MODEL_NVARS_+2] = 0; D[0*_MODEL_NVARS_+3] = 0; \
D[1*_MODEL_NVARS_+0] = 0; D[1*_MODEL_NVARS_+1] = vn+c; D[1*_MODEL_NVARS_+2] = 0; D[1*_MODEL_NVARS_+3] = 0; \
D[2*_MODEL_NVARS_+0] = 0; D[2*_MODEL_NVARS_+1] = 0; D[2*_MODEL_NVARS_+2] = vn; D[2*_MODEL_NVARS_+3] = 0; \
D[3*_MODEL_NVARS_+0] = 0; D[3*_MODEL_NVARS_+1] = 0; D[3*_MODEL_NVARS_+2] = 0; D[3*_MODEL_NVARS_+3] = vn; \
} else if (dir == _YDIR_) { \
D[0*_MODEL_NVARS_+0] = vn-c; D[0*_MODEL_NVARS_+1] = 0; D[0*_MODEL_NVARS_+2] = 0; D[0*_MODEL_NVARS_+3] = 0; \
D[1*_MODEL_NVARS_+0] = 0; D[1*_MODEL_NVARS_+1] = vn; D[1*_MODEL_NVARS_+2] = 0; D[1*_MODEL_NVARS_+3] = 0; \
D[2*_MODEL_NVARS_+0] = 0; D[2*_MODEL_NVARS_+1] = 0; D[2*_MODEL_NVARS_+2] = vn+c; D[2*_MODEL_NVARS_+3] = 0; \
D[3*_MODEL_NVARS_+0] = 0; D[3*_MODEL_NVARS_+1] = 0; D[3*_MODEL_NVARS_+2] = 0; D[3*_MODEL_NVARS_+3] = vn; \
}\
}
#define _XDIR_
#define _MODEL_NVARS_
#define _YDIR_

Compute the eigenvalues, given a solution vector in terms of the conserved variables. The eigenvalues are returned as a matrix D whose diagonal values are the eigenvalues. Admittedly, this is inefficient. The matrix D is stored in a row-major format.

Definition at line 213 of file navierstokes2d.h.

#define _NavierStokes2DLeftEigenvectors_ (   u,
  L,
  ga,
  dir 
)

Compute the left eigenvectors, given a solution vector in terms of the conserved variables. The eigenvectors are returned as a matrix L whose rows correspond to each eigenvector. The matrix L is stored in the row-major format.

Reference:

Definition at line 247 of file navierstokes2d.h.

#define _NavierStokes2DRightEigenvectors_ (   u,
  R,
  ga,
  dir 
)

Compute the right eigenvectors, given a solution vector in terms of the conserved variables. The eigenvectors are returned as a matrix R whose columns correspond to each eigenvector. The matrix R is stored in the row-major format.

Reference:

Definition at line 309 of file navierstokes2d.h.

Function Documentation

int NavierStokes2DInitialize ( void *  s,
void *  m 
)

Initialize the 2D Navier-Stokes (NavierStokes2D) module: Sets the default parameters, read in and set physics-related parameters, and set the physics-related function pointers in HyPar.

This file reads the file "physics.inp" that must have the following format:

begin
    <keyword>   <value>
    <keyword>   <value>
    <keyword>   <value>
    ...
    <keyword>   <value>
end

where the list of keywords are:

Keyword name Type Variable Default value
gamma double NavierStokes2D::gamma 1.4
Pr double NavierStokes2D::Pr 0.72
Re double NavierStokes2D::Re -1
Minf double NavierStokes2D::Minf 1.0
gravity double,doubleNavierStokes2D::grav_x,NavierStokes2D::grav_y 0.0,0.0
rho_ref double NavierStokes2D::rho0 1.0
p_ref double NavierStokes2D::p0 1.0
HB int NavierStokes2D::HB 1
R double NavierStokes2D::R 1.0
upwinding char[] NavierStokes2D::upw_choice "roe" (_ROE_)

Note: "physics.inp" is optional; if absent, default values will be used.

Parameters
sSolver object of type HyPar
mMPI object of type MPIVariables

Definition at line 103 of file NavierStokes2DInitialize.c.

107 {
108  HyPar *solver = (HyPar*) s;
109  MPIVariables *mpi = (MPIVariables*) m;
110  NavierStokes2D *physics = (NavierStokes2D*) solver->physics;
111  int ferr = 0;
112 
113  static int count = 0;
114 
115  if (solver->nvars != _MODEL_NVARS_) {
116  fprintf(stderr,"Error in NavierStokes2DInitialize(): nvars has to be %d.\n",_MODEL_NVARS_);
117  return(1);
118  }
119  if (solver->ndims != _MODEL_NDIMS_) {
120  fprintf(stderr,"Error in NavierStokes2DInitialize(): ndims has to be %d.\n",_MODEL_NDIMS_);
121  return(1);
122  }
123 
124  /* default values */
125  physics->gamma = 1.4;
126  physics->Pr = 0.72;
127  physics->Re = -1;
128  physics->Minf = 1.0;
129  physics->C1 = 1.458e-6;
130  physics->C2 = 110.4;
131  physics->grav_x = 0.0;
132  physics->grav_y = 0.0;
133  physics->rho0 = 1.0;
134  physics->p0 = 1.0;
135  physics->HB = 1;
136  physics->R = 1.0;
137  physics->N_bv = 0.0;
138  strcpy(physics->upw_choice,"roe");
139 
140  /* reading physical model specific inputs - all processes */
141  if (!mpi->rank) {
142  FILE *in;
143  if (!count) printf("Reading physical model inputs from file \"physics.inp\".\n");
144  in = fopen("physics.inp","r");
145  if (!in) printf("Warning: File \"physics.inp\" not found. Using default values.\n");
146  else {
147  char word[_MAX_STRING_SIZE_];
148  ferr = fscanf(in,"%s",word); if (ferr != 1) return(1);
149  if (!strcmp(word, "begin")){
150  while (strcmp(word, "end")){
151  ferr = fscanf(in,"%s",word); if (ferr != 1) return(1);
152  if (!strcmp(word, "gamma")) {
153  ferr = fscanf(in,"%lf",&physics->gamma); if (ferr != 1) return(1);
154  } else if (!strcmp(word,"upwinding")) {
155  ferr = fscanf(in,"%s",physics->upw_choice); if (ferr != 1) return(1);
156  } else if (!strcmp(word,"Pr")) {
157  ferr = fscanf(in,"%lf",&physics->Pr); if (ferr != 1) return(1);
158  } else if (!strcmp(word,"Re")) {
159  ferr = fscanf(in,"%lf",&physics->Re); if (ferr != 1) return(1);
160  } else if (!strcmp(word,"Minf")) {
161  ferr = fscanf(in,"%lf",&physics->Minf); if (ferr != 1) return(1);
162  } else if (!strcmp(word,"gravity")) {
163  ferr = fscanf(in,"%lf",&physics->grav_x); if (ferr != 1) return(1);
164  ferr = fscanf(in,"%lf",&physics->grav_y); if (ferr != 1) return(1);
165  } else if (!strcmp(word,"rho_ref")) {
166  ferr = fscanf(in,"%lf",&physics->rho0); if (ferr != 1) return(1);
167  } else if (!strcmp(word,"p_ref")) {
168  ferr = fscanf(in,"%lf",&physics->p0); if (ferr != 1) return(1);
169  } else if (!strcmp(word,"HB")) {
170  ferr = fscanf(in,"%d",&physics->HB); if (ferr != 1) return(1);
171  if (physics->HB==3) {
172  ferr = fscanf(in,"%lf",&physics->N_bv); if (ferr != 1) return(1);
173  }
174  } else if (!strcmp(word,"R")) {
175  ferr = fscanf(in,"%lf",&physics->R); if (ferr != 1) return(1);
176  } else if (strcmp(word,"end")) {
177  char useless[_MAX_STRING_SIZE_];
178  ferr = fscanf(in,"%s",useless); if (ferr != 1) return(ferr);
179  printf("Warning: keyword %s in file \"physics.inp\" with value %s not ",word,useless);
180  printf("recognized or extraneous. Ignoring.\n");
181  }
182  }
183  } else {
184  fprintf(stderr,"Error: Illegal format in file \"physics.inp\".\n");
185  return(1);
186  }
187  }
188  fclose(in);
189  }
190 
192  IERR MPIBroadcast_double (&physics->gamma ,1 ,0,&mpi->world); CHECKERR(ierr);
193  IERR MPIBroadcast_double (&physics->Pr ,1 ,0,&mpi->world); CHECKERR(ierr);
194  IERR MPIBroadcast_double (&physics->Re ,1 ,0,&mpi->world); CHECKERR(ierr);
195  IERR MPIBroadcast_double (&physics->Minf ,1 ,0,&mpi->world); CHECKERR(ierr);
196  IERR MPIBroadcast_double (&physics->grav_x ,1 ,0,&mpi->world); CHECKERR(ierr);
197  IERR MPIBroadcast_double (&physics->grav_y ,1 ,0,&mpi->world); CHECKERR(ierr);
198  IERR MPIBroadcast_double (&physics->rho0 ,1 ,0,&mpi->world); CHECKERR(ierr);
199  IERR MPIBroadcast_double (&physics->p0 ,1 ,0,&mpi->world); CHECKERR(ierr);
200  IERR MPIBroadcast_double (&physics->R ,1 ,0,&mpi->world); CHECKERR(ierr);
201  IERR MPIBroadcast_double (&physics->N_bv ,1 ,0,&mpi->world); CHECKERR(ierr);
202  IERR MPIBroadcast_integer (&physics->HB ,1 ,0,&mpi->world); CHECKERR(ierr);
203 
204  /* Scaling the Reynolds number with the M_inf */
205  physics->Re /= physics->Minf;
206 
207  /* check that a well-balanced upwinding scheme is being used for cases with gravity */
208  if ( ((physics->grav_x != 0.0) || (physics->grav_y != 0.0))
209  && (strcmp(physics->upw_choice,_LLF_ ))
210  && (strcmp(physics->upw_choice,_RUSANOV_))
211  && (strcmp(physics->upw_choice,_ROE_ )) ) {
212  if (!mpi->rank) {
213  fprintf(stderr,"Error in NavierStokes2DInitialize: %s, %s or %s upwinding is needed for flows ",_LLF_,_ROE_,_RUSANOV_);
214  fprintf(stderr,"with gravitational forces.\n");
215  }
216  return(1);
217  }
218  /* check that solver has the correct choice of diffusion formulation */
219  if (strcmp(solver->spatial_type_par,_NC_2STAGE_) && (physics->Re > 0)) {
220  if (!mpi->rank)
221  fprintf(stderr,"Error in NavierStokes2DInitialize(): Parabolic term spatial discretization must be \"%s\"\n",_NC_2STAGE_);
222  return(1);
223  }
224 
225  /* initializing physical model-specific functions */
226 #if defined(HAVE_CUDA) && defined(CUDA_VAR_ORDERDING_AOS)
227  if (solver->use_gpu) {
228  solver->FFunction = gpuNavierStokes2DFlux;
229  solver->SFunction = gpuNavierStokes2DSource;
230  solver->UFunction = gpuNavierStokes2DModifiedSolution;
231  } else {
232 #endif
233  solver->PreStep = NavierStokes2DPreStep;
235  solver->FFunction = NavierStokes2DFlux;
241 #if defined(HAVE_CUDA) && defined(CUDA_VAR_ORDERDING_AOS)
242  }
243 #endif
244 
245 #if defined(HAVE_CUDA) && defined(CUDA_VAR_ORDERDING_AOS)
246  if (solver->use_gpu) {
247  if (!strcmp(solver->SplitHyperbolicFlux,"yes")) {
248  if (!mpi->rank) {
249  fprintf(stderr,"Error in NavierStokes2DInitialize(): Not yet implemented on GPU!");
250  }
251  return 1;
252  } else {
254  if (!strcmp(physics->upw_choice,_RUSANOV_)) solver->Upwind = gpuNavierStokes2DUpwindRusanov;
255  else {
256  if (!mpi->rank) {
257  fprintf(stderr,"Error in NavierStokes2DInitialize(): %s is not implemented on GPU. ",
258  physics->upw_choice);
259  fprintf(stderr,"Only choice is %s.\n",_RUSANOV_);
260  }
261  return 1;
262  }
263  }
264  } else {
265 #endif
266  if (!strcmp(solver->SplitHyperbolicFlux,"yes")) {
270  if (!strcmp(physics->upw_choice,_ROE_)) {
274  } else if (!strcmp(physics->upw_choice,_RF_)) {
275  solver->Upwind = NavierStokes2DUpwindRF;
278  } else if (!strcmp(physics->upw_choice,_LLF_)) {
282  } else if (!strcmp(physics->upw_choice,_RUSANOV_)) {
286  } else {
287  if (!mpi->rank) {
288  fprintf(stderr,"Error in NavierStokes2DInitialize(): %s is not a valid upwinding scheme ",
289  physics->upw_choice);
290  fprintf(stderr,"for use with split hyperbolic flux form. Use %s, %s, %s, or %s.\n",
292  }
293  return(1);
294  }
295  } else {
297  if (!strcmp(physics->upw_choice,_ROE_ )) solver->Upwind = NavierStokes2DUpwindRoe;
298  else if (!strcmp(physics->upw_choice,_RF_ )) solver->Upwind = NavierStokes2DUpwindRF;
299  else if (!strcmp(physics->upw_choice,_LLF_ )) solver->Upwind = NavierStokes2DUpwindLLF;
300  else if (!strcmp(physics->upw_choice,_SWFS_ )) solver->Upwind = NavierStokes2DUpwindSWFS;
301  else if (!strcmp(physics->upw_choice,_RUSANOV_)) solver->Upwind = NavierStokes2DUpwindRusanov;
302  else {
303  if (!mpi->rank) {
304  fprintf(stderr,"Error in NavierStokes2DInitialize(): %s is not a valid upwinding scheme. ",
305  physics->upw_choice);
306  fprintf(stderr,"Choices are %s, %s, %s, %s, and %s.\n",_ROE_,_RF_,_LLF_,_SWFS_,_RUSANOV_);
307  }
308  return(1);
309  }
310  }
311 #if defined(HAVE_CUDA) && defined(CUDA_VAR_ORDERDING_AOS)
312  }
313 #endif
314 
315  /* set the value of gamma in all the boundary objects */
316  int n;
317  DomainBoundary *boundary = (DomainBoundary*) solver->boundary;
318  for (n = 0; n < solver->nBoundaryZones; n++) boundary[n].gamma = physics->gamma;
319 
320  /* hijack the main solver's dissipation function pointer
321  * to this model's own function, since it's difficult to express
322  * the dissipation terms in the general form */
323 #if defined(HAVE_CUDA) && defined(CUDA_VAR_ORDERDING_AOS)
324  if (solver->use_gpu) {
326  } else {
327 #endif
329 #if defined(HAVE_CUDA) && defined(CUDA_VAR_ORDERDING_AOS)
330  }
331 #endif
332 
333  /* allocate array to hold the gravity field */
334  int *dim = solver->dim_local;
335  int ghosts = solver->ghosts;
336  int d, size = 1; for (d=0; d<_MODEL_NDIMS_; d++) size *= (dim[d] + 2*ghosts);
337  physics->grav_field_f = (double*) calloc (size, sizeof(double));
338  physics->grav_field_g = (double*) calloc (size, sizeof(double));
339  /* allocate arrays to hold the fast Jacobian for split form of the hyperbolic flux */
340  physics->fast_jac = (double*) calloc (2*size*_MODEL_NVARS_*_MODEL_NVARS_,sizeof(double));
341  physics->solution = (double*) calloc (size*_MODEL_NVARS_,sizeof(double));
342  /* initialize the gravity fields */
343  IERR NavierStokes2DGravityField(solver,mpi); CHECKERR(ierr);
344 
345 #if defined(HAVE_CUDA) && defined(CUDA_VAR_ORDERDING_AOS)
346  if (solver->use_gpu) gpuNavierStokes2DInitialize(s,m);
347 #endif
348 
349  count++;
350  return(0);
351 }
int NavierStokes2DPreStep(double *, void *, void *, double)
int(* JFunction)(double *, double *, void *, int, int, int)
Definition: hypar.h:326
int NavierStokes2DStiffJacobian(double *, double *, void *, int, int, int)
int NavierStokes2DUpwindRF(double *, double *, double *, double *, double *, double *, int, void *, double)
int(* dFFunction)(double *, double *, int, void *, double)
Definition: hypar.h:280
double(* ComputeCFL)(void *, void *, double, double)
Definition: hypar.h:269
int MPIBroadcast_double(double *, int, int, void *)
Definition: MPIBroadcast.c:9
#define _MODEL_NVARS_
Definition: euler1d.h:58
int NavierStokes2DUpwinddFRusanovModified(double *, double *, double *, double *, double *, double *, int, void *, double)
char spatial_type_par[_MAX_STRING_SIZE_]
Definition: hypar.h:96
void * physics
Definition: hypar.h:266
int(* Upwind)(double *, double *, double *, double *, double *, double *, int, void *, double)
Definition: hypar.h:295
int(* GetRightEigenvectors)(double *, double *, void *, int)
Definition: hypar.h:359
#define _LLF_
Definition: euler1d.h:66
int(* UpwinddF)(double *, double *, double *, double *, double *, double *, int, void *, double)
Definition: hypar.h:300
int * dim_local
Definition: hypar.h:37
int(* UpwindFdF)(double *, double *, double *, double *, double *, double *, int, void *, double)
Definition: hypar.h:307
int NavierStokes2DUpwinddFLLF(double *, double *, double *, double *, double *, double *, int, void *, double)
int(* FdFFunction)(double *, double *, int, void *, double)
Definition: hypar.h:286
int(* FFunction)(double *, double *, int, void *, double)
Definition: hypar.h:276
#define _MODEL_NDIMS_
Definition: euler1d.h:56
int NavierStokes2DParabolicFunction(double *, double *, void *, void *, double)
int NavierStokes2DLeftEigenvectors(double *u, double *L, void *p, int dir)
int NavierStokes2DFlux(double *f, double *u, int dir, void *s, double t)
int(* ParabolicFunction)(double *, double *, void *, void *, double)
Definition: hypar.h:256
int ghosts
Definition: hypar.h:52
#define _SWFS_
Definition: euler1d.h:68
int NavierStokes2DUpwindRusanov(double *, double *, double *, double *, double *, double *, int, void *, double)
int NavierStokes2DUpwindFdFLLF(double *, double *, double *, double *, double *, double *, int, void *, double)
MPI_Comm world
int NavierStokes2DUpwindFdFRF(double *, double *, double *, double *, double *, double *, int, void *, double)
#define _NC_2STAGE_
Definition: hypar.h:477
#define _MAX_STRING_SIZE_
Definition: basic.h:14
int(* AveragingFunction)(double *, double *, double *, void *)
Definition: hypar.h:354
int NavierStokes2DUpwindFdFRoe(double *, double *, double *, double *, double *, double *, int, void *, double)
void * boundary
Definition: hypar.h:159
Structure containing the variables and function pointers defining a boundary.
double * grav_field_f
int gpuNavierStokes2DInitialize(void *s, void *m)
int NavierStokes2DSource(double *, double *, void *, void *, double)
#define _RUSANOV_
Definition: euler1d.h:70
double * grav_field_g
int(* GetLeftEigenvectors)(double *, double *, void *, int)
Definition: hypar.h:357
double NavierStokes2DComputeCFL(void *s, void *m, double dt, double t)
#define _ROE_
Definition: euler1d.h:62
int nBoundaryZones
Definition: hypar.h:157
int NavierStokes2DRightEigenvectors(double *u, double *R, void *p, int dir)
int(* SFunction)(double *, double *, void *, void *, double)
Definition: hypar.h:317
int NavierStokes2DStiffFlux(double *f, double *u, int dir, void *s, double t)
int nvars
Definition: hypar.h:29
int NavierStokes2DUpwindRusanovModified(double *, double *, double *, double *, double *, double *, int, void *, double)
#define CHECKERR(ierr)
Definition: basic.h:18
int NavierStokes2DUpwinddFRoe(double *, double *, double *, double *, double *, double *, int, void *, double)
int NavierStokes2DUpwindRoe(double *, double *, double *, double *, double *, double *, int, void *, double)
char upw_choice[_MAX_STRING_SIZE_]
int(* UFunction)(double *, double *, int, void *, void *, double)
Definition: hypar.h:321
int use_gpu
Definition: hypar.h:449
int MPIBroadcast_character(char *, int, int, void *)
Definition: MPIBroadcast.c:37
int NavierStokes2DRoeAverage(double *uavg, double *uL, double *uR, void *p)
int NavierStokes2DModifiedSolution(double *, double *, int, void *, void *, double)
char SplitHyperbolicFlux[_MAX_STRING_SIZE_]
Definition: hypar.h:92
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.
int ndims
Definition: hypar.h:26
int NavierStokes2DUpwinddFRF(double *, double *, double *, double *, double *, double *, int, void *, double)
int NavierStokes2DJacobian(double *, double *, void *, int, int, int)
#define _RF_
Definition: euler1d.h:64
#define IERR
Definition: basic.h:16
int(* PreStep)(double *, void *, void *, double)
Definition: hypar.h:339
int NavierStokes2DGravityField(void *s, void *m)
Structure of MPI-related variables.
int NavierStokes2DUpwindSWFS(double *, double *, double *, double *, double *, double *, int, void *, double)
int NavierStokes2DNonStiffFlux(double *f, double *u, int dir, void *s, double t)
Structure containing all solver-specific variables and functions.
Definition: hypar.h:23
int NavierStokes2DUpwindFdFRusanovModified(double *, double *, double *, double *, double *, double *, int, void *, double)
int MPIBroadcast_integer(int *, int, int, void *)
Definition: MPIBroadcast.c:23
int NavierStokes2DUpwindLLF(double *, double *, double *, double *, double *, double *, int, void *, double)
int NavierStokes2DCleanup ( void *  s)

Function to clean up all allocations in the 2D Navier Stokes module.

Parameters
sObject of type NavierStokes2D

Definition at line 11 of file NavierStokes2DCleanup.c.

12 {
13  NavierStokes2D *param = (NavierStokes2D*) s;
14 
15  free(param->grav_field_f);
16  free(param->grav_field_g);
17  free(param->fast_jac);
18  free(param->solution);
19  return(0);
20 }
double * grav_field_f
double * grav_field_g
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.