PLUTO  4.4-patch2
Functions | Variables
glm.c File Reference

GLM module implementation. More...

#include "pluto.h"

Functions

void GLM_Solve (const Sweep *sweep, int beg, int end, Grid *grid)
 
void GLM_Source (const Data *data, double dt, Grid *grid)
 
void GLM_ExtendedSource (const Sweep *sweep, double dt, int beg, int end, Grid *grid)
 
void GLM_Init (const Data *d, const timeStep *Dts, Grid *grid)
 
void GLM_ComputeDivE (const Sweep *sweep, Grid *grid)
 

Variables

double glm_ch = -1.0
 

Detailed Description

Contains functions for the GLM module.

Authors
A. Mignone (migno.nosp@m.ne@t.nosp@m.o.inf.nosp@m.n.it)
P. Tzeferacos (petro.nosp@m.s.tz.nosp@m.efera.nosp@m.cos@.nosp@m.ph.un.nosp@m.ito..nosp@m.it)
Date
Sep 03, 2020

Reference

Function Documentation

◆ GLM_ComputeDivE()

void GLM_ComputeDivE ( const Sweep sweep,
Grid grid 
)

Compute the divergence of E using Godunov fluxes previously obtained at cell interfaces. This function may be used in Resistive RMHD.

◆ GLM_ExtendedSource()

void GLM_ExtendedSource ( const Sweep sweep,
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 equations given by Eq. (24a)–(24c) in

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

◆ GLM_Init()

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

Initialize the maximum propagation speed glm_ch at the beginning of integration cycle.

◆ GLM_Solve()

void GLM_Solve ( const Sweep sweep,
int  beg,
int  end,
Grid grid 
)

Solve the 2x2 linear hyperbolic GLM-MHD system given by the divergence cleaning approach. Modify inteface states (Bx and psi components) for input to full Riemann problem. We use Eq. (42) of Dedner et al (2002)

Parameters
[in,out]sweeppointer to a Sweep structure
[in]begstarting index of computation
[in]endfinal index of computation
[in]gridpointer to Grid structure

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));

◆ GLM_Source()

void GLM_Source ( const Data data,
double  dt,
Grid grid 
)

Include the damping 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

◆ glm_ch

double glm_ch = -1.0

The propagation speed of divergence error.