PLUTO  4.0
 All Data Structures Files Functions Variables Enumerations Macros Pages
Functions | Variables
glm.c File Reference

GLM module implementation. More...

#include "pluto.h"

Functions

void GLM_Solve (const State_1D *state, double **VL, double **VR, int beg, int end, Grid *grid)
 
void GLM_Source (const Data_Arr Q, double dt, Grid *grid)
 
void EGLM_Source (const State_1D *state, double dt, int beg, int end, Grid *grid)
 
void GLM_Init (const Data *d, const Time_Step *Dts, Grid *grid)
 

Variables

double glm_ch = -1.0
 

Detailed Description

Contains functions for the GLM module.

Authors
A. Mignone (migno.nosp@m.ne@p.nosp@m.h.uni.nosp@m.to.i.nosp@m.t)
P. Tzeferacos (petro.nosp@m.s.tz.nosp@m.efera.nosp@m.cos@.nosp@m.ph.un.nosp@m.ito..nosp@m.it)
Date
Aug 16, 2012

Function Documentation

void EGLM_Source ( const State_1D state,
double  dt,
int  beg,
int  end,
Grid grid 
)

Add source terms to the right hand side of the conservative equations, momentum and energy equations only. This yields the extended GLM (EGLM) given by Eq. (24a)–(24c) in

"Hyperbolic Divergence cleaning for the MHD Equations" Dedner et al. (2002), JcP, 175, 645

void GLM_Init ( const Data d,
const Time_Step Dts,
Grid grid 
)

Initialize the maximum propagation speed glm_ch.

void GLM_Solve ( const State_1D state,
double **  VL,
double **  VR,
int  beg,
int  end,
Grid grid 
)

Solve the 2x2 linear hyperbolic GLM-MHD system given by the divergence cleaning approach. Build new states VL and VR for Riemann problem. We use Eq. (42) of

Dedner et al. "Hyperbolic Divergence Cleaning for the MHD equations", JCP 175, 645 (2002)

Parameters
[in,out]statepointer to a State_1D structure
[out]VLleft-interface state to be passed to the Riemann solver
[out]VRright-interface state to be passed to the Riemann solver
[in]begstarting index of computation
[in]endfinal index of computation
[in]gridpointer to array of Grid structures

The purpose of this function is two-fold:

  1. assign a unique value to the normal component of magnetic field and to te scalar psi before the actual Riemann solver is called.
  2. compute GLM fluxes in Bn and psi.

The following MAPLE script has been used

restart;
with(linalg);
A := matrix(2,2,[ 0, 1, c^2, 0]);
R := matrix(2,2,[ 1, 1, c, -c]);
L := inverse(R);
E := matrix(2,2,[c,0,0,-c]);
multiply(R,multiply(E,L));
void GLM_Source ( const Data_Arr  Q,
double  dt,
Grid grid 
)

Include the parabolic source term of the Lagrangian multiplier equation in a split fashion for the mixed GLM formulation. Ref. Mignone & Tzeferacos, JCP (2010) 229, 2117, Equation (27).

Variable Documentation

double glm_ch = -1.0

The propagation speed of divergence error.