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

Computes the HLL Riemann solver. More...

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

Functions

Vector SchemeHLL (const Cell &Cell1, const Cell &Cell2, const Cell &Cell3, const Cell &Cell4, int AxisNo)
 Returns the HLL numerical flux for MHD equations. The scheme uses four cells to estimate the flux at the interface. Cell2 and Cell3 are the first neighbours on the left and right sides. Cell1 and Cell4 are the second neighbours on the left and right sides. More...
 

Detailed Description

Computes the HLL Riemann solver.

Author
Anna Karina Fontes Gomes
Version
4.0
Date
July-2016

Function Documentation

Vector SchemeHLL ( const Cell Cell1,
const Cell Cell2,
const Cell Cell3,
const Cell Cell4,
const int  AxisNo 
)

Returns the HLL numerical flux for MHD equations. The scheme uses four cells to estimate the flux at the interface. Cell2 and Cell3 are the first neighbours on the left and right sides. Cell1 and Cell4 are the second neighbours on the left and right sides.

Parameters
Cell1second neighbour on the left side
Cell2first neighbour on the left side
Cell3first neighbour on the right side
Cell4second neighbour on the right side
AxisNoAxis of interest.
Returns
Vector
12 {
13 
14  // General variables
15 
16  Vector LeftAverage(QuantityNb); //
17  Vector RightAverage(QuantityNb); // Conservative quantities
18  Vector Result(QuantityNb); // MHD numerical flux
19  int aux=0;
20  // Variables for the HLL scheme
21  Vector FL(QuantityNb), FR(QuantityNb); //Left and right physical fluxes
22  Vector VL(3), VR(3); // Left and right velocities
23  Vector BL(3), BR(3); // Left and right velocities
24  real rhoL=0., rhoR=0.; // Left and right densities
25  real eL=0., eR=0.; // Left and right energies
26  real preL=0., preR=0.;
27  real bkL=0., bkR=0.;
28  real aL=0., aR=0.;
29  real bL=0., bR=0.;
30  real cfL=0., cfR=0.;
31  real SL=0., SR=0.;
32  real dx=0.;
33  dx = Cell2.size(AxisNo);
34  real r, Limit, LeftSlope = 0., RightSlope = 0.; // Left and right slopes
35  int i;
36 
37 // --- Limiter function ---------------------------------------------------------
38 
39  for (i=1; i<=QuantityNb; i++)
40  {
41  // --- Compute left cell-average value ---
42 
43  if (Cell2.average(i) != Cell1.average(i))
44  {
45  RightSlope = Cell3.average(i)-Cell2.average(i);
46  LeftSlope = Cell2.average(i)-Cell1.average(i);
47  r = RightSlope/LeftSlope;
48  Limit = Limiter(r);
49  LeftAverage.setValue(i, Cell2.average(i) + 0.5*Limit*LeftSlope);
50  aux = 1;
51  }
52  else
53  LeftAverage.setValue(i, Cell2.average(i));
54 
55  // --- Compute right cell-average value ---
56 
57  if (Cell3.average(i) != Cell2.average(i))
58  {
59  RightSlope = Cell4.average(i)-Cell3.average(i);
60  LeftSlope = Cell3.average(i)-Cell2.average(i);
61  r = RightSlope/LeftSlope;
62  Limit = Limiter(r);
63  RightAverage.setValue(i, Cell3.average(i) - 0.5*Limit*LeftSlope);
64  aux = 1;
65  }
66  else
67  RightAverage.setValue(i, Cell3.average(i));
68  }
69 
70 
71  // --- HLL scheme -------------------------------------------------------------
72 
73  // --- Conservative variables ---
74 
75  // Left and right densities
76  rhoL = LeftAverage.value(1);
77  rhoR = RightAverage.value(1);
78 
79  // Left and right momentum and magnetic field
80  for (int i=1;i<=3;i++)
81  {
82  VL.setValue( i, LeftAverage.value(i+1));
83  VR.setValue( i, RightAverage.value(i+1));
84  BL.setValue( i, LeftAverage.value(i+6));
85  BR.setValue( i, RightAverage.value(i+6));
86  }
87 
88  // Left and right energies
89  eL = LeftAverage.value(5);
90  eR = RightAverage.value(5);
91 
92  // Left and right pressures
93  preL = (Gamma -1.)*(eL - 0.5*(VL*VL)/rhoL - 0.5*(BL*BL));
94  preR = (Gamma -1.)*(eR - 0.5*(VR*VR)/rhoR - 0.5*(BR*BR));
95 
96  // --- Magnetoacoustic waves calculations --
97 
98  bkL = power2(BL.value(AxisNo))/rhoL;
99  bkR = power2(BR.value(AxisNo))/rhoR;
100 
101  aL = Gamma*preL/rhoL;
102  aR = Gamma*preR/rhoR;
103 
104  bL = (BL*BL)/rhoL;
105  bR = (BR*BR)/rhoR;
106 
107  // Left and Right fast speeds
108  cfL = sqrt(0.5*(aL + bL + sqrt(power2(aL + bL) - 4.0*aL*bkL)));
109  cfR = sqrt(0.5*(aR + bR + sqrt(power2(aR + bR) - 4.0*aR*bkR)));
110 
111  // Left and right slopes
112  SL = Min(Min(VL.value(AxisNo)/rhoL - cfL, VR.value(AxisNo)/rhoR - cfR),0.0);
113  SR = Max(Max(VL.value(AxisNo)/rhoL + cfL, VR.value(AxisNo)/rhoR + cfR),0.0);
114 
115  // --- Physical flux ---
116  if(AxisNo ==1){
117  EigenvalueX = Max(Max(Abs(SL),Abs(SR)),EigenvalueX);
118  FL = FluxX(LeftAverage);
119  FR = FluxX(RightAverage);
120  }else if(AxisNo ==2){
121  EigenvalueY = Max(Max(Abs(SL),Abs(SR)),EigenvalueY);
122  FL = FluxY(LeftAverage);
123  FR = FluxY(RightAverage);
124  }else{
125  EigenvalueZ = Max(Max(Abs(SL),Abs(SR)),EigenvalueZ);
126  FL = FluxZ(LeftAverage);
127  FR = FluxZ(RightAverage);
128  }
129 
130 
131  // --- HLL Riemann Solver ---
132 
133  for(int i=1;i<=QuantityNb;i++)
134  {
135  Result.setValue(i, (SR*FL.value(i) - SL*FR.value(i) + SR*SL*(RightAverage.value(i) - LeftAverage.value(i)))/(SR-SL));
136  }
137 
138  // Parabolic-Hyperbolic divergence Cleaning (Dedner, 2002)
139  fluxCorrection(Result, LeftAverage, RightAverage, AxisNo);
140 
141  // Artificial diffusion terms
142  if(Diffusivity && aux==1) Result = Result - ArtificialViscosity(LeftAverage,RightAverage,dx,AxisNo);
143 
144  return Result;
145 
146 }
#define power2(x)
Definition: Carmen.h:70
real EigenvalueZ
Definition: Parameters.cpp:160
Vector FluxY(const Vector &Avg)
Returns the physical flux of MHD equations in Y direction.
Definition: PhysicalFluxMHD.cpp:58
int QuantityNb
Definition: Parameters.cpp:171
Standard class for every vector in Carmen.
Definition: Vector.h:29
Vector ArtificialViscosity(const Vector &Cell1, const Vector &Cell2, real dx, int AxisNo)
Returns the artificial diffusion source terms in the cell UserCell.
Definition: ArtificialViscosity.cpp:11
#define Min(x, y)
Definition: Carmen.h:62
real Gamma
Definition: Parameters.cpp:109
Vector Limiter(const Vector u, const Vector v)
Returns the value of the slope limiter between the slopes u and v.
Definition: Limiter.cpp:56
bool Diffusivity
Definition: Parameters.cpp:121
Vector FluxX(const Vector &Avg)
Returns the physical flux of MHD equations in X direction.
Definition: PhysicalFluxMHD.cpp:8
real size(const int AxisNo) const
Returns the cell size in the direction AxisNo.
Definition: Cell.h:1095
void fluxCorrection(Vector &Flux, const Vector &AvgL, const Vector &AvgR, int AxisNo)
This function apply the divergence-free correction to the numerical flux.
Definition: FluxCorrection.cpp:9
real average(const int QuantityNo) const
Returns the component no. QuantityNo of the cell-average value.
Definition: Cell.h:1128
Vector FluxZ(const Vector &Avg)
Returns the physical flux of MHD equations in Z direction.
Definition: PhysicalFluxMHD.cpp:106
#define Abs(x)
Definition: Carmen.h:78
real EigenvalueY
Definition: Parameters.cpp:159
real EigenvalueX
Definition: Parameters.cpp:158
#define Max(x, y)
Definition: Carmen.h:54
#define real
Definition: PreProcessor.h:31

Here is the caller graph for this function: