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

Supersonic outflow boundary conditions (specific to Euler/Navier-Stokes systems) More...

#include <stdlib.h>
#include <basic.h>
#include <arrayfunctions.h>
#include <boundaryconditions.h>
#include <physicalmodels/euler2d.h>
#include <physicalmodels/navierstokes3d.h>

Go to the source code of this file.

Functions

int BCSupersonicOutflowU (void *b, void *m, int ndims, int nvars, int *size, int ghosts, double *phi, double waqt)
 

Detailed Description

Supersonic outflow boundary conditions (specific to Euler/Navier-Stokes systems)

Author
Debojyoti Ghosh

Definition in file BCSupersonicOutflow.c.

Function Documentation

◆ BCSupersonicOutflowU()

int BCSupersonicOutflowU ( void *  b,
void *  m,
int  ndims,
int  nvars,
int *  size,
int  ghosts,
double *  phi,
double  waqt 
)

Applies the supersonic outflow boundary condition: All flow variables (density, pressure, velocity) are extrapolated from the interior since the outflow is supersonic. This boundary condition is specific to two and three dimensional Euler/Navier-Stokes systems (Euler2D, NavierStokes2D, NavierStokes3D).

Note: The extrapolate boundary condition (_EXTRAPOLATE_) can be used as well for this boundary. I am not entirely sure why I wrote the code for this boundary in such a complicated fashion.

Parameters
bBoundary object of type DomainBoundary
mMPI object of type MPIVariables
ndimsNumber of spatial dimensions
nvarsNumber of variables/DoFs per grid point
sizeInteger array with the number of grid points in each spatial dimension
ghostsNumber of ghost points
phiThe solution array on which to apply the boundary condition
waqtCurrent solution time

Definition at line 24 of file BCSupersonicOutflow.c.

34 {
35  DomainBoundary *boundary = (DomainBoundary*) b;
36 
37  int dim = boundary->dim;
38  int face = boundary->face;
39 
40  if (ndims == 2) {
41 
42  /* create a fake physics object */
43  Euler2D physics;
44  double gamma;
45  gamma = physics.gamma = boundary->gamma;
46  double inv_gamma_m1 = 1.0/(gamma-1.0);
47 
48  if (boundary->on_this_proc) {
49  int bounds[ndims], indexb[ndims], indexi[ndims];
50  _ArraySubtract1D_(bounds,boundary->ie,boundary->is,ndims);
51  _ArraySetValue_(indexb,ndims,0);
52  int done = 0;
53  while (!done) {
54  int p1, p2;
55  _ArrayCopy1D_(indexb,indexi,ndims);
56  _ArrayAdd1D_(indexi,indexi,boundary->is,ndims);
57  if (face == 1) indexi[dim] = ghosts-1-indexb[dim];
58  else if (face == -1) indexi[dim] = size[dim]-indexb[dim]-1;
59  else return(1);
60  _ArrayIndex1DWO_(ndims,size,indexb,boundary->is,ghosts,p1);
61  _ArrayIndex1D_(ndims,size,indexi,ghosts,p2);
62 
63  /* flow variables in the interior */
64  double rho, uvel, vvel, energy, pressure;
65  double rho_gpt, uvel_gpt, vvel_gpt, energy_gpt, pressure_gpt;
66  _Euler2DGetFlowVar_((phi+nvars*p2),rho,uvel,vvel,energy,pressure,(&physics));
67  /* set the ghost point values */
68  rho_gpt = rho;
69  pressure_gpt = pressure;
70  uvel_gpt = uvel;
71  vvel_gpt = vvel;
72  energy_gpt = inv_gamma_m1*pressure_gpt
73  + 0.5 * rho_gpt * (uvel_gpt*uvel_gpt + vvel_gpt*vvel_gpt);
74 
75  phi[nvars*p1+0] = rho_gpt;
76  phi[nvars*p1+1] = rho_gpt * uvel_gpt;
77  phi[nvars*p1+2] = rho_gpt * vvel_gpt;
78  phi[nvars*p1+3] = energy_gpt;
79 
80  _ArrayIncrementIndex_(ndims,bounds,indexb,done);
81  }
82  }
83 
84  } else if (ndims == 3) {
85 
86  /* create a fake physics object */
87  double gamma;
88  gamma = boundary->gamma;
89  double inv_gamma_m1 = 1.0/(gamma-1.0);
90 
91  if (boundary->on_this_proc) {
92  int bounds[ndims], indexb[ndims], indexi[ndims];
93  _ArraySubtract1D_(bounds,boundary->ie,boundary->is,ndims);
94  _ArraySetValue_(indexb,ndims,0);
95  int done = 0;
96  while (!done) {
97  int p1, p2;
98  _ArrayCopy1D_(indexb,indexi,ndims);
99  _ArrayAdd1D_(indexi,indexi,boundary->is,ndims);
100  if (face == 1) indexi[dim] = ghosts-1-indexb[dim];
101  else if (face == -1) indexi[dim] = size[dim]-indexb[dim]-1;
102  else return(1);
103  _ArrayIndex1DWO_(ndims,size,indexb,boundary->is,ghosts,p1);
104  _ArrayIndex1D_(ndims,size,indexi,ghosts,p2);
105 
106  /* flow variables in the interior */
107  double rho, uvel, vvel, wvel, energy, pressure;
108  double rho_gpt, uvel_gpt, vvel_gpt, wvel_gpt, energy_gpt, pressure_gpt;
109  _NavierStokes3DGetFlowVar_((phi+nvars*p2),_NavierStokes3D_stride_,rho,uvel,vvel,wvel,energy,pressure,gamma);
110  /* set the ghost point values */
111  rho_gpt = rho;
112  pressure_gpt = pressure;
113  uvel_gpt = uvel;
114  vvel_gpt = vvel;
115  wvel_gpt = wvel;
116  energy_gpt = inv_gamma_m1*pressure_gpt
117  + 0.5 * rho_gpt
118  * (uvel_gpt*uvel_gpt + vvel_gpt*vvel_gpt + wvel_gpt*wvel_gpt);
119 
120  phi[nvars*p1+0] = rho_gpt;
121  phi[nvars*p1+1] = rho_gpt * uvel_gpt;
122  phi[nvars*p1+2] = rho_gpt * vvel_gpt;
123  phi[nvars*p1+3] = rho_gpt * wvel_gpt;
124  phi[nvars*p1+4] = energy_gpt;
125 
126  _ArrayIncrementIndex_(ndims,bounds,indexb,done);
127  }
128  }
129 
130  }
131  return(0);
132 }
Structure containing the variables and function pointers defining a boundary.
#define _Euler2DGetFlowVar_(u, rho, vx, vy, e, P, p)
Definition: euler2d.h:44
#define _ArraySubtract1D_(x, a, b, size)
#define _ArrayAdd1D_(x, a, b, size)
double gamma
Definition: euler2d.h:245
#define _ArrayIndex1D_(N, imax, i, ghost, index)
#define _ArrayIndex1DWO_(N, imax, i, offset, ghost, index)
static const int _NavierStokes3D_stride_
#define _ArraySetValue_(x, size, value)
#define _ArrayIncrementIndex_(N, imax, i, done)
#define _ArrayCopy1D_(x, y, size)
#define _NavierStokes3DGetFlowVar_(u, stride, rho, vx, vy, vz, e, P, gamma)