Pre-calculations¶
preCalc.f90
¶
Description
This module is used to handle some pre-calculations of constants (moment of inertia, mass, volumes), determine the timesteps for I/O and fix boundary values for temperature/entropy and chemical composition
Quick access
- Variables:
- Routines:
Needed modules
iso_fortran_env
(output_unit()
)constants
: module containing constants and parameters used in the code.num_param
: Module containing numerical and control parametersoutput_data
: This module contains the parameters for output controlprecision_mod
: This module controls the precision used in MagICtruncation
(n_r_max()
,l_max()
,minc()
,n_r_ic_max()
,nalias()
,n_cheb_ic_max()
,m_max()
,n_cheb_max()
,m_min()
,lm_max()
,n_phi_max()
,n_theta_max()
): This module defines the grid points and the truncationinit_fields
(bots()
,tops()
,s_bot()
,s_top()
,n_s_bounds()
,l_reset_t()
,topxi()
,botxi()
,xi_bot()
,xi_top()
,n_xi_bounds()
,omega_ma1()
,omegaosz_ma1()
,omega_ic1()
,omegaosz_ic1()
): This module is used to construct the initial solution.parallel_mod
(rank()
,n_procs()
,rank_with_r_lcr()
): This module contains the blocking informationlogic
(l_mag()
,l_cond_ic()
,l_non_rot()
,l_mag_lf()
,l_newmap()
,l_anel()
,l_heat()
,l_anelastic_liquid()
,l_cmb_field()
,l_save_out()
,l_to()
,l_tomovie()
,l_r_field()
,l_movie()
,l_lcr()
,l_dt_cmb_field()
,l_non_adia()
,l_temperature_diff()
,l_chemical_conv()
,l_probe()
,l_precession()
,l_finite_diff()
,l_full_sphere()
): Module containing the logicals that control the runradial_data
(radial_balance()
): This module defines the MPI decomposition in the radial direction.radial_functions
(rscheme_oc()
,temp0()
,r_cmb()
,ogrun()
,r_surface()
,visc()
,or2()
,r()
,r_icb()
,dltemp0()
,beta()
,rho0()
,rgrav()
,dbeta()
,alpha0()
,dentropy0()
,sigma()
,lambda()
,dlkappa()
,kappa()
,dlvisc()
,dllambda()
,divktemp0()
,radial()
,transportproperties()
,l_r()
): This module initiates all the radial functions (transport properties, density, temperature, cheb transforms, etc.)physical_parameters
(nvareps()
,pr()
,prmag()
,ra()
,rascaled()
,ek()
,ekscaled()
,opr()
,opm()
,o_sr()
,radratio()
,sigma_ratio()
,corfac()
,lffac()
,buofac()
,polind()
,nvarcond()
,nvardiff()
,nvarvisc()
,rho_ratio_ic()
,rho_ratio_ma()
,epsc()
,epsc0()
,ktops()
,kbots()
,interior_model()
,r_lcr()
,n_r_lcr()
,mode()
,tmagcon()
,oek()
,bn()
,imagcon()
,ktopxi()
,kbotxi()
,epscxi()
,epscxi0()
,sc()
,osc()
,chemfac()
,raxi()
,po()
,prec_angle()
): Module containing the physical parametershorizontal_data
(horizontal()
): Module containing functions depending on longitude and latitude plus help arrays depending on degree and orderintegration
(rint_r()
): Radial and geostrophic integration functionsuseful
(logwrite()
,abortrun()
): This module contains several useful routines.special
(l_curr()
,fac_loop()
,loopradratio()
,amp_curr()
,le()
,n_imp()
,l_imp()
,l_radial_flow_bc()
,ellipticity_cmb()
,ellip_fac_cmb()
,ellipticity_icb()
,ellip_fac_icb()
,tide_fac20()
,tide_fac22p()
,tide_fac22n()
,amp_tide()
): This module contains all variables for the case of an imposed IC dipole, an imposed external magnetic field and a special boundary forcing to excite inertial modestime_schemes
(type_tscheme()
): This module defines an abstract class type_tscheme which is employed for the time advance of the code.
Variables
- precalculations/get_hit_times [private]¶
Subroutines and functions
- subroutine precalculations/precalc(tscheme)¶
Purpose of this subroutine is to initialize the calc values, arrays, constants that are used all over the code.
- Parameters:
tscheme [real ]
- Called from:
- Call to:
radial()
,transportproperties()
,horizontal()
,rint_r()
,abortrun()
,logwrite()
- subroutine precalculations/precalctimes(time, n_time_step)¶
Precalc. after time, time and dt has been read from startfile.
- Parameters:
time [real ,inout] :: Current time
n_time_step [integer ,inout] :: Index of time loop
- Called from:
- subroutine precalculations/writeinfo(n_out)¶
Purpose of this subroutine is to write the namelist to file unit n_out. This file has to be open before calling this routine.
- Parameters:
n_out [integer ,in] :: Normalized OC surface :’’,ES14.6)’) surf_cmb
- Called from:
- Call to:
radial.f90
¶
Description
This module initiates all the radial functions (transport properties, density, temperature, cheb transforms, etc.)
Quick access
- Variables:
alpha0
,alpha1
,alpha2
,beta
,cheb_ic
,cheb_int
,cheb_int_ic
,cheb_norm_ic
,chebt_ic
,chebt_ic_even
,d2cheb_ic
,d2temp0
,dbeta
,dcheb_ic
,ddbeta
,ddlalpha0
,ddltemp0
,ddlvisc
,dentropy0
,divktemp0
,dlalpha0
,dlkappa
,dllambda
,dltemp0
,dlvisc
,dr_fac_ic
,dr_top_ic
,dxicond
,epscprof
,getbackground
,getentropygradient
,jvarcon
,kappa
,l_r
,lambda
,ndd_costf1
,ndd_costf1_ic
,ndd_costf2_ic
,ndi_costf2_ic
,o_r_ic
,o_r_ic2
,ogrun
,or1
,or2
,or3
,or4
,orho1
,orho2
,otemp1
,polynomialbackground
,r
,r_cmb
,r_ic
,r_icb
,r_surface
,rgrav
,rho0
,sigma
,temp0
,visc
- Routines:
finalize_radial_functions()
,initialize_radial_functions()
,radial()
,transportproperties()
Needed modules
iso_fortran_env
(output_unit()
)truncation
(n_r_max()
,n_cheb_max()
,n_r_ic_max()
,fd_ratio()
,rcut_l()
,fd_stretch()
,fd_order()
,fd_order_bound()
,l_max()
,n_cheb_ic_max()
): This module defines the grid points and the truncationconstants
(sq4pi()
,one()
,two()
,three()
,four()
,half()
,pi()
): module containing constants and parameters used in the code.physical_parameters
: Module containing the physical parameterslogic
(l_mag()
,l_cond_ic()
,l_heat()
,l_anelastic_liquid()
,l_isothermal()
,l_anel()
,l_non_adia()
,l_centrifuge()
,l_temperature_diff()
,l_single_matrix()
,l_var_l()
,l_finite_diff()
,l_newmap()
,l_full_sphere()
,l_chemical_conv()
): Module containing the logicals that control the runradial_data
(nrstart()
,nrstop()
): This module defines the MPI decomposition in the radial direction.cosine_transform_odd
: This module contains the built-in type I discrete Cosine Transforms. This implementation is based on Numerical Recipes and FFTPACK. This only works forn_r_max-1 = 2**a 3**b 5**c
, with a,b,c integers….radial_scheme
(type_rscheme()
): This is an abstract type that defines the radial scheme used in MagICfinite_differences
(type_fd()
): This module is used to calculate the radial grid when finite differences are requestedradial_der
(get_dr()
): Radial derivatives functionsmem_alloc
(bytes_allocated()
): This little module is used to estimate the global memory allocation used in MagICuseful
(logwrite()
,abortrun()
): This module contains several useful routines.parallel_mod
(rank()
): This module contains the blocking informationoutput_data
(tag()
): This module contains the parameters for output controlnum_param
(alph1()
,alph2()
): Module containing numerical and control parameters
Variables
- radial_functions/alpha0 (*) [real,allocatable/public]¶
d log (alpha) / dr
- radial_functions/alpha1 [real,public]¶
Input parameter for non-linear map to define degree of spacing (0.0:2.0)
- radial_functions/alpha2 [real,public]¶
Input parameter for non-linear map to define central point of different spacing (-1.0:1.0)
- radial_functions/beta (*) [real,allocatable/public]¶
Inverse of density scale height drho0/rho0
- radial_functions/cheb_ic (*,*) [real,allocatable/public]¶
Chebyshev polynomials for IC
- radial_functions/cheb_int (*) [real,allocatable/public]¶
Array for cheb integrals
- radial_functions/cheb_int_ic (*) [real,allocatable/public]¶
Array for integrals of cheb for IC
- radial_functions/cheb_norm_ic [real,public]¶
Chebyshev normalisation for IC
- radial_functions/chebt_ic [costf_odd_t,public]¶
- radial_functions/chebt_ic_even [costf_even_t,public]¶
- radial_functions/d2cheb_ic (*,*) [real,allocatable/public]¶
Second radial derivative cheb_ic
- radial_functions/d2temp0 (*) [real,allocatable/private]¶
Second rad. derivative of background temperature
- radial_functions/dbeta (*) [real,allocatable/public]¶
Radial gradient of beta
- radial_functions/dcheb_ic (*,*) [real,allocatable/public]¶
First radial derivative of cheb_ic
- radial_functions/ddbeta (*) [real,allocatable/public]¶
2nd derivative of beta
- radial_functions/ddlalpha0 (*) [real,allocatable/public]¶
\(d/dr(1/alpha d\alpha/dr)\)
- radial_functions/ddltemp0 (*) [real,allocatable/public]¶
\(d/dr(1/T dT/dr)\)
- radial_functions/ddlvisc (*) [real,allocatable/public]¶
2nd derivative of kinematic viscosity
- radial_functions/dentropy0 (*) [real,allocatable/public]¶
Radial gradient of background entropy
- radial_functions/divktemp0 (*) [real,allocatable/public]¶
Term for liquid anelastic approximation
- radial_functions/dlalpha0 (*) [real,allocatable/public]¶
\(1/\alpha d\alpha/dr\)
- radial_functions/dlkappa (*) [real,allocatable/public]¶
Derivative of thermal diffusivity
- radial_functions/dllambda (*) [real,allocatable/public]¶
Derivative of magnetic diffusivity
- radial_functions/dltemp0 (*) [real,allocatable/public]¶
Inverse of temperature scale height
- radial_functions/dlvisc (*) [real,allocatable/public]¶
Derivative of kinematic viscosity
- radial_functions/dr_fac_ic [real,public]¶
For IC: \(2/(2 r_i)\)
- radial_functions/dr_top_ic (*) [real,allocatable/public]¶
Derivative in real space for r=r_i
- radial_functions/dxicond (*) [real,allocatable/public]¶
Radial gradient of chemical composition
- radial_functions/epscprof (*) [real,allocatable/public]¶
Sources in heat equations
- radial_functions/getbackground [private]¶
‘)
- radial_functions/getentropygradient [private]¶
- radial_functions/jvarcon (*) [real,allocatable/public]¶
Analytical solution for toroidal field potential aj (see init_fields.f90)
- radial_functions/kappa (*) [real,allocatable/public]¶
Thermal diffusivity
- radial_functions/l_r (*) [integer,allocatable/public]¶
Variable degree with radius
- radial_functions/lambda (*) [real,allocatable/public]¶
Array of magnetic diffusivity
- radial_functions/ndd_costf1 [integer,public]¶
Radii for transform
- radial_functions/ndd_costf1_ic [integer,private]¶
Radii for transform
- radial_functions/ndd_costf2_ic [integer,private]¶
Radii for transform
- radial_functions/ndi_costf2_ic [integer,private]¶
Radii for transform
- radial_functions/o_r_ic (*) [real,allocatable/public]¶
Inverse of IC radii
- radial_functions/o_r_ic2 (*) [real,allocatable/public]¶
Inverse of square of IC radii
- radial_functions/ogrun (*) [real,allocatable/public]¶
\(1/\Gamma\)
- radial_functions/or1 (*) [real,allocatable/public]¶
\(1/r\)
- radial_functions/or2 (*) [real,allocatable/public]¶
\(1/r^2\)
- radial_functions/or3 (*) [real,allocatable/public]¶
\(1/r^3\)
- radial_functions/or4 (*) [real,allocatable/public]¶
\(1/r^4\)
- radial_functions/orho1 (*) [real,allocatable/public]¶
\(1/\tilde{\rho}\)
- radial_functions/orho2 (*) [real,allocatable/public]¶
\(1/\tilde{\rho}^2\)
- radial_functions/otemp1 (*) [real,allocatable/public]¶
Inverse of background temperature
- radial_functions/polynomialbackground [private]¶
- radial_functions/r (*) [real,allocatable/public]¶
radii
- radial_functions/r_cmb [real,public]¶
OC radius
- radial_functions/r_ic (*) [real,allocatable/public]¶
IC radii
- radial_functions/r_icb [real,public]¶
IC radius
- radial_functions/r_surface [real,public]¶
Surface radius for extrapolation in units of (r_cmb-r_icb)
- radial_functions/rgrav (*) [real,allocatable/public]¶
Buoyancy term dtemp0/Di
- radial_functions/rho0 (*) [real,allocatable/public]¶
Inverse of background density
- radial_functions/sigma (*) [real,allocatable/public]¶
Electrical conductivity
- radial_functions/temp0 (*) [real,allocatable/public]¶
Background temperature
- radial_functions/visc (*) [real,allocatable/public]¶
Kinematic viscosity
Subroutines and functions
- subroutine radial_functions/initialize_radial_functions()¶
Initial memory allocation
- Called from:
- subroutine radial_functions/finalize_radial_functions()¶
Memory deallocation of radial functions
- Called from:
- subroutine radial_functions/radial()¶
Calculates everything needed for radial functions, transforms etc.
- Called from:
- Call to:
- subroutine radial_functions/transportproperties()¶
Calculates the transport properties: electrical conductivity, kinematic viscosity and thermal conductivity.
- Called from:
- Call to:
horizontal.f90
¶
Description
Module containing functions depending on longitude and latitude plus help arrays depending on degree and order
Quick access
- Variables:
cosn_theta_e2
,costheta
,dlh
,dphi
,dpl0eq
,dtheta1a
,dtheta1s
,dtheta2a
,dtheta2s
,dtheta3a
,dtheta3s
,dtheta4a
,dtheta4s
,gauss
,hdif_b
,hdif_s
,hdif_v
,hdif_xi
,n_theta_cal2ord
,n_theta_ord2cal
,o_sin_theta
,o_sin_theta_e2
,phi
,sintheta
,sintheta_e2
,theta_ord
- Routines:
finalize_horizontal_data()
,gauleg()
,horizontal()
,initialize_horizontal_data()
Needed modules
truncation
(l_max()
,n_theta_max()
,n_phi_max()
,nlat_padded()
,lm_max()
,n_m_max()
,minc()
,m_min()
,m_max()
,l_axi()
): This module defines the grid points and the truncationradial_functions
(r_cmb()
): This module initiates all the radial functions (transport properties, density, temperature, cheb transforms, etc.)physical_parameters
(ek()
): Module containing the physical parametersnum_param
(difeta()
,difnu()
,difkap()
,ldif()
,ldifexp()
,difchem()
): Module containing numerical and control parametersblocking
(lm2l()
,lm2m()
): Module containing blocking informationlogic
(l_non_rot()
,l_scramble_theta()
): Module containing the logicals that control the runfft
: This module contains the native subroutines used to compute FFT’s. They are based on the FFT99 package from Temperton: http://www.cesm.ucar.edu/models/cesm1.2/cesm/cesmBbrowser/html_code/cam/fft99.F90.htmlconstants
(pi()
,zero()
,one()
,two()
,half()
): module containing constants and parameters used in the code.precision_mod
: This module controls the precision used in MagICmem_alloc
(bytes_allocated()
): This little module is used to estimate the global memory allocation used in MagIC
Variables
- horizontal_data/cosn_theta_e2 (*) [real,allocatable/public]¶
\(\cos\theta/\sin^2\theta\)
- horizontal_data/costheta (*) [real,allocatable/public]¶
\(\cos\theta\)
- horizontal_data/dlh (*) [real,allocatable/public]¶
\(\ell(\ell+1)\)
- horizontal_data/dphi (*) [complex,allocatable/public]¶
\(\mathrm{i} m\)
- horizontal_data/dpl0eq (*) [real,allocatable/public]¶
- horizontal_data/dtheta1a (*) [real,allocatable/public]¶
- horizontal_data/dtheta1s (*) [real,allocatable/public]¶
- horizontal_data/dtheta2a (*) [real,allocatable/public]¶
- horizontal_data/dtheta2s (*) [real,allocatable/public]¶
- horizontal_data/dtheta3a (*) [real,allocatable/public]¶
- horizontal_data/dtheta3s (*) [real,allocatable/public]¶
- horizontal_data/dtheta4a (*) [real,allocatable/public]¶
- horizontal_data/dtheta4s (*) [real,allocatable/public]¶
- horizontal_data/gauss (*) [real,allocatable/public]¶
- horizontal_data/hdif_b (*) [real,allocatable/public]¶
- horizontal_data/hdif_s (*) [real,allocatable/public]¶
- horizontal_data/hdif_v (*) [real,allocatable/public]¶
- horizontal_data/hdif_xi (*) [real,allocatable/public]¶
- horizontal_data/n_theta_cal2ord (*) [integer,allocatable/public]¶
- horizontal_data/n_theta_ord2cal (*) [integer,allocatable/public]¶
- horizontal_data/o_sin_theta (*) [real,allocatable/public]¶
\(1/\sin\theta\)
- horizontal_data/o_sin_theta_e2 (*) [real,allocatable/public]¶
\(1/\sin^2\theta\)
- horizontal_data/phi (*) [real,allocatable/public]¶
- horizontal_data/sintheta (*) [real,allocatable/public]¶
\(\sin\theta\)
- horizontal_data/sintheta_e2 (*) [real,allocatable/public]¶
\(\sin^2\theta\)
- horizontal_data/theta_ord (*) [real,allocatable/public]¶
Gauss points (unscrambled)
Subroutines and functions
- subroutine horizontal_data/initialize_horizontal_data()¶
Memory allocation of horizontal functions
- Called from:
- subroutine horizontal_data/finalize_horizontal_data()¶
Memory deallocation of horizontal functions
- Called from:
- Call to:
- subroutine horizontal_data/horizontal()¶
Calculates functions of \(\theta\) and \(\phi\), for example the Legendre functions, and functions of degree \(\ell\) and order \(m\) of the legendres.
- Called from:
- Call to:
- subroutine horizontal_data/gauleg(sinthmin, sinthmax, theta_ord, gauss)¶
Subroutine is based on a NR code. Calculates N zeros of legendre polynomial P(l=N) in the interval [sinThMin,sinThMax]. Zeros are returned in radiants theta_ord(i) The respective weights for Gauss-integration are given in gauss(i).
- Parameters:
sinthmin [real ,in] :: lower bound in radiants
sinthmax [real ,in] :: upper bound in radiants
theta_ord (*) [real ,out] :: zeros cos(theta)
gauss (*) [real ,out] :: associated Gauss-Legendre weights
- Called from: