cosmicray_density

Module

Description

$Id: cosmicray_nolog.f90 19193 2012-06-30 12:55:46Z wdobler $

This modules solves the cosmic ray energy density equation. It follows the description of Hanasz & Lesch (2002,2003) as used in their ZEUS 3D implementation.

** AUTOMATIC CPARAM.INC GENERATION ************************ Declare (for generation of cparam.inc) the number of f array variables and auxiliary variables added by this module

CPARAM logical, parameter :: lcosmicray = .true.

MVAR CONTRIBUTION 1 MAUX CONTRIBUTION 0

PENCILS PROVIDED ecr; ecr1; gecr(3); ugecr; bgecr; bglnrho


Quick access

Variables:

ampl_qcr, ampl_qcr2, amplecr, amplecr2, blimiter_cr, cosmicray_diff, ecr_const, ecr_floor, ecr_floor_log, ecr_min, epsilon_ecr, gammacr, gammacr1, idiag_ecrdivum, idiag_ecrm, idiag_ecrmax, idiag_ecrmin, idiag_ecrmz, idiag_ecrph1mz, idiag_ecrph2mz, idiag_ecrph3mz, idiag_ecrpt, idiag_kmax, initecr, initecr2, jcr_param, kx_ecr, ky_ecr, kz_ecr, lalfven_advect, lcheck_negative_ecr, limiter_cr, lnegl, lupw_ecr, luse_diff_constants, lvariable_tensor_diff, radius_ecr, simplified_cosmicray_tensor, widthecr, x_pos_cr, x_pos_cr2, y_pos_cr, y_pos_cr2, z_pos_cr, z_pos_cr2

Routines:

calc_diagnostics_cosmicray(), calc_pencils_cosmicray(), decr_dt(), get_slices_cosmicray(), impose_ecr_floor(), init_ecr(), initialize_cosmicray(), pencil_criteria_cosmicray(), pencil_interdep_cosmicray(), read_cosmicray_init_pars(), read_cosmicray_run_pars(), register_cosmicray(), rprint_cosmicray(), write_cosmicray_init_pars(), write_cosmicray_run_pars()

Needed modules

Variables

  • cosmicray/ampl_qcr [real,optional/default=0.0]
  • cosmicray/ampl_qcr2 [real,optional/default=0.0]
  • cosmicray/amplecr [real,optional/default=0.1]
  • cosmicray/amplecr2 [real,optional/default=0.0]
  • cosmicray/blimiter_cr [real,optional/default=0.0]
  • cosmicray/cosmicray_diff [real,optional/default=0.0]
  • cosmicray/ecr_const [real,optional/default=0.0]
  • cosmicray/ecr_min [real,optional/default=0.0]
  • cosmicray/epsilon_ecr [real,optional/default=0.0]
  • cosmicray/gammacr [real,optional/default=1.3333333333333333]
  • cosmicray/gammacr1 [real]
  • cosmicray/idiag_ecrdivum [integer,optional/default=0]
  • cosmicray/idiag_ecrm [integer,optional/default=0]
  • cosmicray/idiag_ecrmax [integer,optional/default=0]
  • cosmicray/idiag_ecrpt [integer,optional/default=0]
  • cosmicray/idiag_kmax [integer,optional/default=0]
  • cosmicray/initecr [character,optional/default='zero']
  • cosmicray/initecr2 [character,optional/default='zero']
  • cosmicray/jcr_param [real,optional/default=0.0]
  • cosmicray/kpara [real,optional/default=0.0]
  • cosmicray/kperp [real,optional/default=0.0]
  • cosmicray/kx_ecr [real,optional/default=1.0]
  • cosmicray/ky_ecr [real,optional/default=1.0]
  • cosmicray/kz_ecr [real,optional/default=1.0]
  • cosmicray/lalfven_advect [logical,optional/default=.false.]
  • cosmicray/limiter_cr [real,optional/default=1.0]
  • cosmicray/lnegl [logical,optional/default=.false.]
  • cosmicray/luse_diff_constants [logical,optional/default=.false.]
  • cosmicray/lvariable_tensor_diff [logical,optional/default=.false.]
  • cosmicray/radius_ecr [real,optional/default=1.0]
  • cosmicray/simplified_cosmicray_tensor [logical,optional/default=.false.]
  • cosmicray/widthecr [real,optional/default=0.5]
  • cosmicray/x_pos_cr [real,optional/default=0.0]
  • cosmicray/x_pos_cr2 [real,optional/default=0.0]
  • cosmicray/y_pos_cr [real,optional/default=0.0]
  • cosmicray/y_pos_cr2 [real,optional/default=0.0]
  • cosmicray/z_pos_cr [real,optional/default=0.0]
  • cosmicray/z_pos_cr2 [real,optional/default=0.0]

