chiral

Module

Description

$Id$

This modules solves two reactive scalar advection equations This is used for modeling the spatial evolution of left and right handed aminoacids.

** 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 :: lchiral = .true.

MVAR CONTRIBUTION 2 MAUX CONTRIBUTION 0

** AUTOMATIC REFERENCE-LINK.TEX GENERATION **************** Declare relevant citations from pencil-code/doc/citations/ref.bib for this module. The entries are taken from pencil-code/doc/citations/notes.tex

2025PhRvD.111d3541V,%Vachaspati+Brandenburg “Spectra of magnetic fields from electroweak symmetry breaking” 2004IJAsB…3..209B,% Brandenburg & Multamaki “How long can left and right handed life forms coexist?”


Quick access

Variables:

amplxx_chiral, amplyy_chiral, amplzz_chiral, chiral_crossinhibition, chiral_diff, chiral_diffxx, chiral_diffzz, chiral_fidelity, chiral_fisherh, chiral_fisherk, chiral_fishermu, chiral_fishernu, chiral_fisherr, chiral_fisherr2, chiral_fisherr2_tdep, chiral_fisherr2_tend, chiral_fisherr2_tstart, chiral_list, chiral_reaction, cutoffyy_chiral, del2xx_chiral, del2yy_chiral, del2zz_chiral, enum_chiral_reaction, gradx0, grady0, gxx_chiral, gyy_chiral, gzz_chiral, idiag_bmaxep, idiag_brmsep, idiag_jbmep, idiag_jmaxep, idiag_jrmsep, idiag_qq21m_chiral, idiag_qq21qqm_chiral, idiag_qqm_chiral, idiag_r2tdep, idiag_xx_chiralm, idiag_xx_chiralmax, idiag_yy_chiralm, idiag_yy_chiralmax, idiag_zz_chiralm, idiag_zz_chiralmax, initpoweryy_chiral, initxx_chiral, inityy_chiral, initzz_chiral, ixx2_chiral, ixx_chiral, iyy2_chiral, iyy_chiral, izz_chiral, kx_xx_chiral, kx_yy_chiral, ky_xx_chiral, ky_yy_chiral, kz_xx_chiral, kz_yy_chiral, limposed_gradient, linitialize_aa_from_ep, linitialize_aa_from_ep_alpgrad, linitialize_aa_from_ep_betgrad, llorentzforceep, lupw_chiral, lzz_chiral, qq, radiusxx_chiral, radiusyy_chiral, tinitialize_aa_from_ep, widthxx_chiral, widthyy_chiral, widthzz_chiral, xposxx_chiral, xposyy_chiral, yposxx_chiral, yposyy_chiral, zposxx_chiral, zposyy_chiral

Routines:

calc_diagnostics_chiral(), calc_pencils_chiral(), chiral_before_boundary(), dxy_chiral_dt(), get_slices_chiral(), init_chiral(), initialize_chiral(), pencil_criteria_chiral(), pencil_interdep_chiral(), read_chiral_init_pars(), read_chiral_run_pars(), register_chiral(), rprint_chiral(), write_chiral_init_pars(), write_chiral_run_pars()

Needed modules

