Multidip  1.0
Multi-photon matrix elements
multidip_params.F90
Go to the documentation of this file.
1 ! Copyright 2020
2 !
3 ! For a comprehensive list of the developers that contributed to these codes
4 ! see the UK-AMOR website.
5 !
6 ! This file is part of UKRmol-out (UKRmol+ suite).
7 !
8 ! UKRmol-out is free software: you can redistribute it and/or modify
9 ! it under the terms of the GNU General Public License as published by
10 ! the Free Software Foundation, either version 3 of the License, or
11 ! (at your option) any later version.
12 !
13 ! UKRmol-out is distributed in the hope that it will be useful,
14 ! but WITHOUT ANY WARRANTY; without even the implied warranty of
15 ! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16 ! GNU General Public License for more details.
17 !
18 ! You should have received a copy of the GNU General Public License
19 ! along with UKRmol-out (in source/COPYING). Alternatively, you can also visit
20 ! <https://www.gnu.org/licenses/>.
21 !
30 
31  use precisn_gbl, only: wp
32 
33  implicit none
34 
35  integer, parameter :: ntermsasy = 5
36  integer, parameter :: ntermsppt = 3
37  integer, parameter :: nmaxphotons = bit_size(0)
38  integer, parameter :: max_romb_level = 20
39  integer, parameter :: max_levin_level = 20
40  integer, parameter :: cheb_order = 5
41 
42  integer :: maxtarg = 0
43  logical :: cache_integrals = .false.
44  logical :: check_dipoles = .true.
45  logical :: closed_interm = .false.
46  logical :: coulomb_check = .true.
47  logical :: print_warnings = .true.
48  logical :: custom_kmatrices = .false.
49  logical :: extend_istate = .false.
50  integer :: num_integ_algo = 2
51  real(wp) :: epsrel = 1e-6
52  real(wp) :: coultol = 1e-4
53  real(wp) :: closed_range = 5.0
54 
55  real(wp), parameter :: alpha = 1/137.03599907_wp
56  real(wp), parameter :: rzero = 0
57  real(wp), parameter :: rone = 1
58  real(wp), parameter :: rhalf = 0.5
59  real(wp), parameter :: pi = 4*atan(1.0_wp)
60 
61  complex(wp), parameter :: czero = 0
62  complex(wp), parameter :: cone = 1
63  complex(wp), parameter :: imu = (0, 1)
64 
65  character(len=1), parameter :: compt(3) = ['x', 'y', 'z']
66 
67  integer, parameter :: carti(3) = [3, 1, 2]
68  integer, parameter :: cartm(3) = [+1, -1, 0]
69 
70 contains
71 
72  subroutine read_input_namelist (input, order, lusct, lukmt, lubnd, nkset, polar, omega, verbose, rmt_data, first_IP, rasym, &
73  raw, erange, mpiio)
74 
75  integer, intent(in) :: input
76  integer, intent(inout) :: order, lusct(8), lukmt(8), lubnd, nkset(8), erange(2)
77  logical, intent(inout) :: verbose, mpiio
78  real(wp), intent(inout) :: omega(nMaxPhotons), first_IP, rasym
79  complex(wp), intent(inout) :: polar(3, nMaxPhotons)
80  character(len=256), intent(inout) :: rmt_data, raw
81 
82  namelist /mdip/ order, lusct, lukmt, lubnd, nkset, polar, omega, verbose, rmt_data, first_ip, rasym, raw, erange, mpiio, &
85 
86  order = 1
87  lusct = 0
88  lukmt = 0
89  lubnd = 0
90  nkset = 1
91  polar = 0
92  omega = 0
93  rasym = 0
94  first_ip = 0
95  erange = 0
96  verbose = .false.
97  rmt_data = 'molecular_data'
98  raw = ''
99  mpiio = .false.
100 
101  read (input, nml = mdip)
102 
103  if (order < 1) then
104  print '(A)', 'Order must be positive'
105  stop 1
106  end if
107 
108  if (order > nmaxphotons) then
109  print '(A,I0)', 'Order must be <= ', nmaxphotons
110  stop 1
111  end if
112 
113  if (raw /= '' .and. raw /= 'sph' .and. raw /= 'xyz' .and. raw /= 'both') then
114  print '(A)', 'The parameter "raw" must be empty or one of: "sph", "xyz", "both".'
115  stop 1
116  end if
117 
118 #ifndef WITH_SCALAPACK
119  if (mpiio) then
120  print '(A)', 'MPI-IO mode requires compiling MULTIDIP with ScaLAPACK support'
121  stop 1
122  end if
123 #endif
124 
125  polar(:, order + 1:nmaxphotons) = 0
126  omega(order + 1:nmaxphotons) = 0
127 
128  end subroutine read_input_namelist
129 
130 end module multidip_params
Hard-coded parameters of MULTIDIP.
integer, parameter ntermsppt
Number of terms in expansion of the exponential integral.
subroutine read_input_namelist(input, order, lusct, lukmt, lubnd, nkset, polar, omega, verbose, rmt_data, first_IP, rasym, raw, erange, mpiio)
logical extend_istate
Continue initial state from the known inner region expansion outwards.
real(wp), parameter pi
integer, parameter max_romb_level
Maximal nesting level for Romberg integration.
logical closed_interm
Consider weakly closed channel in intermediate states (unfinished!)
logical custom_kmatrices
Ignore RSOLVE K-matrices and calculate them from scratch.
integer num_integ_algo
Numerical integration mode (1 = Romberg, 2 = Levin)
real(wp), parameter rhalf
integer, dimension(3), parameter carti
Position of a Cartesian coordinate in the real spherical basis (y, z, x).
real(wp) epsrel
Romberg integration relative tolerance for convergence.
character(len=1), dimension(3), parameter compt
logical cache_integrals
Cache Coulomb-Green integrals in memory (but disables some threading)
complex(wp), parameter imu
real(wp), parameter rone
complex(wp), parameter cone
logical coulomb_check
Skip integrals that cannot be asymptotically approximated well.
real(wp) coultol
Coulomb matching relative tolerance (decides whether to use asym. forms)
logical print_warnings
Print warnings about non-converged integrals.
integer, parameter max_levin_level
Maximal nesting level for Levin integration.
real(wp) closed_range
Energy interval (a.u.) before threshold for inclusion of closed channels.
logical check_dipoles
Check that all dipole matrices in molecular_data are nonzero.
integer, parameter cheb_order
Order of Chebyshev interpolation used in Levin quadrature.
real(wp), parameter rzero
integer, parameter ntermsasy
Number of terms in expansion of the Coulomb-Hankel functions.
integer maxtarg
Maximal number of target states to calculate dipoles for (0 = all)
real(wp), parameter alpha
Fine structure constant.
integer, dimension(3), parameter cartm
Real spherical harmonic m-value corresponding to given a Cartesian coord.
complex(wp), parameter czero
integer, parameter nmaxphotons
The limit on number of photons is nested_cgreen_integ