Subroutines and functions

subroutine  cosmicray/register_cosmicray()

Initialise variables which should know that we solve for active scalar: iecr - the cosmic ray energy density; increase nvar accordingly

09-oct-03/tony: coded

Use :

farraymanager

Call to:

svn_id(), blob(), gaussian(), parabola(), wave(), wave_uu(), htube2(), stop_it(), initial_condition_ecr(), dot_mn(), identify_bcs(), fatal_error(), max_mn_name(), save_name()

subroutine  cosmicray/initialize_cosmicray(f)

Perform any necessary post-parameter read initialization Dummy routine

09-oct-03/tony: coded

Parameters:

f (mx,my,mz,mfarray) [real]

Call to:

blob(), gaussian(), parabola(), wave(), wave_uu(), htube2(), stop_it(), initial_condition_ecr(), dot_mn(), identify_bcs(), fatal_error(), max_mn_name(), save_name()

subroutine  cosmicray/init_ecr(f)

initialise cosmic ray energy density field; called from start.f90

09-oct-03/tony: coded

Parameters:

f (mx,my,mz,mfarray) [real]

Use :

mpicomm, sub, initcond, initialcondition (initial_condition_ecr())

Call to:

blob(), gaussian(), parabola(), wave(), wave_uu(), htube2(), stop_it(), initial_condition_ecr(), dot_mn(), identify_bcs(), fatal_error(), max_mn_name(), save_name()

subroutine  cosmicray/pencil_criteria_cosmicray()

The current module’s criteria for which pencils are needed

20-nov-04/anders: coded

Call to:

dot_mn(), identify_bcs(), fatal_error(), max_mn_name(), save_name()

subroutine  cosmicray/pencil_interdep_cosmicray(lpencil_in)

Set all pencils that input pencil array depends on. The most basic pencils should come last (since other pencils chosen here because of dependence may depend on more basic pencils).

20-11-04/anders: coded

Parameters:

lpencil_in (npencils) [logical]

Call to:

dot_mn(), identify_bcs(), fatal_error(), max_mn_name(), save_name()

subroutine  cosmicray/calc_pencils_cosmicray(f, p)

Calculate cosmicray pencils

20-11-04/anders: coded

Parameters:
Use :

sub

Call to:

dot_mn(), identify_bcs(), fatal_error(), max_mn_name(), save_name()

subroutine  cosmicray/decr_dt(f, df, p)

cosmic ray evolution calculate decr/dt + div(u.ecr - flux) = -pcr*divu = -(gammacr-1)*ecr*divu solve as decr/dt + u.grad(ecr) = -gammacr*ecr*divu + div(flux) add du = … - (1/rho)*grad(pcr) to momentum equation Alternatively, use w

09-oct-03/tony: coded

Parameters:
Use :

diagnostics, sub

Call to:

identify_bcs(), fatal_error(), max_mn_name(), save_name()

subroutine  cosmicray/read_cosmicray_init_pars(iomsg)
Parameters:

iomsg [character,out]

Use :

file_io (parallel_unit())

subroutine  cosmicray/write_cosmicray_init_pars(unit)
Parameters:

unit [integer,in]

subroutine  cosmicray/read_cosmicray_run_pars(iomsg)
Parameters:

iomsg [character,out]

Use :

file_io (parallel_unit())

subroutine  cosmicray/write_cosmicray_run_pars(unit)
Parameters:

unit [integer,in]

subroutine  cosmicray/rprint_cosmicray(lreset[, lwrite])

reads and registers print parameters relevant for cosmic rays

6-jul-02/axel: coded

Parameters:
  • lreset [logical]

  • lwrite [logical]

Use :

diagnostics, farraymanager (farray_index_append())

subroutine  cosmicray/get_slices_cosmicray(f, slices)

Write slices for animation of Cosmicray variables.

26-jul-06/tony: coded

Parameters:
Use :

slices_methods (assign_slices_scal())

subroutine  cosmicray/impose_ecr_floor(f)

dummy routine

29-jun-15/axel: adapted from cosmicray.f90

Parameters:

f (mx,my,mz,mfarray) [real,inout]