HyPar  1.0
Finite-Difference Hyperbolic-Parabolic PDE Solver on Cartesian Grids
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
NavierStokes3DPreStep.c
Go to the documentation of this file.
1 
5 #include <arrayfunctions.h>
6 #include <mathfunctions.h>
7 #include <matmult_native.h>
9 #include <hypar.h>
10 
33  double *u,
34  void *s,
35  void *m,
36  double waqt
37  )
38 {
39  HyPar *solver = (HyPar*) s;
40  NavierStokes3D *param = (NavierStokes3D*) solver->physics;
41  int *dim = solver->dim_local;
42  int ghosts = solver->ghosts, dir, p;
43  double *A;
44 
45  static const int ndims = _MODEL_NDIMS_;
46  static const int JacSize = _MODEL_NVARS_*_MODEL_NVARS_;
47  static int index[_MODEL_NDIMS_], bounds[_MODEL_NDIMS_], offset[_MODEL_NDIMS_];
48  static double D[_MODEL_NVARS_*_MODEL_NVARS_],L[_MODEL_NVARS_*_MODEL_NVARS_],
49  R[_MODEL_NVARS_*_MODEL_NVARS_],DL[_MODEL_NVARS_*_MODEL_NVARS_];
50 
51  /* set bounds for array index to include ghost points */
52  _ArrayAddCopy1D_(dim,(2*ghosts),bounds,ndims);
53  /* set offset such that index is compatible with ghost point arrangement */
54  _ArraySetValue_(offset,ndims,-ghosts);
55  /* copy the solution to act as a reference for linearization */
56  _ArrayCopy1D_(u,param->solution,(solver->npoints_local_wghosts*_MODEL_NVARS_));
57 
58  int done = 0; _ArraySetValue_(index,ndims,0);
59  while (!done) {
60  _ArrayIndex1DWO_(ndims,dim,index,offset,ghosts,p);
61  int q = _MODEL_NVARS_*p;
62 
63  dir = _XDIR_;
64  A = (param->fast_jac + _MODEL_NDIMS_*JacSize*p + dir*JacSize);
65  /* get the eigenvalues, and left & right eigenvectors */
69  /* remove the entropy modes (corresponding to eigenvalues u) */
70  D[0] = D[12] = D[18] = 0.0;
71  /* assemble the Jacobian */
72  MatMult5(_MODEL_NVARS_,DL,D,L );
73  MatMult5(_MODEL_NVARS_,A ,R,DL);
74 
75  dir = _YDIR_;
76  A = (param->fast_jac + _MODEL_NDIMS_*JacSize*p + dir*JacSize);
77  /* get the eigenvalues, and left & right eigenvectors */
81  /* remove the entropy modes (corresponding to eigenvalues v) */
82  D[0] = D[6] = D[18] = 0.0;
83  /* assemble the Jacobian */
84  MatMult5(_MODEL_NVARS_,DL,D,L );
85  MatMult5(_MODEL_NVARS_,A ,R,DL);
86 
87  dir = _ZDIR_;
88  A = (param->fast_jac + _MODEL_NDIMS_*JacSize*p + dir*JacSize);
89  /* get the eigenvalues, and left & right eigenvectors */
93  /* remove the entropy modes (corresponding to eigenvalues v) */
94  D[0] = D[6] = D[12] = 0.0;
95  /* assemble the Jacobian */
96  MatMult5(_MODEL_NVARS_,DL,D,L );
97  MatMult5(_MODEL_NVARS_,A ,R,DL);
98 
99  _ArrayIncrementIndex_(ndims,bounds,index,done);
100  }
101  return(0);
102 }
#define _YDIR_
Definition: euler2d.h:41
int npoints_local_wghosts
Definition: hypar.h:42
#define _ArraySetValue_(x, size, value)
3D Navier Stokes equations (compressible flows)
#define _NavierStokes3DLeftEigenvectors_(u, stride, L, ga, dir)
#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
#define _MODEL_NDIMS_
Definition: euler1d.h:56
int ghosts
Definition: hypar.h:52
#define _ZDIR_
#define _NavierStokes3DEigenvalues_(u, stride, D, gamma, dir)
int NavierStokes3DPreStep(double *, void *, void *, double)
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.
#define MatMult5(N, A, X, Y)
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.
#define _NavierStokes3DRightEigenvectors_(u, stride, R, ga, dir)
Contains macros and function definitions for common array operations.
#define _ArrayAddCopy1D_(x, a, y, size)
Structure containing all solver-specific variables and functions.
Definition: hypar.h:23
#define _XDIR_
Definition: euler1d.h:75
Contains macros and function definitions for common matrix multiplication.
static const int _NavierStokes3D_stride_