Radial derivatives and integration¶
radial_derivatives.f90¶
Description
Radial derivatives functions
Quick access
- Variables:
bulk_to_ghost,exch_ghosts,get_dcheb,get_dcheb_complex,get_dcheb_real_1d,get_ddcheb,get_dddcheb,get_ddddr_ghost,get_dddr,get_ddr,get_ddr_ghost,get_ddr_rloc,get_dr,get_dr_real_1d,get_dr_rloc,work_1d_real- Routines:
Needed modules
constants(zero(),one(),three()): module containing constants and parameters used in the code.finite_differences(type_fd()): This module is used to calculate the radial grid when finite differences are requestedprecision_mod: This module controls the precision used in MagICmem_alloc: This little module is used to estimate the global memory allocation used in MagICcosine_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 MagIClogic(l_finite_diff()): Module containing the logicals that control the runparallel_mod: This module contains the blocking informationuseful(abortrun()): This module contains several useful routines.
Variables
- radial_der/bulk_to_ghost [public]¶
- radial_der/exch_ghosts [public]¶
- radial_der/get_dcheb [private]¶
- radial_der/get_dcheb_complex [private]¶
- radial_der/get_dcheb_real_1d [private]¶
- radial_der/get_ddcheb [private]¶
- radial_der/get_dddcheb [private]¶
- radial_der/get_ddddr_ghost [public]¶
- radial_der/get_dddr [public]¶
- radial_der/get_ddr [public]¶
- radial_der/get_ddr_ghost [public]¶
- radial_der/get_ddr_rloc [public]¶
- radial_der/get_dr [public]¶
- radial_der/get_dr_real_1d [private]¶
- radial_der/get_dr_rloc [public]¶
- radial_der/work (*,*) [complex,private/allocatable]¶
- radial_der/work_1d_real (*) [real,private/allocatable]¶
Subroutines and functions
- subroutine radial_der/initialize_der_arrays(n_r_max, llm, ulm)¶
Allocate work arrays to compute derivatives
- Parameters:
n_r_max [integer ,in]
llm [integer ,in]
ulm [integer ,in]
- Called from:
radial_derivatives_even.f90¶
Quick access
- Variables:
- Routines:
Needed modules
constants(zero()): module containing constants and parameters used in the code.precision_mod: This module controls the precision used in MagICcosine_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….
Variables
Subroutines and functions
- subroutine radial_der_even/get_ddr_even(f, df, ddf, n_f_max, n_f_start, n_f_stop, n_r_max, n_cheb_max, dr_fac, work1, work2, chebt_odd, chebt_even)¶
Returns first radial derivative df and second radial derivative ddf of the input function f. Array f(n_f_max,*) may contain several functions numbered by the first index. The subroutine calculates the derivatives of the functions f(n_f_start,*) to f(n_f_stop) by transforming to a Chebychev representation using n_r_max radial grid points. The cheb transforms have to be initialized by calling init_costf1 and init_costf2.
- Parameters:
f (n_f_max,n_r_max) [complex ,in]
df (n_f_max,n_r_max) [complex ,out] :: first derivative of f
ddf (n_f_max,n_r_max) [complex ,out] :: second derivative of f
n_f_max [integer ,in,] :: first dim of f
n_f_start [integer ,in] :: first function to be treated
n_f_stop [integer ,in] :: last function to be treated
n_r_max [integer ,in,] :: number of radial grid points
n_cheb_max [integer ,in] :: number of cheb modes
dr_fac [real ,in] :: mapping factor
work1 (n_f_max,n_r_max) [complex ,out] :: work array needed for costf
work2 (n_f_max,n_r_max) [complex ,out] :: work array needed for costf
chebt_odd [costf_odd_t ,in]
chebt_even [costf_even_t ,in]
- Called from:
- subroutine radial_der_even/get_drns_even(f, df, n_f_max, n_f_start, n_f_stop, n_r_max, n_cheb_max, dr_fac, work1, chebt_odd, chebt_even)¶
Returns first rarial derivative df. Array f(n_f_max,*) may contain several functions numbered by the first index. The subroutine calculates the derivative of the functions f(n_f_start,*) to f(n_f_stop) by transforming to a Chebychev representation using n_r_max radial grid points. The cheb transforms have to be initialized by calling init_costf1 and init_costf2.
- Parameters:
f (n_f_max,n_r_max) [complex ,inout]
df (n_f_max,n_r_max) [complex ,out] :: first derivative of f
n_f_max [integer ,in,] :: first dim of f
n_f_start [integer ,in] :: first function to be treated
n_f_stop [integer ,in] :: last function to be treated
n_r_max [integer ,in,] :: number of radial grid points
n_cheb_max [integer ,in] :: number of cheb modes
dr_fac [real ,in] :: mapping factor
work1 (n_f_max,n_r_max) [complex ,out] :: work array needed for costf
chebt_odd [costf_odd_t ,in]
chebt_even [costf_even_t ,in]
- Called from:
- subroutine radial_der_even/get_ddrns_even(f, df, ddf, n_f_max, n_f_start, n_f_stop, n_r_max, n_cheb_max, dr_fac, work1, chebt_odd, chebt_even)¶
Returns first radial derivative df and second radial derivative ddf of the input function f. Array f(n_f_max,*) may contain several functions numbered by the first index. The subroutine calculates the derivatives of the functions f(n_f_start,*) to f(n_f_stop) by transforming to a Chebychev representation using n_r_max radial grid points. The cheb transforms have to be initialized by calling init_costf1 and init_costf2.
- Parameters:
f (n_f_max,n_r_max) [complex ,inout]
df (n_f_max,n_r_max) [complex ,out] :: first derivative of f
ddf (n_f_max,n_r_max) [complex ,out] :: second derivative of f
n_f_max [integer ,in,] :: first dim of f
n_f_start [integer ,in] :: first function to be treated
n_f_stop [integer ,in] :: last function to be treated
n_r_max [integer ,in,] :: number of radial grid points
n_cheb_max [integer ,in] :: number of cheb modes
dr_fac [real ,in] :: mapping factor
work1 (n_f_max,n_r_max) [complex ,out] :: work array needed for costf
chebt_odd [costf_odd_t ,in]
chebt_even [costf_even_t ,in]
- Called from:
integration.f90¶
Description
Radial and geostrophic integration functions
Quick access
- Variables:
- Routines:
Needed modules
precision_mod: This module controls the precision used in MagICconstants(half(),one(),two(),pi()): module containing constants and parameters used in the code.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 requestedcosine_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….
Variables
Subroutines and functions
- function integration/rintic(f, nrmax, drfac, chebt)¶
This function performs the radial integral over a function f that is given on the appropriate nRmax radial Chebyshev grid points.
- Parameters:
f (nrmax) [real ,inout] :: This is zero order contribution
nrmax [integer ,in,] :: Only even chebs for IC
drfac [real ,in]
chebt [costf_odd_t ,in]
- Return:
rintic [real ]
- Called from:
- function integration/rint_r(f, r, r_scheme)¶
Same as function rInt but for a radial dependent mapping function dr_fac2.
- Parameters:
f (*) [real ,in] :: Input function
r (*) [real ,in] :: Radius
r_scheme [real ] :: Radial scheme (FD or Cheb)
- Return:
rint [real ]
- Called from:
get_poltorrms(),getdlm(),get_e_kin(),get_u_square(),get_e_mag(),outhemi(),outhelicity(),outheat(),outphase(),get_onset(),outperppar(),get_angular_moment(),get_angular_moment_rloc(),output(),get_power(),precalc(),get_amplitude(),updatewp()