Radial derivatives and integration¶
radial_derivatives.f90
¶
Description
Radial derivatives functions
Quick access
- Variables
- Routines
bulk_to_ghost()
,exch_ghosts()
,finalize_der_arrays()
,get_bound_vals()
,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_complex()
,get_dr_real_1d()
,get_dr_rloc()
,initialize_der_arrays()
Needed modules
constants
(zero()
,one()
,three()
): module containing constants and parameters used in the code.precision_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/
get_dcheb
[public]¶
-
-
radial_der/
get_dr
[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
-
subroutine
radial_der/
get_dcheb_complex
(f, df, n_f_max, n_f_start, n_f_stop, n_r_max, n_cheb_max, d_fac)¶ Returns Chebyshev coeffitients of first derivative df and second derivative ddf for a function whose cheb-coeff. are given as columns in array f(n_f_max,n_r_max).
- Parameters
f (n_f_max,n_r_max) [complex ,in]
df (n_f_max,n_r_max) [complex ,out]
n_f_max [integer ,in,] :: Max no of functions
n_f_start [integer ,in] :: No of function to start with
n_f_stop [integer ,in] :: No of function to stop with
n_r_max [integer ,in,] :: second dimension of f,df,ddf
n_cheb_max [integer ,in] :: Number of cheb modes
d_fac [real ,in] :: factor for interval mapping
-
subroutine
radial_der/
get_dcheb_real_1d
(f, df, n_r_max, n_cheb_max, d_fac)¶
-
subroutine
radial_der/
get_ddcheb
(f, df, ddf, n_f_max, n_f_start, n_f_stop, n_r_max, n_cheb_max, d_fac)¶ Returns Chebyshev coefficients of first derivative df and second derivative ddf for a function whose cheb-coeff. are given as columns in array f(n_c_tot,n_r_max).
- Parameters
f (n_f_max,n_r_max) [complex ,in]
df (n_f_max,n_r_max) [complex ,out]
ddf (n_f_max,n_r_max) [complex ,out]
n_f_max [integer ,in,] :: First dimension of f,df,ddf
n_f_start [integer ,in] :: No of column to start with
n_f_stop [integer ,in] :: No of column to stop with
n_r_max [integer ,in,] :: second dimension of f,df,ddf
n_cheb_max [integer ,in] :: Number of cheb modes
d_fac [real ,in] :: factor for interval mapping
- Called from
-
subroutine
radial_der/
get_dddcheb
(f, df, ddf, dddf, n_f_max, n_f_start, n_f_stop, n_r_max, n_cheb_max, d_fac)¶ Returns chebychev coeffitiens of first derivative df and second derivative ddf for a function whose cheb-coeff. are given as columns in array f(n_c_tot,n_r_max).
- Parameters
f (n_f_max,n_r_max) [complex ,in]
df (n_f_max,n_r_max) [complex ,out]
ddf (n_f_max,n_r_max) [complex ,out]
dddf (n_f_max,n_r_max) [complex ,out]
n_f_max [integer ,in,] :: First dimension of f,df,ddf
n_f_start [integer ,in] :: No of column to start with
n_f_stop [integer ,in] :: No of column to stop with
n_r_max [integer ,in,] :: second dimension of f,df,ddf
n_cheb_max [integer ,in] :: Number of cheb modes
d_fac [real ,in] :: factor for interval mapping
- Called from
-
subroutine
radial_der/
get_dr_real_1d
(f, df, n_r_max, r_scheme)¶
-
subroutine
radial_der/
get_dr_complex
(f, df, n_f_max, n_f_start, n_f_stop, n_r_max, r_scheme[, nocopy[, l_dct_in]])¶ Returns first radial derivative df of the input function f. Array f(n_f_max,*) may contain several functions numbered by the first index. The subroutine calculates the derivaties 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 .
- 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
r_scheme [real ]
- Options
nocopy [logical ,in,]
l_dct_in [logical ,in,]
-
subroutine
radial_der/
get_ddr
(f, df, ddf, n_f_max, n_f_start, n_f_stop, n_r_max, r_scheme[, l_dct_in])¶ 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.
- 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
r_scheme [real ]
- Options
l_dct_in [logical ,in,]
- Called from
get_mag_rhs_imp()
,get_phase_rhs_imp()
,get_entropy_rhs_imp()
,get_pol_rhs_imp()
,assemble_pol()
,assemble_single()
,get_single_rhs_imp()
,get_comp_rhs_imp()
,get_tor_rhs_imp()
- Call to
-
subroutine
radial_der/
get_dddr
(f, df, ddf, dddf, n_f_max, n_f_start, n_f_stop, n_r_max, r_scheme[, l_dct_in])¶ Returns first radial derivative df, the second radial deriv. ddf, and the third radial derivative dddf of the input function f. Array f(n_f_max,*) may contain several functions numbered by the first index. The subroutine calculates the derivaties 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.
- 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
dddf (n_f_max,n_r_max) [complex ,out] :: third 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
r_scheme [real ]
- Options
l_dct_in [logical ,in,]
- Called from
- Call to
-
subroutine
radial_der/
get_dr_rloc
(f_rloc, df_rloc, lm_max, nrstart, nrstop, n_r_max, r_scheme)¶ Purpose of this subroutine is to take the first radial derivative of an input complex array distributed over radius. This can only be used with finite differences.
- Parameters
- Called from
dtvrms()
,dtbrms()
,transp_lmloc_to_rloc()
,transp_rloc_to_lmloc_io()
,finish_exp_mag_rdist()
,finish_exp_entropy_rdist()
,finish_exp_pol_rdist()
,get_pol_rhs_imp_ghost()
,finish_exp_smat_rdist()
,finish_exp_comp_rdist()
- Call to
abortrun()
,get_openmp_blocks()
,exch_ghosts()
,get_bound_vals()
-
subroutine
radial_der/
get_ddr_rloc
(f_rloc, df_rloc, ddf_rloc, lm_max, nrstart, nrstop, n_r_max, r_scheme)¶ Purpose of this subroutine is to take the first and second radial derivatives of an input complex array distributed over radius. This can only be used with finite differences.
- Parameters
- Called from
- Call to
abortrun()
,get_openmp_blocks()
,exch_ghosts()
,get_bound_vals()
-
subroutine
radial_der/
get_ddr_ghost
(f_rloc, df_rloc, ddf_rloc, lm_max, start_lm, stop_lm, nrstart, nrstop, r_scheme)¶ Purpose of this subroutine is to take the first and second radial derivatives of an input complex array distributed over radius that has the ghost zones properly filled.
- Parameters
- Called from
get_mag_rhs_imp_ghost()
,get_phase_rhs_imp_ghost()
,get_entropy_rhs_imp_ghost()
,get_comp_rhs_imp_ghost()
,get_tor_rhs_imp_ghost()
- Call to
-
subroutine
radial_der/
get_ddddr_ghost
(f_rloc, df_rloc, ddf_rloc, dddf_rloc, ddddf_rloc, lm_max, start_lm, stop_lm, nrstart, nrstop, r_scheme)¶ Purpose of this subroutine is to take the first and second radial derivatives of an input complex array distributed over radius that has the ghost zones properly filled.
- Parameters
lm_max [integer ,in,]
start_lm [integer ,in]
stop_lm [integer ,in]
nrstart [integer ,in]
nrstop [integer ,in]
r_scheme [real ]
- Called from
- Call to
-
subroutine
radial_der/
exch_ghosts
(f, lm_max, nrstart, nrstop, nghosts)¶
-
subroutine
radial_der/
get_bound_vals
(fbot, ftop, lm_max, nrstart, nrstop, n_r_max, nbounds)¶ - Parameters
- Called from
-
subroutine
radial_der/
bulk_to_ghost
(x, x_g, ng, nrstart, nrstop, lm_max, start_lm, stop_lm)¶ This subroutine is used to copy an array that is defined from nRstart to nRstop to an array that is defined from nRstart-1 to nRstop+1
radial_derivatives_even.f90
¶
Quick access
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 rarial 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 derivaties 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
- Call to
-
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 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 derivaties 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
- Call to
-
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 rarial 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 derivaties 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
- Call to
-
subroutine
radial_der_even/
get_dcheb_even
(f, df, n_f_max, n_f_start, n_f_stop, n_r_max, n_cheb_max, d_fac)¶ Returns Chebyshev coefficients of first derivative df and second derivative ddf for a function whose cheb-coeff. are given as columns in array f(n_f_max,n_r_max).
- Parameters
f (n_f_max,n_r_max) [complex ,in]
df (n_f_max,n_r_max) [complex ,out]
n_f_max [integer ,in,] :: First dimension of f,df
n_f_start [integer ,in] :: No of function to start with
n_f_stop [integer ,in] :: No of function to stop with
n_r_max [integer ,in,] :: second dimension of f,df
n_cheb_max [integer ,in] :: Number of cheb modes
d_fac [real ,in] :: factor for interval mapping
- Called from
-
subroutine
radial_der_even/
get_ddcheb_even
(f, df, ddf, n_f_max, n_f_start, n_f_stop, n_r_max, n_cheb_max, d_fac)¶ Returns Chebyshev coefficients of first derivative df and second derivative ddf for a function whose cheb-coeff. are given as columns in array f(n_f_max,n_r_max).
- Parameters
f (n_f_max,n_r_max) [complex ,in]
df (n_f_max,n_r_max) [complex ,out]
ddf (n_f_max,n_r_max) [complex ,out]
n_f_max [integer ,in,] :: First dimension of f,df,ddf
n_f_start [integer ,in] :: No of function to start with
n_f_stop [integer ,in] :: No of function to stop with
n_r_max [integer ,in,] :: second dimension of f,df,ddf
n_cheb_max [integer ,in] :: Number of cheb modes
d_fac [real ,in] :: factor for interval mapping
- Called from
integration.f90
¶
Description
Radial integration functions
Quick access
- 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 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
-
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 Chebychev grid points.
- Parameters
f (nrmax) [real ,inout] :: This is zero order contribution
nrmax [integer ,in,optional/default=len(f)] :: 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_force()
,dtbrms()
,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()
,spectrum_scal()
,spectrum_vec()
,get_amplitude()
,updatewp()
- Call to
-
function
integration/
simps
(f, r)¶ Simpson’s method to integrate a function
- Parameters
f (*) [real ,in] :: Input function
r (*) [real ,in] :: Radius
- Return
rint [real ]
- Called from
rint_r()
,initialize_geos()
,outgeos()
,initialize_outto_mod()
,outto()
-
subroutine
integration/
cylmean_otc
(a, v, n_s_max, n_s_otc, r, s, theta[, zdensin])¶ This routine computes a z-averaging using Simpsons’s rule outside T.C.
- Parameters
- Options
zdensin [real ,in,]
- Called from
-
subroutine
integration/
cylmean_itc
(a, vn, vs, n_s_max, n_s_otc, r, s, theta[, zdensin])¶ This routine computes a z-averaging using Simpsons’s rule inside T.C.
- Parameters
- Options
zdensin [real ,in,]
- Called from