HyPar  1.0
Finite-Difference Hyperbolic-Parabolic PDE Solver on Cartesian Grids
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
NavierStokes3DFlux.c
Go to the documentation of this file.
1 
6 #include <stdlib.h>
7 #include <basic.h>
8 #include <arrayfunctions.h>
9 #include <mathfunctions.h>
10 #include <matmult_native.h>
12 #include <hypar.h>
13 
24  double *f,
25  double *u,
26  int dir,
27  void *s,
28  double t
29  )
30 {
31  HyPar *solver = (HyPar*) s;
32  NavierStokes3D *param = (NavierStokes3D*) solver->physics;
33  int i;
34 
35  int *dim = solver->dim_local;
36  int ghosts = solver->ghosts;
37  int ndims = solver->ndims;
38  int index[ndims], bounds[ndims], offset[ndims];
39 
40  /* set bounds for array index to include ghost points */
41  _ArrayCopy1D_(dim,bounds,ndims);
42  for (i=0; i<ndims; i++) bounds[i] += 2*ghosts;
43 
44  /* set offset such that index is compatible with ghost point arrangement */
45  _ArraySetValue_(offset,ndims,-ghosts);
46 
47  int done = 0; _ArraySetValue_(index,ndims,0);
48  while (!done) {
49  int p; _ArrayIndex1DWO_(ndims,dim,index,offset,ghosts,p);
50  double rho, vx, vy, vz, e, P;
53  _ArrayIncrementIndex_(ndims,bounds,index,done);
54  }
55 
56  return(0);
57 }
58 
71  double *f,
72  double *u,
73  int dir,
74  void *s,
75  double t
76  )
77 {
78  HyPar *solver = (HyPar*) s;
79  NavierStokes3D *param = (NavierStokes3D*) solver->physics;
80  int *dim = solver->dim_local;
81  int ghosts = solver->ghosts;
82  static const int JacSize = _MODEL_NVARS_*_MODEL_NVARS_;
83  static int index[_MODEL_NDIMS_], bounds[_MODEL_NDIMS_], offset[_MODEL_NDIMS_];
84 
85  /* set bounds for array index to include ghost points */
86  _ArrayAddCopy1D_(dim,(2*ghosts),bounds,_MODEL_NDIMS_);
87  /* set offset such that index is compatible with ghost point arrangement */
88  _ArraySetValue_(offset,_MODEL_NDIMS_,-ghosts);
89 
90  int done = 0; _ArraySetValue_(index,_MODEL_NDIMS_,0);
91  while (!done) {
92  int p; _ArrayIndex1DWO_(_MODEL_NDIMS_,dim,index,offset,ghosts,p);
93  double *Af = param->fast_jac+(_MODEL_NDIMS_*p+dir)*JacSize;
95  _ArrayIncrementIndex_(_MODEL_NDIMS_,bounds,index,done);
96  }
97 
98  return(0);
99 }
100 
112  double *f,
113  double *u,
114  int dir,
115  void *s,
116  double t
117  )
118 {
119  HyPar *solver = (HyPar*) s;
120  NavierStokes3D *param = (NavierStokes3D*) solver->physics;
121  int *dim = solver->dim_local;
122  int ghosts = solver->ghosts;
123  static const int JacSize = _MODEL_NVARS_*_MODEL_NVARS_;
124  static int index[_MODEL_NDIMS_], bounds[_MODEL_NDIMS_], offset[_MODEL_NDIMS_];
125  static double ftot[_MODEL_NVARS_], fstiff[_MODEL_NVARS_];
126 
127  /* set bounds for array index to include ghost points */
128  _ArrayAddCopy1D_(dim,(2*ghosts),bounds,_MODEL_NDIMS_);
129  /* set offset such that index is compatible with ghost point arrangement */
130  _ArraySetValue_(offset,_MODEL_NDIMS_,-ghosts);
131 
132  int done = 0; _ArraySetValue_(index,_MODEL_NDIMS_,0);
133  while (!done) {
134  int p; _ArrayIndex1DWO_(_MODEL_NDIMS_,dim,index,offset,ghosts,p);
135  /* compute total flux */
136  double rho, vx, vy, vz, e, P;
137  _NavierStokes3DGetFlowVar_((u+_MODEL_NVARS_*p),_NavierStokes3D_stride_,rho,vx,vy,vz,e,P,param->gamma);
138  _NavierStokes3DSetFlux_(ftot,_NavierStokes3D_stride_,rho,vx,vy,vz,e,P,dir);
139  /* compute stiff stuff */
140  double *Af = param->fast_jac+(_MODEL_NDIMS_*p+dir)*JacSize;
141  MatVecMult5(_MODEL_NVARS_,fstiff,Af,(u+_MODEL_NVARS_*p));
142  /* subtract stiff flux from total flux */
143  _ArraySubtract1D_((f+_MODEL_NVARS_*p),ftot,fstiff,_MODEL_NVARS_);
144  /* Done */
145  _ArrayIncrementIndex_(_MODEL_NDIMS_,bounds,index,done);
146  }
147 
148  return(0);
149 }
#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)
#define _MODEL_NVARS_
Definition: euler1d.h:58
void * physics
Definition: hypar.h:266
int * dim_local
Definition: hypar.h:37
int NavierStokes3DNonStiffFlux(double *f, double *u, int dir, void *s, double t)
#define _NavierStokes3DSetFlux_(f, stride, rho, vx, vy, vz, e, P, dir)
#define _MODEL_NDIMS_
Definition: euler1d.h:56
#define _ArraySubtract1D_(x, a, b, size)
int ghosts
Definition: hypar.h:52
#define MatVecMult5(N, y, A, x)
int NavierStokes3DFlux(double *f, double *u, int dir, void *s, double t)
Contains function definitions for common mathematical functions.
#define _ArrayIndex1DWO_(N, imax, i, offset, ghost, index)
#define _ArrayCopy1D_(x, y, size)
Contains structure definition for hypar.
Structure containing variables and parameters specific to the 3D Navier Stokes equations. This structure contains the physical parameters, variables, and function pointers specific to the 3D Navier-Stokes equations.
Some basic definitions and macros.
int ndims
Definition: hypar.h:26
Contains macros and function definitions for common array operations.
#define _ArrayAddCopy1D_(x, a, y, size)
int NavierStokes3DStiffFlux(double *f, double *u, int dir, void *s, double t)
Structure containing all solver-specific variables and functions.
Definition: hypar.h:23
Contains macros and function definitions for common matrix multiplication.
static const int _NavierStokes3D_stride_