Radial scheme¶
radial_scheme.f90
¶
Description
This is an abstract type that defines the radial scheme used in MagIC
Quick access
- Types:
- Variables:
Needed modules
precision_mod
: This module controls the precision used in MagIC
Types
- type radial_scheme/type_rscheme¶
- Type fields:
% boundary_fac [real ]
% d2rmat (*,*) [real ,allocatable]
% d3rmat (*,*) [real ,allocatable]
% d4rmat (*,*) [real ,allocatable]
% ddddr (*,*) [real ,allocatable]
% dddr (*,*) [real ,allocatable]
% ddr (*,*) [real ,allocatable]
% ddr_bot (*,*) [real ,allocatable]
% ddr_top (*,*) [real ,allocatable]
% dr (*,*) [real ,allocatable]
% dr_bot (*,*) [real ,allocatable]
% dr_top (*,*) [real ,allocatable]
% drmat (*,*) [real ,allocatable]
% drx (*) [real ,allocatable]
% n_max [integer ]
% nrmax [integer ]
% order [integer ]
% order_boundary [integer ]
% rmat (*,*) [real ,allocatable]
% rnorm [real ]
% version [character(len=72)]
Variables
chebyshev.f90
¶
Quick access
- Types:
- Variables:
Needed modules
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 MagICconstants
(half()
,one()
,two()
,three()
,four()
,pi()
,ci()
): module containing constants and parameters used in the code.blocking
(llm()
,ulm()
): Module containing blocking informationradial_scheme
(type_rscheme()
): This is an abstract type that defines the radial scheme used in MagICuseful
(factorise()
): This module contains several useful routines.cosine_transform_odd
(costf_odd_t()
): 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….num_param
(map_function()
): Module containing numerical and control parameters
Types
- type chebyshev/type_cheb_odd¶
- Type fields:
% alpha1 [real ]
% alpha2 [real ]
% chebt_oc [costf_odd_t ]
% dddrx (*) [real ,allocatable]
% ddrx (*) [real ,allocatable]
% l_map [logical ]
% work_costf (*,*) [complex ,pointer]
% x_cheb (*) [real ,allocatable]
Variables
finite_differences.f90
¶
Description
This module is used to calculate the radial grid when finite differences are requested
Quick access
- Types:
- Variables:
get_der_mat
,get_fd_coeffs
,get_fd_grid
,populate_fd_weights
,robin_bc
Needed modules
precision_mod
: This module controls the precision used in MagICconstants
(zero()
,one()
,two()
,half()
): module containing constants and parameters used in the code.useful
(logwrite()
,abortrun()
): This module contains several useful routines.mem_alloc
(bytes_allocated()
): This little module is used to estimate the global memory allocation used in MagICradial_scheme
(type_rscheme()
): This is an abstract type that defines the radial scheme used in MagIC
Types
- type finite_differences/type_fd¶
- Type fields:
% ddddr_bot (*,*) [real ,allocatable]
% ddddr_top (*,*) [real ,allocatable]
% dddr_bot (*,*) [real ,allocatable]
% dddr_top (*,*) [real ,allocatable]
Variables
Chebyshev polynomials and cosine transforms¶
chebyshev_polynoms.f90
¶
Quick access
- Routines:
Needed modules
precision_mod
: This module controls the precision used in MagICconstants
(pi()
,half()
,one()
,two()
,four()
): module containing constants and parameters used in the code.num_param
(map_function()
): Module containing numerical and control parameters
Variables
Subroutines and functions
- subroutine chebyshev_polynoms_mod/get_chebs_even(n_r, a, b, y, n_r_max, cheb, dcheb, d2cheb, dim1, dim2)¶
Construct even Chebyshev polynomials and their first, second and third derivative up to degree
2*(n_r/2)
at(n_r/2)
points x in the interval \([a,b]\). Since the polynoms are only defined in \([-1,1]\) we have to map the points x in [a,b] onto points y in the interval \([-1,1]\). For even Chebs we need only half the point of the map, these(n_r/2)
points are in the interval \([1,0[\).- Parameters:
n_r [integer ,in] :: only even chebs stored !
a [real ,in]
b [real ,in] :: interval boundaries \([a,b]\)
y (n_r_max) [real ,in] :: n_r grid points in interval \([a,b]\)
n_r_max [integer ,in,] :: max number of radial points, dims of y
cheb (dim1,dim2) [real ,out] ::
cheb(i,j)
is Chebyshev pol.dcheb (dim1,dim2) [real ,out] :: first derivative of cheb
d2cheb (dim1,dim2) [real ,out] :: second derivative o cheb
dim1 [integer ,in]
dim2 [integer ,in] :: dimensions of cheb,dcheb,……
- Called from:
- subroutine chebyshev_polynoms_mod/cheb_grid(a, b, n, x, y, a1, a2, x0, lbd, l_map)¶
Given the interval \([a,b]\) the routine returns the n+1 points that should be used to support a Chebyshev expansion. These are the n+1 extrema
y(i)
of the Chebyshev polynomial of degree n in the interval \([-1,1]\). The respective points mapped into the interval of question \([a,b]\) are thex(i)
.Note
x(i)
andy(i)
are stored in the reversed order:x(1)=b
,x(n+1)=a
,y(1)=1
,y(n+1)=-1
- Parameters:
a [real ,in]
b [real ,in] :: interval boundaries
n [integer ,in] :: degree of Cheb polynomial to be represented by the grid points
x (*) [real ,out] :: grid points in interval \([a,b]\)
y (*) [real ,out] :: grid points in interval \([-1,1]\)
a1 [real ,in]
a2 [real ,in]
x0 [real ,in]
lbd [real ,in]
l_map [logical ,in] :: Chebyshev mapping
- Called from:
cosine_transform_odd.f90
¶
Description
This module contains the built-in type I discrete Cosine Transforms. This
implementation is based on Numerical Recipes and FFTPACK. This only works
for n_r_max-1 = 2**a 3**b 5**c
, with a,b,c integers.
Quick access
- Types:
- Variables:
Needed modules
iso_fortran_env
(output_unit()
)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 MagICconstants
(half()
,one()
,two()
,pi()
,sin36()
,cos36()
,sin60()
,sin72()
,cos72()
): module containing constants and parameters used in the code.useful
(factorise()
,abortrun()
): This module contains several useful routines.
Types
- type cosine_transform_odd/costf_odd_t¶
- Type fields:
% d_costf_init (*) [real ,allocatable]
% i_costf_init (*) [integer ,allocatable]
Variables
cosine_transform_even.f90
¶
Quick access
- Types:
- Variables:
Needed modules
iso_fortran_env
(output_unit()
)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 MagICtruncation
(lm_max()
): This module defines the grid points and the truncationconstants
(half()
,one()
,two()
,pi()
,sin36()
,cos36()
,sin60()
,sin72()
,cos72()
): module containing constants and parameters used in the code.useful
(factorise()
,abortrun()
): This module contains several useful routines.
Types
- type cosine_transform_even/costf_even_t¶
- Type fields:
% d_costf_init (*) [real ,allocatable]
% i_costf_init (*) [integer ,allocatable]
Variables
fft_fac.f90
¶
Quick access
- Routines:
Needed modules
precision_mod
: This module controls the precision used in MagICconstants
(sin36()
,sin60()
,sin72()
,cos36()
,cos72()
,half()
): module containing constants and parameters used in the code.
Variables
Subroutines and functions
- subroutine fft_fac_mod/fft_fac_real(a, b, c, d, trigs, nv, l1, l2, n, ifac, la)¶
Main part of Fourier / Chebyshev transform called in
costf1
,costf2
- Parameters:
a (*) [real ,in]
b (*) [real ,in]
c (*) [real ,out]
d (*) [real ,out]
trigs (2 * n) [real ,in]
nv [integer ,in]
l1 [integer ,in]
l2 [integer ,in]
n [integer ,in,]
ifac [integer ,in]
la [integer ,in]
- subroutine fft_fac_mod/fft_fac_complex(a, b, c, d, trigs, nv, l1, l2, n, ifac, la)¶
Main part of Fourier / Chebyshev transform called in
costf1
,costf2
- Parameters:
a (*) [complex ,in]
b (*) [complex ,in]
c (*) [complex ,out]
d (*) [complex ,out]
trigs (2 * n) [real ,in]
nv [integer ,in]
l1 [integer ,in]
l2 [integer ,in]
n [integer ,in,]
ifac [integer ,in]
la [integer ,in]
Fourier transforms¶
fft.f90
¶
Description
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.html I simply got rid of the ‘go to’ and Fortran legacy statements Those transforms only work for prime decomposition that involve factors of 2, 3, 5
Quick access
- Variables:
d_fft_init
,fft991
,fft99a
,fft99b
,i_fft_init
,nd
,ni
,vpassm
- Routines:
fft_many()
,fft_to_real()
,finalize_fft()
,ifft_many()
,init_fft()
Needed modules
iso_fortran_env
(output_unit()
)precision_mod
: This module controls the precision used in MagICuseful
(factorise()
,abortrun()
): This module contains several useful routines.constants
(pi()
,sin36()
,sin60()
,sin72()
,cos36()
,cos72()
,one()
,two()
,half()
): module containing constants and parameters used in the code.truncation
(n_phi_max()
,nlat_padded()
): This module defines the grid points and the truncationparallel_mod
: This module contains the blocking information
Variables
- fft/d_fft_init (*) [real,private/allocatable]¶
- fft/fft991 [private]¶
- fft/fft99a [private]¶
- fft/fft99b [private]¶
- fft/i_fft_init (100) [integer,private]¶
- fft/nd [integer,private]¶
- fft/ni [integer,private/parameter/optional/default=100]¶
- fft/vpassm [private]¶
Subroutines and functions
- subroutine fft/init_fft(n)¶
Purpose of this subroutine is to calculate and store several values that will be needed for a fast Fourier transforms.
- Parameters:
n [integer ,in] :: Dimension of problem, number of grid points
- Called from:
- Call to:
- subroutine fft/finalize_fft()¶
Memory deallocation of FFT help arrays
- Called from:
- subroutine fft/fft_to_real(f, ld_f, nrep)¶
- Parameters:
f (ld_f,nrep) [real ,inout]
ld_f [integer ,in,]
nrep [integer ,in,]
- subroutine fft/fft_many(g, f)¶
Fourier transform:
f(nlat,nlon) -> fhat(nlat,nlon/2+1)
- Parameters:
g (*,*) [real ,in]
f (nlat_padded,1 + n_phi_max / 2) [complex ,out]
- Called from:
- Call to:
- subroutine fft/ifft_many(f, g)¶
Inverse Fourier transform:
fhat(nlat, nlon/2+1) -> f(nlat,nlon)
- Parameters:
f (*,*) [complex ,in]
g (nlat_padded,n_phi_max) [real ,out]
- Called from:
native_qst_to_spat()
,native_sphtor_to_spat()
,native_sph_to_spat()
,native_sph_to_grad_spat()
- Call to:
Spherical harmonic transforms¶
plms.f90
¶
Quick access
- Routines:
Needed modules
precision_mod
: This module controls the precision used in MagICconstants
(osq4pi()
,one()
,two()
): module containing constants and parameters used in the code.
Variables
Subroutines and functions
- subroutine plms_theta/plm_theta(theta, max_degree, min_order, max_order, m0, plma, dtheta_plma, norm)¶
This produces the \(P_{\ell}^m\) for all degrees \(\ell\) and orders \(m\) for a given colatitude. \(\theta\), as well as \(\sin \theta d P_{\ell}^m /d\theta\). The \(P_{\ell}^m\) are stored with a single lm index stored with m first.
Several normalisation are supported:
n=0 – surface normalised,
n=1 – Schmidt normalised,
n=2 – fully normalised.
- Parameters:
theta [real ,in] :: angle in radians
max_degree [integer ,in] :: required max degree of plm
min_order [integer ,in] :: required min order of plm
max_order [integer ,in] :: required max order of plm
m0 [integer ,in] :: basic wave number
plma (*) [real ,out] :: associated Legendre polynomials at theta
dtheta_plma (*) [real ,out] :: their theta derivative
norm [integer ,in] :: =0 fully normalised
- Called from:
horizontal()
,initialize_magnetic_energy()
,initialize_transforms()
shtransforms.f90
¶
Description
This module is used when the native built-in SH transforms of MagIC are used. Those are much slower than SHTns, and are not recommanded for production runs!
Quick access
- Variables:
- Routines:
finalize_transforms()
,initialize_transforms()
,native_axi_to_spat()
,native_qst_to_spat()
,native_spat_to_sph()
,native_spat_to_sph_tor()
,native_sph_to_grad_spat()
,native_sph_to_spat()
,native_sphtor_to_spat()
,native_toraxi_to_spat()
Needed modules
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 MagICtruncation
(lm_max()
,n_m_max()
,l_max()
,l_axi()
,n_theta_max()
,minc()
,n_phi_max()
,m_max()
,m_min()
): This module defines the grid points and the truncationhorizontal_data
(gauleg()
,o_sin_theta_e2()
): Module containing functions depending on longitude and latitude plus help arrays depending on degree and orderconstants
(zero()
,half()
,one()
,ci()
,pi()
,two()
): module containing constants and parameters used in the code.fft
(fft_many()
,ifft_many()
): 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.htmlparallel_mod
(get_openmp_blocks()
): This module contains the blocking information
Variables
- shtransforms/d_mc2m (*) [real,allocatable/public]¶
- shtransforms/dplm (*,*) [real,allocatable]¶
- shtransforms/lstart (*) [integer,allocatable]¶
- shtransforms/lstop (*) [integer,allocatable]¶
- shtransforms/plm (*,*) [real,allocatable]¶
Legendre polynomials \(P_{\ell m}\)
- shtransforms/wdplm (*,*) [real,allocatable]¶
- shtransforms/wplm (*,*) [real,allocatable]¶
Subroutines and functions
- subroutine shtransforms/initialize_transforms()¶
This subroutine allocates the arrays needed when the native SH transforms are used. This also defines the Legendre polynomials and the Gauss weights.
- Called from:
- Call to:
- subroutine shtransforms/finalize_transforms()¶
This subroutine handles the memory deallocation of arrays used with the native SH transforms.
- Called from:
- subroutine shtransforms/native_qst_to_spat(qlm, slm, tlm, brc, btc, bpc, lcut)¶
Vector spherical harmonic transform: take Q,S,T and transform them to a vector field
- Parameters:
qlm (lm_max) [complex ,in] :: Poloidal
slm (lm_max) [complex ,in] :: Spheroidal
tlm (lm_max) [complex ,in] :: Toroidal
brc (n_theta_max,n_phi_max) [real ,out]
btc (n_theta_max,n_phi_max) [real ,out]
bpc (n_theta_max,n_phi_max) [real ,out]
lcut [integer ,in]
- Called from:
torpol_to_spat()
,torpol_to_curl_spat_ic()
,torpol_to_spat_ic()
,torpol_to_curl_spat()
- Call to:
- subroutine shtransforms/native_sphtor_to_spat(slm, tlm, btc, bpc, lcut)¶
Use spheroidal and toroidal potentials to transform them to angular vector components btheta and bphi
- Parameters:
slm (lm_max) [complex ,in]
tlm (lm_max) [complex ,in]
btc (n_theta_max,n_phi_max) [real ,out]
bpc (n_theta_max,n_phi_max) [real ,out]
lcut [integer ,in]
- Called from:
- Call to:
- subroutine shtransforms/native_axi_to_spat(slm, sc)¶
- Parameters:
slm (1 + l_max) [complex ,in]
sc (n_theta_max) [real ,out]
- Called from:
- subroutine shtransforms/native_toraxi_to_spat(tlm, btc, bpc)¶
Use spheroidal and toroidal potentials to transform them to angular vector components btheta and bphi
- Parameters:
tlm (1 + l_max) [complex ,in]
btc (n_theta_max) [real ,out]
bpc (n_theta_max) [real ,out]
- Called from:
- subroutine shtransforms/native_sph_to_spat(slm, sc, lcut)¶
Spherical Harmonic Transform for a scalar input field
- Parameters:
slm (lm_max) [complex ,in]
sc (n_theta_max,n_phi_max) [real ,out]
lcut [integer ,in]
- Called from:
- Call to:
- subroutine shtransforms/native_sph_to_grad_spat(slm, gradtc, gradpc, lcut)¶
Transform
s(l)
intodsdt(theta)
anddsdp(theta)
- Parameters:
slm (lm_max) [complex ,in]
gradtc (n_theta_max,n_phi_max) [real ,out]
gradpc (n_theta_max,n_phi_max) [real ,out]
lcut [integer ,in]
- Called from:
- Call to:
- subroutine shtransforms/native_spat_to_sph(scal, f1lm, lcut)¶
Legendre transform (n_r,n_theta,m) to (n_r,l,m)
- Parameters:
scal (*,*) [real ,inout]
f1lm (lm_max) [complex ,out]
lcut [integer ,in]
- Called from:
- Call to:
- subroutine shtransforms/native_spat_to_sph_tor(vt, vp, f1lm, f2lm, lcut)¶
Vector Legendre transform
vt(n_r,n_theta,m)
,vp(n_r,n_theta,m)
toSpheroidal(n_r,l,m)
andToroidal(n_r,l,m)
- Parameters:
- Called from:
- Call to:
sht_native.f90
¶
Quick access
- Routines:
axi_to_spat()
,finalize_sht()
,initialize_sht()
,pol_to_curlr_spat()
,pol_to_grad_spat()
,scal_to_grad_spat()
,scal_to_sh()
,scal_to_spat()
,spat_to_qst()
,spat_to_sphertor()
,sphtor_to_spat()
,toraxi_to_spat()
,torpol_to_curl_spat()
,torpol_to_curl_spat_ic()
,torpol_to_dphspat()
,torpol_to_spat()
,torpol_to_spat_ic()
Needed modules
iso_fortran_env
(output_unit()
)precision_mod
(cp()
): This module controls the precision used in MagICconstants
(ci()
,one()
,zero()
): module containing constants and parameters used in the code.truncation
(m_max()
,l_max()
,n_theta_max()
,n_phi_max()
,minc()
,lm_max()
): This module defines the grid points and the truncationhorizontal_data
(dlh()
,o_sin_theta_e2()
,o_sin_theta()
): Module containing functions depending on longitude and latitude plus help arrays depending on degree and orderparallel_mod
: This module contains the blocking informationshtransforms
: This module is used when the native built-in SH transforms of MagIC are used. Those are much slower than SHTns, and are not recommanded for production runs!
Variables
Subroutines and functions
- subroutine sht/initialize_sht(l_scrambled_theta)¶
- Parameters:
l_scrambled_theta [logical ,out]
- Called from:
- Call to:
- subroutine sht/finalize_sht()¶
- Called from:
- Call to:
- subroutine sht/scal_to_spat(slm, fieldc, lcut)¶
transform a spherical harmonic field into grid space
- Parameters:
slm (*) [complex ,in]
fieldc (*,*) [real ,out]
lcut [integer ,in]
- Called from:
- Call to:
- subroutine sht/scal_to_grad_spat(slm, gradtc, gradpc, lcut)¶
transform a scalar spherical harmonic field into it’s gradient on the grid
- Parameters:
slm (*) [complex ,in]
gradtc (*,*) [real ,out]
gradpc (*,*) [real ,out]
lcut [integer ,in]
- Call to:
- subroutine sht/pol_to_grad_spat(slm, gradtc, gradpc, lcut)¶
- Parameters:
slm (*) [complex ,in]
gradtc (*,*) [real ,out]
gradpc (*,*) [real ,out]
lcut [integer ,in]
- Call to:
- subroutine sht/torpol_to_spat(wlm, dwlm, zlm, vrc, vtc, vpc, lcut)¶
- Parameters:
wlm (*) [complex ,in]
dwlm (*) [complex ,in]
zlm (*) [complex ,in]
vrc (*,*) [real ,out]
vtc (*,*) [real ,out]
vpc (*,*) [real ,out]
lcut [integer ,in]
- Called from:
- Call to:
- subroutine sht/sphtor_to_spat(dwlm, zlm, vtc, vpc, lcut)¶
- Parameters:
dwlm (*) [complex ,in]
zlm (*) [complex ,in]
vtc (*,*) [real ,out]
vpc (*,*) [real ,out]
lcut [integer ,in]
- Call to:
- subroutine sht/torpol_to_curl_spat_ic(r, r_icb, dblm, ddblm, jlm, djlm, cbr, cbt, cbp)¶
This is a QST transform that contains the transform for the inner core to compute the three components of the curl of B.
- Parameters:
r [real ,in]
r_icb [real ,in]
dblm (*) [complex ,in]
ddblm (*) [complex ,in]
jlm (*) [complex ,in]
djlm (*) [complex ,in]
cbr (*,*) [real ,out]
cbt (*,*) [real ,out]
cbp (*,*) [real ,out]
- Call to:
- subroutine sht/torpol_to_spat_ic(r, r_icb, wlm, dwlm, zlm, br, bt, bp)¶
This is a QST transform that contains the transform for the inner core.
- Parameters:
r [real ,in]
r_icb [real ,in]
wlm (*) [complex ,in]
dwlm (*) [complex ,in]
zlm (*) [complex ,in]
br (*,*) [real ,out]
bt (*,*) [real ,out]
bp (*,*) [real ,out]
- Called from:
- Call to:
- subroutine sht/torpol_to_dphspat(dwlm, zlm, dvtdp, dvpdp, lcut)¶
Computes horizontal phi derivative of a toroidal/poloidal field
- Parameters:
dwlm (*) [complex ,in]
zlm (*) [complex ,in]
dvtdp (*,*) [real ,out]
dvpdp (*,*) [real ,out]
lcut [integer ,in]
- Call to:
- subroutine sht/pol_to_curlr_spat(qlm, cvrc, lcut)¶
- Parameters:
qlm (*) [complex ,in]
cvrc (*,*) [real ,out]
lcut [integer ,in]
- Call to:
- subroutine sht/torpol_to_curl_spat(or2, blm, ddblm, jlm, djlm, cvrc, cvtc, cvpc, lcut)¶
– Input variables
- Parameters:
or2 [real ,in]
blm (*) [complex ,in]
ddblm (*) [complex ,in]
jlm (*) [complex ,in]
djlm (*) [complex ,in]
cvrc (*,*) [real ,out]
cvtc (*,*) [real ,out]
cvpc (*,*) [real ,out]
lcut [integer ,in]
- Call to:
- subroutine sht/scal_to_sh(f, flm, lcut)¶
- Parameters:
f (*,*) [real ,inout]
flm (*) [complex ,out]
lcut [integer ,in]
- Called from:
- Call to:
- subroutine sht/spat_to_qst(f, g, h, qlm, slm, tlm, lcut)¶
- Parameters:
f (*,*) [real ,inout]
g (*,*) [real ,inout]
h (*,*) [real ,inout]
qlm (*) [complex ,out]
slm (*) [complex ,out]
tlm (*) [complex ,out]
lcut [integer ,in]
- Call to:
- subroutine sht/spat_to_sphertor(f, g, flm, glm, lcut)¶
- Parameters:
f (*,*) [real ,inout]
g (*,*) [real ,inout]
flm (*) [complex ,out]
glm (*) [complex ,out]
lcut [integer ,in]
- Called from:
- Call to:
- subroutine sht/axi_to_spat(fl_ax, f)¶
- Parameters:
fl_ax (1 + l_max) [complex ,in]
f (*) [real ,out]
- Call to:
Linear algebra¶
matrices.f90
¶
Quick access
- Types:
- Variables:
Needed modules
precision_mod
: This module controls the precision used in MagIC
Types
- type real_matrices/type_realmat¶
- Type fields:
% dat (*,*) [real ,allocatable]
% l_pivot [logical ]
% ncol [integer ]
% nrow [integer ]
% pivot (*) [integer ,allocatable]
Variables
Quick access
- Types:
Needed modules
precision_mod
: This module controls the precision used in MagICmem_alloc
: This little module is used to estimate the global memory allocation used in MagIC
Types
- type dense_matrices/type_densemat¶
Variables
Subroutines and functions
- subroutine dense_matrices/prepare(this, info)¶
- Parameters:
this [real ]
info [integer ,out]
- Call to:
- subroutine dense_matrices/solve_real_single(this, rhs)¶
- Parameters:
this [real ]
rhs (*) [real ,inout]
- subroutine dense_matrices/solve_complex_single(this, rhs)¶
- Parameters:
this [real ]
rhs (*) [complex ,inout]
- subroutine dense_matrices/solve_real_multi(this, rhs, nrhs)¶
- Parameters:
this [real ]
rhs (*,*) [real ,inout]
nrhs [integer ,in]
- subroutine dense_matrices/set_data(this, dat)¶
- Parameters:
this [real ]
dat (*,*) [real ,in]
Quick access
- Types:
Needed modules
precision_mod
: This module controls the precision used in MagICmem_alloc
: This little module is used to estimate the global memory allocation used in MagICalgebra
(solve_tridiag()
,prepare_tridiag()
,prepare_band()
,solve_band()
)
Types
- type band_matrices/type_bandmat¶
- Type fields:
% d (*) [real ,allocatable]
% dl (*) [real ,allocatable]
% du (*) [real ,allocatable]
% du2 (*) [real ,allocatable]
% kl [integer ]
% ku [integer ]
Variables
Subroutines and functions
- subroutine band_matrices/prepare(this, info)¶
- Parameters:
this [real ]
info [integer ,out]
- Call to:
- subroutine band_matrices/solve_real_single(this, rhs)¶
- Parameters:
this [real ]
rhs (*) [real ,inout]
- subroutine band_matrices/solve_complex_single(this, rhs)¶
- Parameters:
this [real ]
rhs (*) [complex ,inout]
- subroutine band_matrices/solve_real_multi(this, rhs, nrhs)¶
- Parameters:
this [real ]
rhs (*,*) [real ,inout]
nrhs [integer ,in]
- subroutine band_matrices/set_data(this, dat)¶
- Parameters:
this [real ]
dat (*,*) [real ,in]
algebra.f90
¶
Quick access
- Variables:
gemm
,gemv
,solve_band
,solve_band_complex_rhs
,solve_band_real_rhs
,solve_band_real_rhs_multi
,solve_bordered
,solve_bordered_complex_rhs
,solve_bordered_real_rhs
,solve_bordered_real_rhs_multi
,solve_mat
,solve_mat_complex_rhs
,solve_mat_real_rhs
,solve_mat_real_rhs_multi
,solve_tridiag
,solve_tridiag_complex_rhs
,solve_tridiag_real_rhs
,solve_tridiag_real_rhs_multi
,zero_tolerance
- Routines:
prepare_band()
,prepare_bordered()
,prepare_mat()
,prepare_tridiag()
Needed modules
omp_lib
precision_mod
: This module controls the precision used in MagICconstants
(one()
,zero()
): module containing constants and parameters used in the code.useful
(abortrun()
): This module contains several useful routines.
Variables
- algebra/gemm [private]¶
- algebra/gemv [private]¶
- algebra/solve_band [public]¶
- algebra/solve_band_complex_rhs [private]¶
- algebra/solve_band_real_rhs [private]¶
- algebra/solve_band_real_rhs_multi [private]¶
- algebra/solve_bordered [public]¶
- algebra/solve_bordered_complex_rhs [private]¶
- algebra/solve_bordered_real_rhs [private]¶
- algebra/solve_bordered_real_rhs_multi [private]¶
- algebra/solve_mat [public]¶
- algebra/solve_mat_complex_rhs [private]¶
- algebra/solve_mat_real_rhs [private]¶
- algebra/solve_mat_real_rhs_multi [private]¶
- algebra/solve_tridiag [public]¶
- algebra/solve_tridiag_complex_rhs [private]¶
- algebra/solve_tridiag_real_rhs [private]¶
- algebra/solve_tridiag_real_rhs_multi [private]¶
- algebra/zero_tolerance [real,private/parameter/optional/default=1.0e-15_cp]¶
Subroutines and functions
- subroutine algebra/prepare_mat(a, ia, n, ip, info)¶
LU decomposition the real matrix a(n,n) via gaussian elimination
- Parameters:
a (ia,*) [real ,inout] :: real nxn matrix on input, lu-decomposed matrix on output
ia [integer ,in,] :: first dimension of a (must be >= n)
n [integer ,in] :: 2nd dimension and rank of a
ip (*) [integer ,out] :: pivot pointer array
info [integer ,out] :: error message when /= 0
- Called from:
- Call to:
- subroutine algebra/prepare_band(abd, n, kl, ku, pivot, info)¶
- Parameters:
abd (1 + 2 * kl + ku,n) [real ,inout]
n [integer ,in,]
kl [integer ,in,required]
ku [integer ,in,required]
pivot (n) [integer ,out]
info [integer ,out]
- Called from:
- subroutine algebra/prepare_tridiag(dl, d, du, du2, n, pivot, info)¶
- Parameters:
dl (*) [real ,inout]
d (*) [real ,inout]
du (*) [real ,inout]
du2 (*) [real ,out]
n [integer ,in]
pivot (*) [integer ,out]
info [integer ,out]
- subroutine algebra/prepare_bordered(a1, a2, a3, a4, lena1, n_boundaries, kl, ku, pivota1, pivota4, info)¶
– Input variables
- Parameters:
a1 (1 + 2 * kl + ku,lena1) [real ,inout]
a2 (lena1,n_boundaries) [real ,inout]
a3 (lena1) [real ,inout]
a4 (n_boundaries,n_boundaries) [real ,inout]
lena1 [integer ,in,]
n_boundaries [integer ,in,]
kl [integer ,in,required]
ku [integer ,in,required]
pivota1 (lena1) [integer ,out]
pivota4 (n_boundaries) [integer ,out]
info [integer ,out]
- Called from:
- Call to: