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

Supersonic inflow 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 BCSupersonicInflowU (void *b, void *m, int ndims, int nvars, int *size, int ghosts, double *phi, double waqt)
 

Detailed Description

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

Author
Debojyoti Ghosh

Definition in file BCSupersonicInflow.c.

Function Documentation

◆ BCSupersonicInflowU()

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

Applies the supersonic (steady) inflow boundary condition: All the flow variables (density, pressure, velocity) are specified at the physical boundary ghost points, since it is supersonic inflow. This boundary condition is specific to two and three dimensional Euler/Navier-Stokes systems (Euler2D, NavierStokes2D, NavierStokes3D).

Note: the Dirichlet boundary condition (_DIRICHLET_) could be used as well for supersonic inflow; however the specified Dirichlet state should be in terms of the conserved variables, while the specified supersonic inflow state here is in terms of the flow variables.

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 BCSupersonicInflow.c.

34 {
35  DomainBoundary *boundary = (DomainBoundary*) b;
36 
37  if (ndims == 2) {
38 
39  double gamma;
40  gamma = boundary->gamma;
41  double inv_gamma_m1 = 1.0/(gamma-1.0);
42 
43  if (boundary->on_this_proc) {
44  int bounds[ndims], indexb[ndims];
45  _ArraySubtract1D_(bounds,boundary->ie,boundary->is,ndims);
46  _ArraySetValue_(indexb,ndims,0);
47  int done = 0;
48  while (!done) {
49  int p1; _ArrayIndex1DWO_(ndims,size,indexb,boundary->is,ghosts,p1);
50 
51  /* set the ghost point values */
52  double rho_gpt, uvel_gpt, vvel_gpt, energy_gpt, pressure_gpt;
53  rho_gpt = boundary->FlowDensity;
54  pressure_gpt = boundary->FlowPressure;
55  uvel_gpt = boundary->FlowVelocity[0];
56  vvel_gpt = boundary->FlowVelocity[1];
57  energy_gpt = inv_gamma_m1*pressure_gpt
58  + 0.5 * rho_gpt * (uvel_gpt*uvel_gpt + vvel_gpt*vvel_gpt);
59 
60  phi[nvars*p1+0] = rho_gpt;
61  phi[nvars*p1+1] = rho_gpt * uvel_gpt;
62  phi[nvars*p1+2] = rho_gpt * vvel_gpt;
63  phi[nvars*p1+3] = energy_gpt;
64 
65  _ArrayIncrementIndex_(ndims,bounds,indexb,done);
66  }
67  }
68 
69  } else if (ndims == 3) {
70 
71  double gamma;
72  gamma = boundary->gamma;
73  double inv_gamma_m1 = 1.0/(gamma-1.0);
74 
75  if (boundary->on_this_proc) {
76  int bounds[ndims], indexb[ndims];
77  _ArraySubtract1D_(bounds,boundary->ie,boundary->is,ndims);
78  _ArraySetValue_(indexb,ndims,0);
79  int done = 0;
80  while (!done) {
81  int p1; _ArrayIndex1DWO_(ndims,size,indexb,boundary->is,ghosts,p1);
82 
83  /* set the ghost point values */
84  double rho_gpt, uvel_gpt, vvel_gpt, wvel_gpt, energy_gpt, pressure_gpt;
85  rho_gpt = boundary->FlowDensity;
86  pressure_gpt = boundary->FlowPressure;
87  uvel_gpt = boundary->FlowVelocity[0];
88  vvel_gpt = boundary->FlowVelocity[1];
89  wvel_gpt = boundary->FlowVelocity[2];
90  energy_gpt = inv_gamma_m1*pressure_gpt
91  + 0.5 * rho_gpt
92  * (uvel_gpt*uvel_gpt + vvel_gpt*vvel_gpt + wvel_gpt*wvel_gpt);
93 
94  phi[nvars*p1+0] = rho_gpt;
95  phi[nvars*p1+1] = rho_gpt * uvel_gpt;
96  phi[nvars*p1+2] = rho_gpt * vvel_gpt;
97  phi[nvars*p1+3] = rho_gpt * wvel_gpt;
98  phi[nvars*p1+4] = energy_gpt;
99 
100  _ArrayIncrementIndex_(ndims,bounds,indexb,done);
101  }
102  }
103 
104  }
105  return(0);
106 }
Structure containing the variables and function pointers defining a boundary.
#define _ArraySubtract1D_(x, a, b, size)
#define _ArrayIndex1DWO_(N, imax, i, offset, ghost, index)
#define _ArraySetValue_(x, size, value)
#define _ArrayIncrementIndex_(N, imax, i, done)