Variables

  • chiral/amplxx_chiral [real,private/optional/default=0.1]
  • chiral/amplyy_chiral [real,private/optional/default=0.1]
  • chiral/amplzz_chiral [real,private/optional/default=0.1]
  • chiral/chiral_crossinhibition [real,private/optional/default=1.0]
  • chiral/chiral_diff [real,private/optional/default=0.0]
  • chiral/chiral_diffxx [real,private/optional/default=impossible]
  • chiral/chiral_diffzz [real,private/optional/default=impossible]
  • chiral/chiral_fidelity [real,private/optional/default=1.0]
  • chiral/chiral_fisherh [real,private/optional/default=1.0]
  • chiral/chiral_fisherk [real,private/optional/default=1.0]
  • chiral/chiral_fishermu [real,private/optional/default=0.0]
  • chiral/chiral_fishernu [real,private/optional/default=0.0]
  • chiral/chiral_fisherr [real,private/optional/default=0.0]
  • chiral/chiral_fisherr2 [real,private/optional/default=0.0]
  • chiral/chiral_fisherr2_tdep [real,private]
  • chiral/chiral_fisherr2_tend [real,private/optional/default=0.0]
  • chiral/chiral_fisherr2_tstart [real,private/optional/default=0.0]
  • chiral/chiral_list [private]
  • chiral/chiral_reaction [character,private/optional/default='bahn_model']
  • chiral/cutoffyy_chiral [real,private/optional/default=0.0]
  • chiral/del2xx_chiral (nx) [real,private]
  • chiral/del2yy_chiral (nx) [real,private]
  • chiral/del2zz_chiral (nx) [real,private]
  • chiral/enum_chiral_reaction [integer,private/optional/default=0]
  • chiral/gradx0 (3) [real,private/optional/default=(/0.0,0.0,0.0/)]
  • chiral/grady0 (3) [real,private/optional/default=(/0.0,0.0,0.0/)]
  • chiral/gxx_chiral (nx,3) [real,private]
  • chiral/gyy_chiral (nx,3) [real,private]
  • chiral/gzz_chiral (nx,3) [real,private]
  • chiral/idiag_bmaxep [integer,private/optional/default=0]
  • chiral/idiag_brmsep [integer,private/optional/default=0]
  • chiral/idiag_jbmep [integer,private/optional/default=0]
  • chiral/idiag_jmaxep [integer,private/optional/default=0]
  • chiral/idiag_jrmsep [integer,private/optional/default=0]
  • chiral/idiag_qq21m_chiral [integer,private/optional/default=0]
  • chiral/idiag_qq21qqm_chiral [integer,private/optional/default=0]
  • chiral/idiag_qqm_chiral [integer,private/optional/default=0]
  • chiral/idiag_r2tdep [integer,private/optional/default=0]
  • chiral/idiag_xx_chiralm [integer,private/optional/default=0]
  • chiral/idiag_xx_chiralmax [integer,private/optional/default=0]
  • chiral/idiag_yy_chiralm [integer,private/optional/default=0]
  • chiral/idiag_yy_chiralmax [integer,private/optional/default=0]
  • chiral/idiag_zz_chiralm [integer,private/optional/default=0]
  • chiral/idiag_zz_chiralmax [integer,private/optional/default=0]
  • chiral/initpoweryy_chiral [real,private/optional/default=2.0]
  • chiral/initxx_chiral [character,private/optional/default='zero']
  • chiral/inityy_chiral [character,private/optional/default='zero']
  • chiral/initzz_chiral [character,private/optional/default='zero']
  • chiral/ixx2_chiral [integer,public/optional/default=0]
  • chiral/ixx_chiral [integer,public/optional/default=0]
  • chiral/iyy2_chiral [integer,public/optional/default=0]
  • chiral/iyy_chiral [integer,public/optional/default=0]
  • chiral/izz_chiral [integer,private/optional/default=0]
  • chiral/kx_xx_chiral [real,private/optional/default=1.0]
  • chiral/kx_yy_chiral [real,private/optional/default=1.0]
  • chiral/ky_xx_chiral [real,private/optional/default=1.0]
  • chiral/ky_yy_chiral [real,private/optional/default=1.0]
  • chiral/kz_xx_chiral [real,private/optional/default=1.0]
  • chiral/kz_yy_chiral [real,private/optional/default=1.0]
  • chiral/limposed_gradient [logical,private/optional/default=.false.]
  • chiral/linitialize_aa_from_ep [logical,private/optional/default=.false.]
  • chiral/linitialize_aa_from_ep_alpgrad [logical,private/optional/default=.false.]
  • chiral/linitialize_aa_from_ep_betgrad [logical,private/optional/default=.false.]
  • chiral/llorentzforceep [logical,private/optional/default=.false.]
  • chiral/lupw_chiral [logical,private/optional/default=.false.]
  • chiral/lzz_chiral [logical,private/optional/default=.false.]
  • chiral/pp [real,private]
  • chiral/qq [real,private]
  • chiral/radiusxx_chiral [real,private/optional/default=0.0]
  • chiral/radiusyy_chiral [real,private/optional/default=0.0]
  • chiral/tinitialize_aa_from_ep [real,private/optional/default=0.0]
  • chiral/widthxx_chiral [real,private/optional/default=0.5]
  • chiral/widthyy_chiral [real,private/optional/default=0.5]
  • chiral/widthzz_chiral [real,private/optional/default=0.5]
  • chiral/xposxx_chiral [real,private/optional/default=0.0]
  • chiral/xposyy_chiral [real,private/optional/default=0.0]
  • chiral/yposxx_chiral [real,private/optional/default=0.0]
  • chiral/yposyy_chiral [real,private/optional/default=0.0]
  • chiral/zposxx_chiral [real,private/optional/default=0.0]
  • chiral/zposyy_chiral [real,private/optional/default=0.0]

Subroutines and functions

subroutine  chiral/register_chiral()

Initialise variables which should know that we solve for passive scalar: iXX_chiral and iYY_chiral; increase nvar accordingly

28-may-04/axel: adapted from pscalar

Use :

farraymanager

Call to:

