HyPar  1.0
Finite-Difference Hyperbolic-Parabolic PDE Solver on Cartesian Grids
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
BCNoslipWall.c
Go to the documentation of this file.
1 
5 #include <stdlib.h>
6 #include <basic.h>
7 #include <arrayfunctions.h>
8 #include <boundaryconditions.h>
9 
12 
21  void *b,
22  void *m,
23  int ndims,
24  int nvars,
25  int *size,
26  int ghosts,
27  double *phi,
28  double waqt
29  )
30 {
31  DomainBoundary *boundary = (DomainBoundary*) b;
32 
33  int dim = boundary->dim;
34  int face = boundary->face;
35 
36  if (ndims == 2) {
37 
38  /* create a fake physics object */
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], indexi[ndims];
45  _ArraySubtract1D_(bounds,boundary->ie,boundary->is,ndims);
46  _ArraySetValue_(indexb,ndims,0);
47  int done = 0;
48  while (!done) {
49  int p1, p2;
50  _ArrayCopy1D_(indexb,indexi,ndims);
51  _ArrayAdd1D_(indexi,indexi,boundary->is,ndims);
52  if (face == 1) indexi[dim] = ghosts-1-indexb[dim];
53  else if (face == -1) indexi[dim] = size[dim]-indexb[dim]-1;
54  else return(1);
55  _ArrayIndex1DWO_(ndims,size,indexb,boundary->is,ghosts,p1);
56  _ArrayIndex1D_(ndims,size,indexi,ghosts,p2);
57 
58  /* flow variables in the interior */
59  double rho, uvel, vvel, energy, pressure;
60  double rho_gpt, uvel_gpt, vvel_gpt, energy_gpt, pressure_gpt;
61  _NavierStokes2DGetFlowVar_((phi+nvars*p2),rho,uvel,vvel,energy,pressure,gamma);
62  /* set the ghost point values */
63  rho_gpt = rho;
64  pressure_gpt = pressure;
65  uvel_gpt = 2.0*boundary->FlowVelocity[0] - uvel;
66  vvel_gpt = 2.0*boundary->FlowVelocity[1] - vvel;
67  energy_gpt = inv_gamma_m1*pressure_gpt
68  + 0.5 * rho_gpt * (uvel_gpt*uvel_gpt + vvel_gpt*vvel_gpt);
69 
70  phi[nvars*p1+0] = rho_gpt;
71  phi[nvars*p1+1] = rho_gpt * uvel_gpt;
72  phi[nvars*p1+2] = rho_gpt * vvel_gpt;
73  phi[nvars*p1+3] = energy_gpt;
74 
75  _ArrayIncrementIndex_(ndims,bounds,indexb,done);
76  }
77  }
78 
79  } else if (ndims == 3) {
80 
81  double gamma;
82  gamma = boundary->gamma;
83  double inv_gamma_m1 = 1.0/(gamma-1.0);
84 
85  if (boundary->on_this_proc) {
86  int bounds[ndims], indexb[ndims], indexi[ndims];
87  _ArraySubtract1D_(bounds,boundary->ie,boundary->is,ndims);
88  _ArraySetValue_(indexb,ndims,0);
89  int done = 0;
90  while (!done) {
91  int p1, p2;
92  _ArrayCopy1D_(indexb,indexi,ndims);
93  _ArrayAdd1D_(indexi,indexi,boundary->is,ndims);
94  if (face == 1) indexi[dim] = ghosts-1-indexb[dim];
95  else if (face == -1) indexi[dim] = size[dim]-indexb[dim]-1;
96  else return(1);
97  _ArrayIndex1DWO_(ndims,size,indexb,boundary->is,ghosts,p1);
98  _ArrayIndex1D_(ndims,size,indexi,ghosts,p2);
99 
100  /* flow variables in the interior */
101  double rho, uvel, vvel, wvel, energy, pressure;
102  double rho_gpt, uvel_gpt, vvel_gpt, wvel_gpt, energy_gpt, pressure_gpt;
103  _NavierStokes3DGetFlowVar_((phi+nvars*p2),_NavierStokes3D_stride_,rho,uvel,vvel,wvel,energy,pressure,gamma);
104  /* set the ghost point values */
105  rho_gpt = rho;
106  pressure_gpt = pressure;
107  uvel_gpt = 2.0*boundary->FlowVelocity[0] - uvel;
108  vvel_gpt = 2.0*boundary->FlowVelocity[1] - vvel;
109  wvel_gpt = 2.0*boundary->FlowVelocity[2] - wvel;
110  energy_gpt = inv_gamma_m1*pressure_gpt
111  + 0.5 * rho_gpt
112  * (uvel_gpt*uvel_gpt + vvel_gpt*vvel_gpt + wvel_gpt*wvel_gpt);
113 
114  phi[nvars*p1+0] = rho_gpt;
115  phi[nvars*p1+1] = rho_gpt * uvel_gpt;
116  phi[nvars*p1+2] = rho_gpt * vvel_gpt;
117  phi[nvars*p1+3] = rho_gpt * wvel_gpt;
118  phi[nvars*p1+4] = energy_gpt;
119 
120  _ArrayIncrementIndex_(ndims,bounds,indexb,done);
121  }
122  }
123 
124  }
125  return(0);
126 }
#define _NavierStokes3DGetFlowVar_(u, stride, rho, vx, vy, vz, e, P, gamma)
#define _ArraySetValue_(x, size, value)
3D Navier Stokes equations (compressible flows)
#define _ArrayIncrementIndex_(N, imax, i, done)
int BCNoslipWallU(void *, void *, int, int, int *, int, double *, double)
Definition: BCNoslipWall.c:20
Containts the structures and definitions for boundary condition implementation.
#define _ArraySubtract1D_(x, a, b, size)
#define _ArrayIndex1D_(N, imax, i, ghost, index)
2D Navier Stokes equations (compressible flows)
Structure containing the variables and function pointers defining a boundary.
#define _ArrayIndex1DWO_(N, imax, i, offset, ghost, index)
#define _ArrayCopy1D_(x, y, size)
#define _NavierStokes2DGetFlowVar_(u, rho, vx, vy, e, P, gamma)
Some basic definitions and macros.
Contains macros and function definitions for common array operations.
#define _ArrayAdd1D_(x, a, b, size)
static const int _NavierStokes3D_stride_