cosmicray

Module

Description

$Id$

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.

NB: This module solves for ln(ecr): ecr is here used for lnecr. The alternative module cosmicray_nolog.f90 works with ecr,

and the _nolog version has been more heavily used/developed.

** 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; gecr(3); ugecr


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/amplecr [real,private/optional/default=0.1]
  • cosmicray/amplecr2 [real,private/optional/default=0.0]
  • cosmicray/cosmicray_diff [real,private/optional/default=0.0]
  • cosmicray/ecr_const [real,private/optional/default=1.0]
  • cosmicray/ecr_floor [real,private/optional/default=-1.0]
  • cosmicray/ecr_floor_log [real,private]
  • cosmicray/epsilon_ecr [real,private/optional/default=0.0]
  • cosmicray/gammacr [real,private/optional/default=1.3333333333333333]
  • cosmicray/gammacr1 [real,private]
  • cosmicray/idiag_ecrm [integer,private/optional/default=0]
  • cosmicray/idiag_ecrmax [integer,private/optional/default=0]
  • cosmicray/idiag_kmax [integer,private/optional/default=0]
  • cosmicray/initecr [character,private/optional/default='zero']
  • cosmicray/initecr2 [character,private/optional/default='zero']
  • cosmicray/k_para [real,private/target/optional/default=0.0]
  • cosmicray/k_perp [real,private/target/optional/default=0.0]
  • cosmicray/kx_ecr [real,private/optional/default=1.0]
  • cosmicray/ky_ecr [real,private/optional/default=1.0]
  • cosmicray/kz_ecr [real,private/optional/default=1.0]
  • cosmicray/limiter_cr [real,private/optional/default=1.0]
  • cosmicray/lnegl [logical,private/optional/default=.false.]
  • cosmicray/lupw_ecr [logical,private/optional/default=.false.]
  • cosmicray/luse_diff_constants [logical,private/optional/default=.false.]
  • cosmicray/lvariable_tensor_diff [logical,private/optional/default=.false.]
  • cosmicray/radius_ecr [real,private/optional/default=0.0]
  • cosmicray/simplified_cosmicray_tensor [logical,private/optional/default=.false.]
  • cosmicray/tensor_diffusion [private]
  • cosmicray/vkperp (nx) [real,private]
  • cosmicray/widthecr [real,private/optional/default=0.5]

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, sharedvariables (put_shared_variable())

Call to:

svn_id(), fatal_error(), blob(), gaussian(), parabola(), wave(), wave_uu(), htube2(), initial_condition_ecr(), identify_bcs(), dot2_mn(), calc_diagnostics_cosmicray(), max_mn_name(), multsv_mn(), dot_mn(), g2ij()

subroutine  cosmicray/initialize_cosmicray(f)

Perform any necessary post-parameter read initialization

Parameters:

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

Use :

messages (fatal_error())

Call to:

fatal_error(), blob(), gaussian(), parabola(), wave(), wave_uu(), htube2(), initial_condition_ecr(), identify_bcs(), dot2_mn(), calc_diagnostics_cosmicray(), max_mn_name(), multsv_mn(), dot_mn(), g2ij()

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 :

sub, initcond, initialcondition (initial_condition_ecr())

Call to:

blob(), gaussian(), parabola(), wave(), wave_uu(), htube2(), fatal_error(), initial_condition_ecr(), identify_bcs(), dot2_mn(), calc_diagnostics_cosmicray(), max_mn_name(), multsv_mn(), dot_mn(), g2ij()

subroutine  cosmicray/pencil_criteria_cosmicray()

All pencils that the Cosmicray module depends on are specified here.

20-11-04/anders: coded

Call to:

identify_bcs(), dot2_mn(), calc_diagnostics_cosmicray(), max_mn_name(), multsv_mn(), dot_mn(), g2ij()

subroutine  cosmicray/pencil_interdep_cosmicray(lpencil_in)

Interdependency among pencils provided by the Cosmicray module. is specified here.

20-11-04/anders: coded

Parameters:

lpencil_in (npencils) [logical]

Call to:

identify_bcs(), dot2_mn(), calc_diagnostics_cosmicray(), max_mn_name(), multsv_mn(), dot_mn(), g2ij()

subroutine  cosmicray/calc_pencils_cosmicray(f, p)

Calculate Cosmicray pencils. Most basic pencils should come first, as others may depend on them.

20-11-04/anders: coded

Parameters:
Use :

sub (u_dot_grad(), grad())

Call to:

identify_bcs(), dot2_mn(), calc_diagnostics_cosmicray(), max_mn_name(), multsv_mn(), dot_mn(), g2ij()

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*divu + div(flux(ecr)) + (K grad(ecr)).(grad(ecr))

add du = … - (1/rho)*grad(pcr) to momentum equation

ecr=ecrn

09-oct-03/tony: coded 04-dec-03/snod: modified for lnecr (=ecr)

Parameters:
Use :

sub

Call to:

identify_bcs(), dot2_mn(), calc_diagnostics_cosmicray(), max_mn_name(), multsv_mn(), dot_mn(), g2ij()

subroutine  cosmicray/calc_diagnostics_cosmicray(p)

diagnostics

Parameters:

p [pencil_case]

Use :

diagnostics

Call to:

max_mn_name(), dot2_mn(), multsv_mn(), dot_mn(), g2ij()

subroutine  cosmicray/read_cosmicray_init_pars(iomsg)
Parameters:

iomsg [character,out]

Use :

file_io (parallel_unit())

Call to:

dot2_mn(), multsv_mn(), dot_mn(), g2ij()

subroutine  cosmicray/write_cosmicray_init_pars(unit)
Parameters:

unit [integer,in]

Call to:

dot2_mn(), multsv_mn(), dot_mn(), g2ij()

subroutine  cosmicray/read_cosmicray_run_pars(iomsg)
Parameters:

iomsg [character,out]

Use :

file_io (parallel_unit())

Call to:

dot2_mn(), multsv_mn(), dot_mn(), g2ij()

subroutine  cosmicray/write_cosmicray_run_pars(unit)
Parameters:

unit [integer,in]

Call to:

dot2_mn(), multsv_mn(), dot_mn(), g2ij()

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

Call to:

dot2_mn(), multsv_mn(), dot_mn(), g2ij()

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

Call to:

dot2_mn(), multsv_mn(), dot_mn(), g2ij()

subroutine  cosmicray/impose_ecr_floor(f)

Impose a minimum cosmic ray energy density by setting all lower densities to the minimum value (ecr_floor).

19-may-15/grsarson: adapted from impose_density_floor

Parameters:

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

subroutine  cosmicray/pushpars2c(p_par)
Parameters:

p_par (100) [integer]

Use :

syscalls (copy_addr()), general (string_to_enum())