Carmen Code
 All Classes Files Functions Variables Macros Pages
Functions
PhysicalFluxMHD.cpp File Reference

Computes the MHD physical flux. More...

#include "Carmen.h"
Include dependency graph for PhysicalFluxMHD.cpp:

Functions

Vector FluxX (const Vector &Avg)
 Returns the physical flux of MHD equations in X direction. More...
 
Vector FluxY (const Vector &Avg)
 Returns the physical flux of MHD equations in Y direction. More...
 
Vector FluxZ (const Vector &Avg)
 Returns the physical flux of MHD equations in Z direction. More...
 

Detailed Description

Computes the MHD physical flux.

Author
Anna Karina Fontes Gomes
Version
2.0

Function Documentation

Vector FluxX ( const Vector Avg)

Returns the physical flux of MHD equations in X direction.

Parameters
AvgAverage vector.
Returns
Vector
9 {
10  real rho;
11  real vx, vy, vz;
12  real pre, e;
13  real Bx, By, Bz;
14  real Bx2, By2, Bz2, B2;
15  real vx2, vy2, vz2, v2;
16  real half = 0.5;
17  Vector F(QuantityNb);
18 
19  //Variables
20  rho = Avg.value(1);
21  vx = Avg.value(2)/rho;
22  vy = Avg.value(3)/rho;
23  vz = Avg.value(4)/rho;
24  e = Avg.value(5);
25  Bx = Avg.value(7);
26  By = Avg.value(8);
27  Bz = Avg.value(9);
28 
29  Bx2 = Bx*Bx;
30  By2 = By*By;
31  Bz2 = Bz*Bz;
32  B2 = half*(Bz2+Bx2+By2);
33 
34  vx2 = vx*vx;
35  vy2 = vy*vy;
36  vz2 = vz*vz;
37  v2 = half*(vz2+vx2+vy2);
38 
39  //pressure
40  pre = (Gamma -1.)*(e - rho*v2 - B2);
41 
42  //Physical flux - x-direction
43  F.setValue(1,rho*vx);
44  F.setValue(2,rho*vx2 + pre + half*(Bz2+By2-Bx2));
45  F.setValue(3,rho*vx*vy - Bx*By);
46  F.setValue(4,rho*vx*vz - Bx*Bz);
47  F.setValue(5,(e + pre + B2)*vx - Bx*(vx*Bx + vy*By + vz*Bz) );
48  F.setValue(6,0.0);
49  F.setValue(7,0.0);
50  F.setValue(8,vx*By - vy*Bx);
51  F.setValue(9,vx*Bz - vz*Bx);
52 
53 
54  return F;
55 }
int QuantityNb
Definition: Parameters.cpp:171
Standard class for every vector in Carmen.
Definition: Vector.h:29
real Gamma
Definition: Parameters.cpp:109
real value(const int n) const
Returns the value of the component n.
Definition: Vector.h:565
#define real
Definition: PreProcessor.h:31

Here is the caller graph for this function:

Vector FluxY ( const Vector Avg)

Returns the physical flux of MHD equations in Y direction.

Parameters
AvgAverage vector.
Returns
Vector
59 {
60  real rho;
61  real vx, vy, vz;
62  real pre, e;
63  real Bx, By, Bz;
64  real Bx2, By2, Bz2, B2;
65  real vx2, vy2, vz2, v2;
66  real half = 0.5;
67 
68  Vector G(QuantityNb);
69 
70  //Variables
71  rho = Avg.value(1);
72  vx = Avg.value(2)/rho;
73  vy = Avg.value(3)/rho;
74  vz = Avg.value(4)/rho;
75  e = Avg.value(5);
76  Bx = Avg.value(7);
77  By = Avg.value(8);
78  Bz = Avg.value(9);
79 
80  Bx2 = Bx*Bx;
81  By2 = By*By;
82  Bz2 = Bz*Bz;
83  B2 = half*(Bz2+Bx2+By2);
84 
85  vx2 = vx*vx;
86  vy2 = vy*vy;
87  vz2 = vz*vz;
88  v2 = half*(vz2+vx2+vy2);
89 
90  //pressure
91  pre = (Gamma -1.)*(e - rho*v2 - B2);
92 
93  //Physical flux - y-direction
94  G.setValue(1,rho*vy);
95  G.setValue(2,rho*vx*vy - Bx*By);
96  G.setValue(3,rho*vy2 + pre + half*(Bx2+Bz2-By2));
97  G.setValue(4,rho*vy*vz - By*Bz);
98  G.setValue(5,(e + pre + B2)*vy - By*(vx*Bx + vy*By + vz*Bz));
99  G.setValue(6,0.0);
100  G.setValue(7,vy*Bx - vx*By);
101  G.setValue(8,0.0);
102  G.setValue(9,vy*Bz - vz*By);
103  return G;
104 }
int QuantityNb
Definition: Parameters.cpp:171
Standard class for every vector in Carmen.
Definition: Vector.h:29
real Gamma
Definition: Parameters.cpp:109
real value(const int n) const
Returns the value of the component n.
Definition: Vector.h:565
#define real
Definition: PreProcessor.h:31

Here is the caller graph for this function:

Vector FluxZ ( const Vector Avg)

Returns the physical flux of MHD equations in Z direction.

Parameters
AvgAverage vector.
Returns
Vector
107 {
108  real rho;
109  real vx, vy, vz;
110  real pre, e;
111  real Bx, By, Bz;
112  real Bx2, By2, Bz2, B2;
113  real vx2, vy2, vz2, v2;
114  real half = 0.5;
115 
116  Vector H(QuantityNb);
117 
118  //Variables
119  rho = Avg.value(1);
120  vx = Avg.value(2)/rho;
121  vy = Avg.value(3)/rho;
122  vz = Avg.value(4)/rho;
123  e = Avg.value(5);
124  Bx = Avg.value(7);
125  By = Avg.value(8);
126  Bz = Avg.value(9);
127 
128  Bx2 = Bx*Bx;
129  By2 = By*By;
130  Bz2 = Bz*Bz;
131  B2 = half*(Bz2+Bx2+By2);
132 
133  vx2 = vx*vx;
134  vy2 = vy*vy;
135  vz2 = vz*vz;
136  v2 = half*(vz2+vx2+vy2);
137 
138  //pressure
139  pre = (Gamma -1.)*(e - rho*v2 - B2);
140 
141  //Physical flux - y-direction
142  H.setValue(1,rho*vz);
143  H.setValue(2,rho*vz*vx - Bz*Bx);
144  H.setValue(3,rho*vz*vy - Bz*By);
145  H.setValue(4,rho*vz2 + pre + half*(Bx2+By2-Bz2));
146  H.setValue(5,(e + pre + B2)*vz - Bz*(vx*Bx + vy*By + vz*Bz));
147  H.setValue(6,0.0);
148  H.setValue(7,vz*Bx - vx*Bz);
149  H.setValue(8,vz*By - vy*Bz);
150  H.setValue(9,0.0);
151 
152  return H;
153 }
int QuantityNb
Definition: Parameters.cpp:171
Standard class for every vector in Carmen.
Definition: Vector.h:29
real Gamma
Definition: Parameters.cpp:109
real value(const int n) const
Returns the value of the component n.
Definition: Vector.h:565
#define real
Definition: PreProcessor.h:31

Here is the caller graph for this function: