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

Thermal slip-wall boundary conditions. More...

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

Go to the source code of this file.

Functions

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

Detailed Description

Thermal slip-wall boundary conditions.

Author
Debojyoti Ghosh

Definition in file BCThermalSlipWall.c.

Function Documentation

◆ BCThermalSlipWallU()

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

Applies the thermal slip-wall boundary condition: This is specific to the 3D Navier-Stokes system (NavierStokes3D). It is used for simulating inviscid walls or symmetric boundaries, where the temperature is specified. The density and the tangential velocity at the ghost points are extrapolated for the interior. The normal velocity at the ghost points is set such that the interpolated velocity at the boundary face is the specified wall velocity. The pressure at the ghost points is set by multiplying the extrapolated density by the specified temperature.

Note: It is assumed that the temperature already contains the gas constant factor, i.e., \( T = P/\rho\).

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 27 of file BCThermalSlipWall.c.

37 {
38  DomainBoundary *boundary = (DomainBoundary*) b;
39 
40  int dim = boundary->dim;
41  int face = boundary->face;
42 
43  if (ndims == 3) {
44 
45  /* create a fake physics object */
46  double gamma;
47  gamma = boundary->gamma;
48  double inv_gamma_m1 = 1.0/(gamma-1.0);
49 
50  if (boundary->on_this_proc) {
51 
52  int *temperature_field_size = boundary->UnsteadyTemperatureSize;
53  int n_time_levels = temperature_field_size[dim];
54  double *time_levels = boundary->UnsteadyTimeLevels;
55  double *temperature_data = boundary->UnsteadyTemperatureData;
56 
57  int it = n_time_levels - 1;
58  while ((time_levels[it] > waqt) && (it > 0)) it--;
59 
60  int bounds[ndims], indexb[ndims], indexi[ndims];
61  _ArraySubtract1D_(bounds,boundary->ie,boundary->is,ndims);
62  _ArraySetValue_(indexb,ndims,0);
63 
64  int done = 0;
65  while (!done) {
66 
67  int p1, p2;
68  _ArrayCopy1D_(indexb,indexi,ndims);
69  _ArrayAdd1D_(indexi,indexi,boundary->is,ndims);
70  if (face == 1) indexi[dim] = ghosts-1-indexb[dim];
71  else if (face == -1) indexi[dim] = size[dim]-indexb[dim]-1;
72  else return(1);
73  _ArrayIndex1DWO_(ndims,size,indexb,boundary->is,ghosts,p1);
74  _ArrayIndex1D_(ndims,size,indexi,ghosts,p2);
75 
76  /* get the specified temperature */
77  int index1[ndims]; _ArrayCopy1D_(indexb,index1,ndims);
78  index1[dim] = it;
79  int q; _ArrayIndex1D_(ndims,temperature_field_size,index1,0,q);
80  double temperature_b = temperature_data[q];
81 
82  /* flow variables in the interior */
83  double rho, uvel, vvel, wvel, energy, pressure;
84  _NavierStokes3DGetFlowVar_((phi+nvars*p2),_NavierStokes3D_stride_,rho,uvel,vvel,wvel,energy,pressure,gamma);
85  /* set the ghost point values */
86  double rho_gpt, uvel_gpt, vvel_gpt, wvel_gpt, energy_gpt, pressure_gpt;
87  rho_gpt = rho;
88  if (dim == _XDIR_) {
89  uvel_gpt = 2.0*boundary->FlowVelocity[_XDIR_] - uvel;
90  vvel_gpt = vvel;
91  wvel_gpt = wvel;
92  } else if (dim == _YDIR_) {
93  uvel_gpt = uvel;
94  vvel_gpt = 2.0*boundary->FlowVelocity[_YDIR_] - vvel;
95  wvel_gpt = wvel;
96  } else if (dim == _ZDIR_) {
97  uvel_gpt = uvel;
98  vvel_gpt = vvel;
99  wvel_gpt = 2.0*boundary->FlowVelocity[_ZDIR_] - wvel;
100  } else {
101  uvel_gpt = 0.0;
102  vvel_gpt = 0.0;
103  wvel_gpt = 0.0;
104  }
105  pressure_gpt = rho_gpt * temperature_b;
106  energy_gpt = inv_gamma_m1*pressure_gpt
107  + 0.5 * rho_gpt
108  * (uvel_gpt*uvel_gpt + vvel_gpt*vvel_gpt + wvel_gpt*wvel_gpt);
109 
110  phi[nvars*p1+0] = rho_gpt;
111  phi[nvars*p1+1] = rho_gpt * uvel_gpt;
112  phi[nvars*p1+2] = rho_gpt * vvel_gpt;
113  phi[nvars*p1+3] = rho_gpt * wvel_gpt;
114  phi[nvars*p1+4] = energy_gpt;
115 
116  _ArrayIncrementIndex_(ndims,bounds,indexb,done);
117  }
118  }
119 
120  }
121  return(0);
122 }
#define _ZDIR_
Structure containing the variables and function pointers defining a boundary.
#define _ArraySubtract1D_(x, a, b, size)
#define _ArrayAdd1D_(x, a, b, size)
#define _YDIR_
Definition: euler2d.h:41
double * UnsteadyTimeLevels
#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 _XDIR_
Definition: euler1d.h:75
double * UnsteadyTemperatureData
#define _ArrayCopy1D_(x, y, size)
#define _NavierStokes3DGetFlowVar_(u, stride, rho, vx, vy, vz, e, P, gamma)