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

Contains the functions compute flux Jacobians for the 1D Euler system. More...

Go to the source code of this file.

Functions

int Euler1DJacobian (double *Jac, double *u, void *p, int dir, int nvars, int upw)
 
int Euler1DStiffJacobian (double *Jac, double *u, void *p, int dir, int nvars, int upw)
 

Detailed Description

Contains the functions compute flux Jacobians for the 1D Euler system.

Author
Debojyoti Ghosh

Definition in file Euler1DJacobian.c.

Function Documentation

◆ Euler1DJacobian()

int Euler1DJacobian ( double *  Jac,
double *  u,
void *  p,
int  dir,
int  nvars,
int  upw 
)

Function to compute the flux Jacobian of the 1D Euler equations, given the solution at a grid point. The Jacobian is square matrix of size nvar=3, and is returned as a 1D array (double) of 9 elements in row-major format.

Parameters
JacJacobian matrix: 1D array of size nvar^2 = 9
usolution at a grid point (array of size nvar = 3)
pobject containing the physics-related parameters
dirdimension (x/y/z) (not used, since this is 1D system)
nvarsnumber of vector components
upw0 -> send back complete Jacobian, 1 -> send back Jacobian of right(+)-moving flux, -1 -> send back Jacobian of left(-)-moving flux

Definition at line 14 of file Euler1DJacobian.c.

22 {
23  Euler1D *param = (Euler1D*) p;
26 
27  /* get the eigenvalues and left,right eigenvectors */
28  _Euler1DEigenvalues_ (u,D,param,0);
29  _Euler1DLeftEigenvectors_ (u,L,param,0);
30  _Euler1DRightEigenvectors_(u,R,param,0);
31 
32  int aupw = absolute(upw), k;
33  k = 0; D[k] = absolute( (1-aupw)*D[k] + 0.5*aupw*(1+upw)*max(0,D[k]) + 0.5*aupw*(1-upw)*min(0,D[k]) );
34  k = 4; D[k] = absolute( (1-aupw)*D[k] + 0.5*aupw*(1+upw)*max(0,D[k]) + 0.5*aupw*(1-upw)*min(0,D[k]) );
35  k = 8; D[k] = absolute( (1-aupw)*D[k] + 0.5*aupw*(1+upw)*max(0,D[k]) + 0.5*aupw*(1-upw)*min(0,D[k]) );
36 
37  MatMult3(3,DL,D,L);
38  MatMult3(3,Jac,R,DL);
39 
40  return(0);
41 }
#define absolute(a)
Definition: math_ops.h:32
Structure containing variables and parameters specific to the 1D Euler equations. This structure cont...
Definition: euler1d.h:273
#define min(a, b)
Definition: math_ops.h:14
#define _Euler1DLeftEigenvectors_(u, L, p, dir)
Definition: euler1d.h:225
#define _Euler1DEigenvalues_(u, D, p, dir)
Definition: euler1d.h:207
#define _Euler1DRightEigenvectors_(u, R, p, dir)
Definition: euler1d.h:249
#define MatMult3(N, A, X, Y)
#define _MODEL_NVARS_
Definition: euler1d.h:58
#define max(a, b)
Definition: math_ops.h:18

◆ Euler1DStiffJacobian()

int Euler1DStiffJacobian ( double *  Jac,
double *  u,
void *  p,
int  dir,
int  nvars,
int  upw 
)

Function to compute the Jacobian of the fast flux (representing the acoustic waves) of the 1D Euler equations, given the solution at a grid point (see _Euler1DSetStiffFlux_, _Euler1DSetStiffJac_). The Jacobian is square matrix of size nvar=3, and is returned as a 1D array (double) of 9 elements in row-major format.

Parameters
JacJacobian matrix: 1D array of size nvar^2 = 9
usolution at a grid point (array of size nvar = 3)
pobject containing the physics-related parameters
dirdimension (x/y/z) (not used, since this is 1D system)
nvarsnumber of vector components
upw0 -> send back complete Jacobian, 1 -> send back Jacobian of right(+)-moving flux, -1 -> send back Jacobian of left(-)-moving flux

Definition at line 48 of file Euler1DJacobian.c.

56 {
57  Euler1D *param = (Euler1D*) p;
60 
61  /* get the eigenvalues and left,right eigenvectors */
62  _Euler1DEigenvalues_ (u,D,param,0);
63  _Euler1DLeftEigenvectors_ (u,L,param,0);
64  _Euler1DRightEigenvectors_(u,R,param,0);
65 
66  int aupw = absolute(upw), k;
67  k = 0; D[k] = 0.0; /* remove the entropy eigenmode */
68  k = 4; D[k] = absolute( (1-aupw)*D[k] + 0.5*aupw*(1+upw)*max(0,D[k]) + 0.5*aupw*(1-upw)*min(0,D[k]) );
69  k = 8; D[k] = absolute( (1-aupw)*D[k] + 0.5*aupw*(1+upw)*max(0,D[k]) + 0.5*aupw*(1-upw)*min(0,D[k]) );
70 
71  MatMult3(3,DL,D,L);
72  MatMult3(3,Jac,R,DL);
73 
74  return(0);
75 }
#define absolute(a)
Definition: math_ops.h:32
Structure containing variables and parameters specific to the 1D Euler equations. This structure cont...
Definition: euler1d.h:273
#define min(a, b)
Definition: math_ops.h:14
#define _Euler1DLeftEigenvectors_(u, L, p, dir)
Definition: euler1d.h:225
#define _Euler1DEigenvalues_(u, D, p, dir)
Definition: euler1d.h:207
#define _Euler1DRightEigenvectors_(u, R, p, dir)
Definition: euler1d.h:249
#define MatMult3(N, A, X, Y)
#define _MODEL_NVARS_
Definition: euler1d.h:58
#define max(a, b)
Definition: math_ops.h:18