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

Routines

get_hit_times(), precalc(), precalctimes(), writeinfo()

Needed modules

Variables

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

magic

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 ,out]

  • n_time_step [integer ,out]

Called from

magic

Call to

get_hit_times()

subroutine precalculations/get_hit_times(t, n_t_max, n_t, l_t, t_start, t_stop, dt, n_tot, n_step, string_bn, time, tscale)

This subroutine checks whether any specific times t(*) are given on input. If so, it returns their number n_r and sets l_t to true. If not, t(*) may also be defined by giving a time step dt or a number n_tot of desired output times and t_stop>t_start.

Parameters
  • t (n_t_max) [real ,inout] :: Times for output

  • n_t_max [integer ,in,optional/default=len(t)] :: Dimension of t(*)

  • n_t [integer ,out] :: No. of output times

  • l_t [logical ,out] :: =.true. if output times are defined

  • t_start [real ,inout] :: Starting time for output

  • t_stop [real ,inout] :: Stop time for output

  • dt [real ,inout] :: Time step for output

  • n_tot [integer ,inout] :: No. of output (times) if no times defined

  • n_step [integer ,inout] :: Ouput step in no. of time steps

  • string_bn [character(len=*),in]

  • time [real ,in] :: Time of start file

  • tscale [real ,in] :: Scale unit for time

Called from

precalctimes()

Call to

abortrun()

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

magic

Call to

abortrun()

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, jvarcon, kappa, l_r, lambda, ndd_costf1, ndd_costf1_ic, ndd_costf2_ic, ndi_costf1_ic, ndi_costf2_ic, o_r_ic, o_r_ic2, ogrun, or1, or2, or3, or4, orho1, orho2, otemp1, r, r_cmb, r_ic, r_icb, r_surface, rgrav, rho0, sigma, temp0, visc

Routines

finalize_radial_functions(), getbackground(), getentropygradient(), initialize_radial_functions(), polynomialbackground(), radial(), transportproperties()

Needed modules

Variables

  • radial_functions/alpha0 (*) [real,allocatable/public]

    Thermal expansion coefficient

  • 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/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,public]

    Radii for transform

  • radial_functions/ndd_costf2_ic [integer,public]

    Radii for transform

  • radial_functions/ndi_costf1_ic [integer,public]

    Radii for transform

  • radial_functions/ndi_costf2_ic [integer,public]

    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/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

magic

subroutine radial_functions/finalize_radial_functions()

Memory deallocation of radial functions

Called from

magic

subroutine radial_functions/radial()

Calculates everything needed for radial functions, transforms etc.

Called from

precalc()

Call to

logwrite(), getentropygradient(), getbackground(), polynomialbackground(), cheb_grid(), get_chebs_even()

subroutine radial_functions/transportproperties()

Calculates the transport properties: electrical conductivity, kinematic viscosity and thermal conductivity.

Called from

precalc()

Call to

abortrun()

subroutine radial_functions/getentropygradient()

This subroutine allows to calculate the background entropy gradient in case stable stratification is required

Called from

radial()

subroutine radial_functions/getbackground(input, boundaryval, output[, coeff])

Linear solver of the form: df/dx +coeff*f= input with f(1)=boundaryVal

Parameters
  • input (n_r_max) [real ,in]

  • boundaryval [real ,in]

  • output (n_r_max) [real ,out]

Options

coeff (n_r_max) [real ,in,]

Called from

radial()

Call to

prepare_mat(), abortrun()

subroutine radial_functions/polynomialbackground(coeffdens, coefftemp[, coeffgrav])

This subroutine allows to calculate a reference state based on an input polynomial function.

Parameters
  • coeffdens (*) [real ,in]

  • coefftemp (*) [real ,in]

Options

coeffgrav (*) [real ,in,]

Called from

radial()

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

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

magic

subroutine horizontal_data/finalize_horizontal_data()

Memory deallocation of horizontal functions

Called from

magic

Call to

finalize_fft()

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

precalc()

Call to

gauleg(), plm_theta(), init_fft()

subroutine horizontal_data/gauleg(sinthmin, sinthmax, theta_ord, gauss, n_th_max)

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 (n_th_max) [real ,out] :: zeros cos(theta)

  • gauss (n_th_max) [real ,out] :: associated Gauss-Legendre weights

  • n_th_max [integer ,in] :: desired maximum degree

Called from

horizontal(), initialize_transforms()