svn_id(), blob(), hat(), gaussian(), wave(), cosxz_cosz(), cosy_sinz(), cosx_cosz(), cosx_cosy_cosz(), cosx_siny_cosz(), fatal_error(), cosyz_sinz(), cosy_cosz(), power_randomphase(), initial_condition_chiral(), identify_bcs(), multsv_mn(), calc_diagnostics_chiral(), max_mn_name(), save_name(), g2ij(), multmv_mn()

subroutine  chiral/initialize_chiral(f)

Perform any necessary post-parameter read initialization Dummy routine

28-may-04/axel: adapted from pscalar

Parameters:

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

Call to:

blob(), hat(), gaussian(), wave(), cosxz_cosz(), cosy_sinz(), cosx_cosz(), cosx_cosy_cosz(), cosx_siny_cosz(), fatal_error(), cosyz_sinz(), cosy_cosz(), power_randomphase(), initial_condition_chiral(), identify_bcs(), multsv_mn(), calc_diagnostics_chiral(), max_mn_name(), save_name(), g2ij(), multmv_mn()

subroutine  chiral/init_chiral(f)

initialise passive scalar field; called from start.f90

28-may-04/axel: adapted from pscalar

Parameters:

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

Use :

sub, general, initcond, initialcondition (initial_condition_chiral())

Call to:

blob(), hat(), gaussian(), wave(), cosxz_cosz(), cosy_sinz(), cosx_cosz(), cosx_cosy_cosz(), cosx_siny_cosz(), fatal_error(), cosyz_sinz(), cosy_cosz(), power_randomphase(), initial_condition_chiral(), identify_bcs(), multsv_mn(), calc_diagnostics_chiral(), max_mn_name(), save_name(), g2ij(), multmv_mn()

subroutine  chiral/pencil_criteria_chiral()

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

21-nov-04/anders: coded

Call to:

identify_bcs(), multsv_mn(), fatal_error(), calc_diagnostics_chiral(), max_mn_name(), save_name(), g2ij(), multmv_mn()

subroutine  chiral/pencil_interdep_chiral(lpencil_in)

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

21-11-04/anders: coded

Parameters:

lpencil_in (npencils) [logical]

Call to:

identify_bcs(), multsv_mn(), fatal_error(), calc_diagnostics_chiral(), max_mn_name(), save_name(), g2ij(), multmv_mn()

subroutine  chiral/calc_pencils_chiral(f, p)

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

21-11-04/anders: coded

Parameters:
Call to:

identify_bcs(), multsv_mn(), fatal_error(), calc_diagnostics_chiral(), max_mn_name(), save_name(), g2ij(), multmv_mn()

subroutine  chiral/dxy_chiral_dt(f, df, p)

passive scalar evolution calculate chirality equations in reduced form; see q-bio/0401036

28-may-04/axel: adapted from pscalar

1-jul-09/axel: included gradX, gradY, and allowed for different reactions

Parameters:
Use :

sub

Call to:

identify_bcs(), multsv_mn(), fatal_error(), calc_diagnostics_chiral(), max_mn_name(), save_name(), g2ij(), multmv_mn()

subroutine  chiral/calc_diagnostics_chiral(f, p)

diagnostics

output for double and triple correlators (assume z-gradient of cc) <u_k u_j d_j c> = <u_k c uu.gradXX_chiral>

Parameters:
Use :

diagnostics, sub

Call to:

max_mn_name(), save_name(), g2ij(), multmv_mn()

subroutine  chiral/chiral_before_boundary(f)

initialize aa from Euler potentials (EP), but only once when t=0

4-jul-09/axel: coded

Parameters:

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

Use :

sub

subroutine  chiral/read_chiral_init_pars(iomsg)
Parameters:

iomsg [character,out]

Use :

file_io (parallel_unit())

subroutine  chiral/write_chiral_init_pars(unit)
Parameters:

unit [integer,in]

subroutine  chiral/read_chiral_run_pars(iomsg)
Parameters:

iomsg [character,out]

Use :

file_io (parallel_unit())

subroutine  chiral/write_chiral_run_pars(unit)
Parameters:

unit [integer,in]

subroutine  chiral/rprint_chiral(lreset[, lwrite])

reads and registers print parameters relevant for chirality

28-may-04/axel: adapted from pscalar

Parameters:
  • lreset [logical]

  • lwrite [logical]

Use :

diagnostics

subroutine  chiral/get_slices_chiral(f, slices)

Write slices for animation of Chiral variables.

26-jul-06/tony: coded

Parameters:
Use :

slices_methods (assign_slices_scal())

subroutine  chiral/pushpars2c(p_par)
Parameters:

p_par (30) [integer]

Use :

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