GBTOlib: library for evaluation of molecular integrals in mixed Gaussian / B-spline basis  111
Functions/Subroutines | Variables
const_gbl Module Reference

Functions/Subroutines

subroutine redirect_stdout_to_file (myrank, master_stdout)
 This routine can be called only once. It redirects standard output into the file with name "stdout_file_name"//myrank where // stands for string concatenation and myrank is integer. The idea is that the mpi_mod module routine mpi_start redirects the standard output for each process into its own file by calling this routine with the rank of the process as argument. See also description of the parameters stdout_unit_base, stdout_unit_step. More...
 
recursive subroutine set_verbosity_level (level)
 Reset the output verbosity. More...
 

Variables

integer, protected stdout = output_unit
 unit for standard output. This value is set by the routine redirect_stdout_to_file but it can be done only once during the course of the program. Typically, the value of stdout will be chosen different for each rank (in case of MPI runs) so that each rank will use its own unit for std. output. This also ensures that no modifications in the code are needed concerning standard output. More...
 
integer, parameter stdin = input_unit
 unit for standard input More...
 
integer, parameter stderr = error_unit
 unit for standard error More...
 
integer devnull = output_unit
 unit for null stream (to be reset in redirect_stdout_to_file) More...
 
integer, parameter char_len = character_storage_size
 The size in bits of a character storage unit. More...
 
integer level1 = output_unit
 Main standard output: representative summary of the most important data and workflow progress. More...
 
integer level2 = output_unit
 Extended standard output: additional data that are mostly not of direct interest. More...
 
integer level3 = output_unit
 Debug standard output: debugging-related output (subroutine entry/exit markers, shell data) More...
 
character(len= *), parameter stdout_file_name = "log_file."
 Leading part of the file-name that will be used to collect standard output for each process. Each process gets its own file for standard output, e.g. "log_file.2" will contain standard output from MPI process with rank = 2. More...
 
integer, parameter stdout_unit_base = 10000
 The base number that will stand for the unit number for each processes's standard output (see mpi_mod). This number must be chosen large enough so it cannot collide with any other file units that may be open during the program run. More...
 
integer, parameter stdout_unit_step = 1
 Parameter used in redirect_stdout_to_file to compute the unit number for each processes's standard output: stdout = stdout_unit_base+myrank*stdout_unit_step This parameter allows to choose e.g. master's stdout to be unit 6 while the other processes will not interfere with any other fort units, e.g.: stdout_unit_base=6, stdout_unit_step=1000. More...
 
integer, parameter line_len = 207
 length of one line More...
 
integer, parameter len_ufmat = 11
 
character(len_ufmat-2), parameter fmat = 'formatted'
 
character(len_ufmat), parameter ufmat = 'unformatted'
 
character(len=19), parameter no_header = 'No header specified'
 Default character string used for the int_input_outputname variable. More...
 
integer, parameter number_of_types_el_ints = 5
 How many different types of 1-electron integrals have been implemented so far: Overlap (1), Kinetic energy (2), Kinetic energy+Nuclear attraction (3), Nuclear attraction (4), property integrals (5). More...
 
character(len= *), parameter one_electron_ints = 'One electron integrals'
 Character parameters used by the method l_mol_basismolecular_integrals to identify the type of molecular integral requested. These headers are used by the user to request the specific integral. More...
 
character(len= *), parameter overlap_ints = 'Overlap integrals'
 
character(len= *), parameter kinetic_ints = 'Kinetic energy integrals'
 
character(len= *), parameter property_ints = 'Property integrals'
 
character(len= *), parameter nuc_rep_att_ints = 'Nuclear attraction integrals'
 
character(len= *), parameter one_elham = 'One electron Hamiltonian integrals'
 
character(len= *), parameter two_el_ints = 'Two-electron integrals'
 
character(len= *), parameter two_el_ints_prefetched = 'Prefetched two-electron integrals'
 
character(len= *), parameter one_p_sym_ints = 'One-particle AO->MO transformed integrals for symmetric operators'
 
character(len= *), parameter two_p_sym_ints = 'Two-particle AO->MO transformed integrals for symmetric operators'
 
character(len= *), parameter ijkl_indices_header = 'ijkl indices'
 
character(line_len), parameter data_file_obj_id_10 = 'DATA FILE version 1.0'
 Character string present on the first line of a file which identifies the file as readable by the data_file_obj object. The length of this string must be line_len since that length is used in the data_file module when reading the header from the file. More...
 
character(line_len), parameter data_file_obj_id_11 = 'DATA FILE version 1.1'
 
character(line_len), parameter data_file_obj_id = data_file_obj_id_11
 
integer, parameter orbs_line = 6
 Number of molecular orbitals for which orbital coefficients will be printed on one line by the routine molecular_orbital_basis_objprint_orbitals. More...
 
integer, parameter max_contr_len = 20
 Maximum allowed number of contraction coefficients defining the contracted GTO function. This parameter can be increased/decreased as needed. It is used throught the program to define dimensions of some arrays. More...
 
integer, parameter fit_order = 2
 The coefficients below are coefficients in the fit: mmax(T) = -18.206_cfp + 0.42734_cfp*T + 0.00075795_cfp*T**2 This fit was obtained using Mathematica and estimates, for each T, the maximum value of m for which Fm(T) can be calculated using the ASYMPTOTIC formula to full double precision accuracy rel. error .le. 10^(-15)). Fm(T) is the Boys function. The values are stored in single precision since that is enough to calculate mmax(T) with high enough accuracy sufficient for the estimate. More...
 
real(kind=cfp), dimension(1:fit_order+1), parameter wp_fit_terms = (/-18.206_cfp,0.42734_cfp,0.00075795_cfp/)
 
real(kind=cfp), dimension(1:fit_order+1), parameter ep_fit_terms = (/-21.408_cfp,0.19460_cfp,0.00093852_cfp/)
 The same as above (ASYMPTOTIC formula limit) but this time for full quadruple precision of Fm(T). More...
 
real(kind=cfp), dimension(1:fit_order+1), parameter wp_fit_terms_up = (/-18.561_cfp,0.3748_cfp,0.00088246_cfp/)
 The same as above but this time for the UPWARD recursion to full double precision accuracy. More...
 
real(kind=cfp), dimension(1:fit_order+1), parameter ep_fit_terms_up = (/-19.331_cfp,0.14211_cfp,0.00106170_cfp/)
 The same as above (UPWARD recursion formula limit) but this time for full quadruple precision of Fm(T). More...
 
real(kind=cfp), parameter boys_tol = epsilon(1.0_cfp)
 Convergence parameter for the calculation of the Boys function using boys_function. The convergence is chosen as the machine epsilon for real(kind=cfp) reals. More...
 
real(kind=ep1), parameter boys_tol_qprec = epsilon(1.0_ep1)
 Convergence parameter for the calculation of the Boys function using boys_function. We need to define the quad precision separately due to the routine boys_function_quad which is explicitly written in quad precision. More...
 
real(kind=cfp), parameter boys_f_dprec_asym_thr = 60.0_cfp
 Default value of Tmax which is used in init_boys to precalculate Fm(T) on the grid of T and m values. T=60 corresponds to the value for which, in double precision, the Boys function is calculated using the series expansion to full double precision for mmax=24. The corresponding value of T for quad precision is T=140. This parameter can be varied at will - it has only a minor effect on efficiency of the Taylor-based evaluation of Fm(T). More...
 
integer, parameter mmax = 24
 Default value of mmax for which is used in cgto_hgp to precalculate Fm(T) on the grid of T and m values. The minimum value of mmax is 4*L, where L is the maximum GTO L used in the calculations. The rest of the Boys function parameters have been set-up to give fully accurate results for m up to mmax=100. The value 24 below is just a sensible value that will be enough for most calculations. If the integral algorithm fails complaining that mmax is too small then set it here to the value at least 4*L where is the maximum L in your GTO atomic basis. More...
 
real(kind=cfp), parameter boys_f_grid_step_wp = 0.1_cfp
 Step in T in the grid of values of Fm(T) calculated for Taylor-expansion based evaluation of the Boys function for double precision calculations. More...
 
real(kind=cfp), parameter boys_f_grid_step_ep = 0.01_cfp
 Step in T in the grid of values of Fm(T) calculated for Taylor-expansion based evaluation of the Boys function for quad precision calculations. More...
 
integer, parameter imax_wp = 140
 Empirically found value for the largest power of T required to calculate the Boys function \(F_{m}(T)\) for \(T=0,\dots,60\) and \(m=0,\dots,24\) with the accuracy epsilon(1.0_wp). More...
 
integer, parameter imax_ep = 340
 The same as imax_wp but for quad precision result. More...
 
integer, parameter taylor_k_wp = 8
 Order of the Taylor expansion for Taylor-expansion-based evaluation of the Boys function Fm(T) for double precision calculations. More...
 
integer, parameter taylor_k_ep = 10
 Order of the Taylor expansion for Taylor-expansion-based evaluation of the Boys function Fm(T) for quad precision calculations. More...
 
real(kind=cfp), parameter int_rel_prec = 10**(-(precision(cfp_dummy)-5.0_cfp))
 Molecular integrals with relative precision last than or equal to this value will trigger an error. This parameter is used only as a default in integral_optionsprec, so this parameter can be effectively adjusted on run-time. More...
 
real(kind=cfp), parameter int_del_thr = 10**(-(precision(cfp_dummy)-4.0_cfp))
 Molecular integrals (contracted) smaller than this value will be neglected. Similarly to int_rel_prec this value is used only as a sensible default in integral_optionstol. More...
 
real(kind=cfp), parameter thrs_symm_ortho = 10**(-(precision(cfp_dummy)-9.0_cfp))
 Maximum allowed value for a cross overlap of two orbitals after the symmetric orthogonalization has been performed. More...
 
real(kind=cfp), parameter thrs_gs_ortho = 10**(-(precision(cfp_dummy)-5.0_cfp))
 Maximum allowed value for a cross overlap of two orbitals after the Gramm-Schmidt orthogonalization has been performed. More...
 
real(kind=cfp), parameter thrs_cf_sym_ortho_trans = 10**(-(precision(cfp_dummy)-5.0_cfp))
 Coefficients in the transformation matrix for the symmetric orthogonalization smaller than this value will be neglected. More...
 
real(kind=cfp), parameter thrs_lin_dep_gs_ortho = 10**(-(precision(cfp_dummy)-8.0_cfp))
 Threshold value for self-overlap of an orbital (before normalization) for the Gramm-Schmidt orthogonalization. If an orbital has a self-overlap (before normalization) smaller than this value then we assume that linear dependency in the orbital basis is present. In other words we decide that self-overlaps smaller than this value would lead to numerical problems. More...
 
real(kind=cfp), parameter thrs_orb_cf = 10**(-(precision(cfp_dummy)-3.0_cfp))
 Threshold for magnitude of final orbital coefficients as used by molecular_orbital_basis_objdelete_small_coefficients. More...
 
integer, parameter i_min_sparse = 10000
 The minimum number of elements to allocate for the gather part in the routine finalize_two_electron_integrals_sparse. More...
 
real(kind=cfp), parameter epsabs = 10**(-(precision(cfp_dummy)-5.0_cfp))
 Absolute precision for the numerical quadrature routine dqags. More...
 
real(kind=cfp), parameter epsrel = 10**(-(precision(cfp_dummy)-5.0_cfp))
 Relative precision for the numerical quadrature routine dqags. More...
 
integer, parameter limit = 1000
 Determines the maximum number of subintervals in the partition of the given integration interval in dqags. More...
 
integer, parameter lenw = 4*LIMIT
 Dimensioning parameter for dqags. More...
 
integer, parameter max_epstab = 152
 Dimensioning parameter in DQELG which controls the maximum length of the vector that can be extrapolated. The original SLATEC routine used 52. More...
 
integer, dimension(8, 8), parameter abel_prod_tab =RESHAPE( (/ 1,2,3,4,5,6,7,8, 2,1,4,3,6,5,8,7, 3,4,1,2,7,8,5,6, 4,3,2,1,8,7,6,5, 5,6,7,8,1,2,3,4, 6,5,8,7,2,1,4,3, 7,8,5,6,3,4,1,2, 8,7,6,5,4,3,2,1 /), (/8,8/) )
 Multiplication table for Abelian point groups. More...
 
integer, parameter c1_id = 1
 Numerical identifier of the C1 point group-symmetry. More...
 
integer, parameter cs_id = 2
 Numerical identifier of the Cs point group-symmetry. More...
 
integer, parameter c2_id = 3
 Numerical identifier of the C2 point group-symmetry. More...
 
integer, parameter ci_id = 4
 Numerical identifier of the Ci point group-symmetry. More...
 
integer, parameter c2v_id = 5
 Numerical identifier of the C2v point group-symmetry. More...
 
integer, parameter c2h_id = 6
 Numerical identifier of the C2h point group-symmetry. More...
 
integer, parameter d2_id = 7
 Numerical identifier of the D2 point group-symmetry. More...
 
integer, parameter d2h_id = 8
 Numerical identifier of the D2h point group-symmetry. More...
 
integer, parameter nuc_nam_len = 2
 Length of the character variable 'name' in the nucleus_type object. More...
 
character(len=nuc_nam_lennam_scattering_centre = 'sc'
 Name of the scattering centre. More...
 
integer, parameter id_scattering_centre = 0
 ID assigned to the scattering centre. Do not change this value. More...
 
integer, parameter sym_op_nam_len = 3
 Length of the character variable sym_op specifying the symmetry operation in the geometry_obj object. More...
 
character(len=sym_op_nam_len), dimension(4), parameter c2v_names = (/'A1','B1','B2','A2'/)
 Names and order of IRR for C2v point group. Molpro order. More...
 
integer, dimension(4, 2), parameter c2v_char_tab =RESHAPE( (/ 1, 1,-1,-1, 1,-1, 1,-1 /), (/4,2/) )
 s_v, s_vp symmetry elements. Atkins, Tables for Group theory. More...
 
character(len=sym_op_nam_len), dimension(8), parameter d2h_names = (/'Ag ','B3u','B2u','B1g','B1u','B2g','B3g','Au '/)
 Names and order of IRR for D2h point group. Molpro order. More...
 
integer, dimension(8, 3), parameter d2h_char_tab =RESHAPE( (/ 1, 1, 1, 1,-1,-1,-1,-1, 1, 1,-1,-1, 1, 1,-1,-1, 1,-1, 1,-1, 1,-1, 1,-1 /), (/8,3/) )
 s_3, s_2, s_1 symmetry elements. Atkins, Tables for Group theory. More...
 
character(len=sym_op_nam_len), dimension(4), parameter c2h_names = (/'Ag','Au','Bu','Bg'/)
 Names and order of IRR for C2h point group. Molpro order. More...
 
integer, dimension(4, 3), parameter c2h_char_tab =RESHAPE( (/ 1, 1,-1,-1, 1,-1,-1, 1, 1,-1, 1,-1 /), (/4,3/) )
 c2, c_i, s_h symmetry elements. Atkins, Tables for Group theory. More...
 
character(len=sym_op_nam_len), dimension(4), parameter d2_names = (/'A ','B3','B2','B1'/)
 Names and order of IRR for D2 point group. Molpro order. More...
 
integer, dimension(4, 2), parameter d2_char_tab =RESHAPE( (/ 1,-1, 1,-1, 1, 1,-1,-1 /), (/4,2/) )
 c2_2, c2_1 symmetry elements. Atkins, Tables for Group theory. More...
 
character(len=sym_op_nam_len), dimension(2), parameter c2_names = (/'A','B'/)
 Names and order of IRR for C2 point group. Molpro order. More...
 
integer, dimension(2, 1), parameter c2_char_tab =RESHAPE( (/ 1,-1 /), (/2,1/) )
 c2 symmetry element. Atkins, Tables for Group theory. More...
 
character(len=sym_op_nam_len), dimension(2), parameter cs_names = (/'Ap ','App'/)
 Names and order of IRR for Cs point group. Molpro order. More...
 
integer, dimension(2, 1), parameter cs_char_tab =RESHAPE( (/ 1,-1 /), (/2,1/) )
 s_h symmetry element. Atkins, Tables for Group theory. More...
 
character(len=sym_op_nam_len), dimension(2), parameter ci_names = (/'Ag','Au'/)
 Names and order of IRR for Ci point group. Molpro order. More...
 
integer, dimension(2, 1), parameter ci_char_tab =RESHAPE( (/ 1,-1 /), (/2,1/) )
 c_i symmetry element. Atkins, Tables for Group theory. More...
 
character(len=sym_op_nam_len), dimension(1), parameter c1_names = (/'A'/)
 IRR for Ci point group. Molpro order. More...
 
integer, dimension(1, 1), parameter c1_char_tab =RESHAPE( (/ 1 /), (/1,1/) )
 Idetity symmetry element. Atkins, Tables for Group theory. More...
 
character(len=sym_op_nam_len), dimension(8, 8), parameter pg_irr_names = RESHAPE ((/'A ',' ',' ',' ',' ',' ',' ',' ', "A' ",'A" ',' ',' ',' ',' ',' ',' ', 'A ','B ',' ',' ',' ',' ',' ',' ', 'Ag ','Au ',' ',' ',' ',' ',' ',' ', 'A1 ','B1 ','B2 ','A2 ',' ',' ',' ',' ', 'Ag ','Au ','Bu ','Bg ',' ',' ',' ',' ', 'A ','B3 ','B2 ','B1 ',' ',' ',' ',' ', 'Ag ','B3u','B2u','B1g','B1u','B2g','B3g','Au '/), (/8,8/))
 Combined irreducible representation labels used in Psi4 Molden files. More...
 
real(kind=cfp), parameter y_lm_size_threshold = 100.0_cfp
 If it is requested that the Y_lm functions (evaluated by eval_BTO_CGTO_Y_lm) are saved to disk then all Y_lm with size (in MiB) greater than this value will be saved to disk. If you want to save all Y_lm regardless of their size then set this parameter to 0.0_cfp. More...
 
integer, parameter mib = 1048576
 Number of bytes in one Mib. More...
 
integer, parameter cache_line_size = 64
 Cache line size (bytes). More...
 
integer, parameter tile = 64
 Tile size used for cache blocking. This is derived from cacheline size. More...
 
real(kind=cfp), dimension(103), parameter amass = (/ 1.0078246_cfp, 4.002601_cfp, 7.01600_cfp, 9.01218_cfp, 11.009307_cfp, 12.000000_cfp, 14.0030738_cfp, 15.9949141_cfp, 18.9984022_cfp, 19.992441_cfp, 22.9898_cfp, 23.98504_cfp, 26.98153_cfp, 27.976929_cfp, 30.973764_cfp, 31.9720727_cfp, 34.9688531_cfp, 39.962386_cfp, 38.96371_cfp, 39.96259_cfp, 44.95592_cfp, 48._cfp, 50.9440_cfp, 51.9405_cfp, 54.9380_cfp, 55.9349_cfp, 58.9332_cfp, 57.9353_cfp, 62.9296_cfp, 63.9291_cfp, 68.9257_cfp, 73.9219_cfp, 74.9216_cfp, 79.9165_cfp, 78.91839_cfp, 83.91151_cfp, 84.9117_cfp, 87.9056_cfp, 88.9059_cfp, 89.9043_cfp, 92.9060_cfp, 97.9055_cfp, 98._cfp, 101.9037_cfp, 102.9048_cfp, 107.90389_cfp, 106.90509_cfp, 113.9036_cfp, 114.9041_cfp, 120._cfp, 120.9038_cfp, 129.9067_cfp, 126.90466_cfp, 131.90416_cfp, 132.9051_cfp, 137.9050_cfp, 138.9061_cfp, 139.9053_cfp, 140.9074_cfp, 141.9075_cfp, 145._cfp, 151.9195_cfp, 152.9209_cfp, 157.9241_cfp, 159.9250_cfp, 163.9288_cfp, 164.9303_cfp, 165.9304_cfp, 168.9344_cfp, 173.9390_cfp, 174.9409_cfp, 179.9468_cfp, 180.9480_cfp, 183.9510_cfp, 186.9560_cfp, 192._cfp, 192.9633_cfp, 194.9648_cfp, 196.9666_cfp, 201.970625_cfp, 204.9745_cfp, 207.9766_cfp, 208.9804_cfp, 209._cfp, 210._cfp, 222._cfp, 223._cfp, 226._cfp, 227._cfp, 232._cfp, 231._cfp, 238._cfp, 237._cfp, 244._cfp, 243._cfp, 247._cfp, 247._cfp, 251._cfp, 252._cfp, 257._cfp, 258._cfp, 259._cfp, 260._cfp /)
 Masses of the elements in the periodic table. Taken from DENPROP. More...
 

Function/Subroutine Documentation

◆ redirect_stdout_to_file()

subroutine const_gbl::redirect_stdout_to_file ( integer, intent(in)  myrank,
logical, intent(in), optional  master_stdout 
)

This routine can be called only once. It redirects standard output into the file with name "stdout_file_name"//myrank where // stands for string concatenation and myrank is integer. The idea is that the mpi_mod module routine mpi_start redirects the standard output for each process into its own file by calling this routine with the rank of the process as argument. See also description of the parameters stdout_unit_base, stdout_unit_step.

Here is the call graph for this function:

◆ set_verbosity_level()

recursive subroutine const_gbl::set_verbosity_level ( integer, intent(in), optional  level)

Reset the output verbosity.

Author
J Benda
Date
2020

Adjust the output verbosity of the library. The known levels are 1 to 4. Setting the level to 0 will turn off all optional output. Setting the level to 4 and higher will make the library print maximum information. Specifying the level is optional; without any argument given, this subroutine will reset the verbosity to its default level, which is 1.

Variable Documentation

◆ abel_prod_tab

integer, dimension(8,8), parameter const_gbl::abel_prod_tab =RESHAPE( (/ 1,2,3,4,5,6,7,8, 2,1,4,3,6,5,8,7, 3,4,1,2,7,8,5,6, 4,3,2,1,8,7,6,5, 5,6,7,8,1,2,3,4, 6,5,8,7,2,1,4,3, 7,8,5,6,3,4,1,2, 8,7,6,5,4,3,2,1 /), (/8,8/) )

Multiplication table for Abelian point groups.

◆ amass

real(kind=cfp), dimension(103), parameter const_gbl::amass = (/ 1.0078246_cfp, 4.002601_cfp, 7.01600_cfp, 9.01218_cfp, 11.009307_cfp, 12.000000_cfp, 14.0030738_cfp, 15.9949141_cfp, 18.9984022_cfp, 19.992441_cfp, 22.9898_cfp, 23.98504_cfp, 26.98153_cfp, 27.976929_cfp, 30.973764_cfp, 31.9720727_cfp, 34.9688531_cfp, 39.962386_cfp, 38.96371_cfp, 39.96259_cfp, 44.95592_cfp, 48._cfp, 50.9440_cfp, 51.9405_cfp, 54.9380_cfp, 55.9349_cfp, 58.9332_cfp, 57.9353_cfp, 62.9296_cfp, 63.9291_cfp, 68.9257_cfp, 73.9219_cfp, 74.9216_cfp, 79.9165_cfp, 78.91839_cfp, 83.91151_cfp, 84.9117_cfp, 87.9056_cfp, 88.9059_cfp, 89.9043_cfp, 92.9060_cfp, 97.9055_cfp, 98._cfp, 101.9037_cfp, 102.9048_cfp, 107.90389_cfp, 106.90509_cfp, 113.9036_cfp, 114.9041_cfp, 120._cfp, 120.9038_cfp, 129.9067_cfp, 126.90466_cfp, 131.90416_cfp, 132.9051_cfp, 137.9050_cfp, 138.9061_cfp, 139.9053_cfp, 140.9074_cfp, 141.9075_cfp, 145._cfp, 151.9195_cfp, 152.9209_cfp, 157.9241_cfp, 159.9250_cfp, 163.9288_cfp, 164.9303_cfp, 165.9304_cfp, 168.9344_cfp, 173.9390_cfp, 174.9409_cfp, 179.9468_cfp, 180.9480_cfp, 183.9510_cfp, 186.9560_cfp, 192._cfp, 192.9633_cfp, 194.9648_cfp, 196.9666_cfp, 201.970625_cfp, 204.9745_cfp, 207.9766_cfp, 208.9804_cfp, 209._cfp, 210._cfp, 222._cfp, 223._cfp, 226._cfp, 227._cfp, 232._cfp, 231._cfp, 238._cfp, 237._cfp, 244._cfp, 243._cfp, 247._cfp, 247._cfp, 251._cfp, 252._cfp, 257._cfp, 258._cfp, 259._cfp, 260._cfp /)

Masses of the elements in the periodic table. Taken from DENPROP.

Warning
The masses of some elements (e.g. Boron) are way off from the NIST values!

◆ boys_f_dprec_asym_thr

real(kind=cfp), parameter const_gbl::boys_f_dprec_asym_thr = 60.0_cfp

Default value of Tmax which is used in init_boys to precalculate Fm(T) on the grid of T and m values. T=60 corresponds to the value for which, in double precision, the Boys function is calculated using the series expansion to full double precision for mmax=24. The corresponding value of T for quad precision is T=140. This parameter can be varied at will - it has only a minor effect on efficiency of the Taylor-based evaluation of Fm(T).

◆ boys_f_grid_step_ep

real(kind=cfp), parameter const_gbl::boys_f_grid_step_ep = 0.01_cfp

Step in T in the grid of values of Fm(T) calculated for Taylor-expansion based evaluation of the Boys function for quad precision calculations.

◆ boys_f_grid_step_wp

real(kind=cfp), parameter const_gbl::boys_f_grid_step_wp = 0.1_cfp

Step in T in the grid of values of Fm(T) calculated for Taylor-expansion based evaluation of the Boys function for double precision calculations.

◆ boys_tol

real(kind=cfp), parameter const_gbl::boys_tol = epsilon(1.0_cfp)

Convergence parameter for the calculation of the Boys function using boys_function. The convergence is chosen as the machine epsilon for real(kind=cfp) reals.

◆ boys_tol_qprec

real(kind=ep1), parameter const_gbl::boys_tol_qprec = epsilon(1.0_ep1)

Convergence parameter for the calculation of the Boys function using boys_function. We need to define the quad precision separately due to the routine boys_function_quad which is explicitly written in quad precision.

◆ c1_char_tab

integer, dimension(1,1), parameter const_gbl::c1_char_tab =RESHAPE( (/ 1 /), (/1,1/) )

Idetity symmetry element. Atkins, Tables for Group theory.

◆ c1_id

integer, parameter const_gbl::c1_id = 1

Numerical identifier of the C1 point group-symmetry.

◆ c1_names

character(len=sym_op_nam_len), dimension(1), parameter const_gbl::c1_names = (/'A'/)

IRR for Ci point group. Molpro order.

◆ c2_char_tab

integer, dimension(2,1), parameter const_gbl::c2_char_tab =RESHAPE( (/ 1,-1 /), (/2,1/) )

c2 symmetry element. Atkins, Tables for Group theory.

◆ c2_id

integer, parameter const_gbl::c2_id = 3

Numerical identifier of the C2 point group-symmetry.

◆ c2_names

character(len=sym_op_nam_len), dimension(2), parameter const_gbl::c2_names = (/'A','B'/)

Names and order of IRR for C2 point group. Molpro order.

◆ c2h_char_tab

integer, dimension(4,3), parameter const_gbl::c2h_char_tab =RESHAPE( (/ 1, 1,-1,-1, 1,-1,-1, 1, 1,-1, 1,-1 /), (/4,3/) )

c2, c_i, s_h symmetry elements. Atkins, Tables for Group theory.

◆ c2h_id

integer, parameter const_gbl::c2h_id = 6

Numerical identifier of the C2h point group-symmetry.

◆ c2h_names

character(len=sym_op_nam_len), dimension(4), parameter const_gbl::c2h_names = (/'Ag','Au','Bu','Bg'/)

Names and order of IRR for C2h point group. Molpro order.

◆ c2v_char_tab

integer, dimension(4,2), parameter const_gbl::c2v_char_tab =RESHAPE( (/ 1, 1,-1,-1, 1,-1, 1,-1 /), (/4,2/) )

s_v, s_vp symmetry elements. Atkins, Tables for Group theory.

◆ c2v_id

integer, parameter const_gbl::c2v_id = 5

Numerical identifier of the C2v point group-symmetry.

◆ c2v_names

character(len=sym_op_nam_len), dimension(4), parameter const_gbl::c2v_names = (/'A1','B1','B2','A2'/)

Names and order of IRR for C2v point group. Molpro order.

◆ cache_line_size

integer, parameter const_gbl::cache_line_size = 64

Cache line size (bytes).

◆ char_len

integer, parameter const_gbl::char_len = character_storage_size

The size in bits of a character storage unit.

◆ ci_char_tab

integer, dimension(2,1), parameter const_gbl::ci_char_tab =RESHAPE( (/ 1,-1 /), (/2,1/) )

c_i symmetry element. Atkins, Tables for Group theory.

◆ ci_id

integer, parameter const_gbl::ci_id = 4

Numerical identifier of the Ci point group-symmetry.

◆ ci_names

character(len=sym_op_nam_len), dimension(2), parameter const_gbl::ci_names = (/'Ag','Au'/)

Names and order of IRR for Ci point group. Molpro order.

◆ cs_char_tab

integer, dimension(2,1), parameter const_gbl::cs_char_tab =RESHAPE( (/ 1,-1 /), (/2,1/) )

s_h symmetry element. Atkins, Tables for Group theory.

◆ cs_id

integer, parameter const_gbl::cs_id = 2

Numerical identifier of the Cs point group-symmetry.

◆ cs_names

character(len=sym_op_nam_len), dimension(2), parameter const_gbl::cs_names = (/'Ap ','App'/)

Names and order of IRR for Cs point group. Molpro order.

◆ d2_char_tab

integer, dimension(4,2), parameter const_gbl::d2_char_tab =RESHAPE( (/ 1,-1, 1,-1, 1, 1,-1,-1 /), (/4,2/) )

c2_2, c2_1 symmetry elements. Atkins, Tables for Group theory.

◆ d2_id

integer, parameter const_gbl::d2_id = 7

Numerical identifier of the D2 point group-symmetry.

◆ d2_names

character(len=sym_op_nam_len), dimension(4), parameter const_gbl::d2_names = (/'A ','B3','B2','B1'/)

Names and order of IRR for D2 point group. Molpro order.

◆ d2h_char_tab

integer, dimension(8,3), parameter const_gbl::d2h_char_tab =RESHAPE( (/ 1, 1, 1, 1,-1,-1,-1,-1, 1, 1,-1,-1, 1, 1,-1,-1, 1,-1, 1,-1, 1,-1, 1,-1 /), (/8,3/) )

s_3, s_2, s_1 symmetry elements. Atkins, Tables for Group theory.

◆ d2h_id

integer, parameter const_gbl::d2h_id = 8

Numerical identifier of the D2h point group-symmetry.

◆ d2h_names

character(len=sym_op_nam_len), dimension(8), parameter const_gbl::d2h_names = (/'Ag ','B3u','B2u','B1g','B1u','B2g','B3g','Au '/)

Names and order of IRR for D2h point group. Molpro order.

◆ data_file_obj_id

character(line_len), parameter const_gbl::data_file_obj_id = data_file_obj_id_11

◆ data_file_obj_id_10

character(line_len), parameter const_gbl::data_file_obj_id_10 = 'DATA FILE version 1.0'

Character string present on the first line of a file which identifies the file as readable by the data_file_obj object. The length of this string must be line_len since that length is used in the data_file module when reading the header from the file.

◆ data_file_obj_id_11

character(line_len), parameter const_gbl::data_file_obj_id_11 = 'DATA FILE version 1.1'

◆ devnull

integer const_gbl::devnull = output_unit

unit for null stream (to be reset in redirect_stdout_to_file)

◆ ep_fit_terms

real(kind=cfp), dimension(1:fit_order+1), parameter const_gbl::ep_fit_terms = (/-21.408_cfp,0.19460_cfp,0.00093852_cfp/)

The same as above (ASYMPTOTIC formula limit) but this time for full quadruple precision of Fm(T).

◆ ep_fit_terms_up

real(kind=cfp), dimension(1:fit_order+1), parameter const_gbl::ep_fit_terms_up = (/-19.331_cfp,0.14211_cfp,0.00106170_cfp/)

The same as above (UPWARD recursion formula limit) but this time for full quadruple precision of Fm(T).

◆ epsabs

real(kind=cfp), parameter const_gbl::epsabs = 10**(-(precision(cfp_dummy)-5.0_cfp))

Absolute precision for the numerical quadrature routine dqags.

◆ epsrel

real(kind=cfp), parameter const_gbl::epsrel = 10**(-(precision(cfp_dummy)-5.0_cfp))

Relative precision for the numerical quadrature routine dqags.

◆ fit_order

integer, parameter const_gbl::fit_order = 2

The coefficients below are coefficients in the fit: mmax(T) = -18.206_cfp + 0.42734_cfp*T + 0.00075795_cfp*T**2 This fit was obtained using Mathematica and estimates, for each T, the maximum value of m for which Fm(T) can be calculated using the ASYMPTOTIC formula to full double precision accuracy rel. error .le. 10^(-15)). Fm(T) is the Boys function. The values are stored in single precision since that is enough to calculate mmax(T) with high enough accuracy sufficient for the estimate.

◆ fmat

character(len_ufmat-2), parameter const_gbl::fmat = 'formatted'

◆ i_min_sparse

integer, parameter const_gbl::i_min_sparse = 10000

The minimum number of elements to allocate for the gather part in the routine finalize_two_electron_integrals_sparse.

◆ id_scattering_centre

integer, parameter const_gbl::id_scattering_centre = 0

ID assigned to the scattering centre. Do not change this value.

◆ ijkl_indices_header

character(len=*), parameter const_gbl::ijkl_indices_header = 'ijkl indices'

◆ imax_ep

integer, parameter const_gbl::imax_ep = 340

The same as imax_wp but for quad precision result.

◆ imax_wp

integer, parameter const_gbl::imax_wp = 140

Empirically found value for the largest power of T required to calculate the Boys function \(F_{m}(T)\) for \(T=0,\dots,60\) and \(m=0,\dots,24\) with the accuracy epsilon(1.0_wp).

◆ int_del_thr

real(kind=cfp), parameter const_gbl::int_del_thr = 10**(-(precision(cfp_dummy)-4.0_cfp))

Molecular integrals (contracted) smaller than this value will be neglected. Similarly to int_rel_prec this value is used only as a sensible default in integral_optionstol.

◆ int_rel_prec

real(kind=cfp), parameter const_gbl::int_rel_prec = 10**(-(precision(cfp_dummy)-5.0_cfp))

Molecular integrals with relative precision last than or equal to this value will trigger an error. This parameter is used only as a default in integral_optionsprec, so this parameter can be effectively adjusted on run-time.

Warning
Not all integral calculation routines are necessarily using this parameter.

◆ kinetic_ints

character(len=*), parameter const_gbl::kinetic_ints = 'Kinetic energy integrals'

◆ len_ufmat

integer, parameter const_gbl::len_ufmat = 11

◆ lenw

integer, parameter const_gbl::lenw = 4*LIMIT

Dimensioning parameter for dqags.

◆ level1

integer const_gbl::level1 = output_unit

Main standard output: representative summary of the most important data and workflow progress.

◆ level2

integer const_gbl::level2 = output_unit

Extended standard output: additional data that are mostly not of direct interest.

◆ level3

integer const_gbl::level3 = output_unit

Debug standard output: debugging-related output (subroutine entry/exit markers, shell data)

◆ limit

integer, parameter const_gbl::limit = 1000

Determines the maximum number of subintervals in the partition of the given integration interval in dqags.

◆ line_len

integer, parameter const_gbl::line_len = 207

length of one line

◆ max_contr_len

integer, parameter const_gbl::max_contr_len = 20

Maximum allowed number of contraction coefficients defining the contracted GTO function. This parameter can be increased/decreased as needed. It is used throught the program to define dimensions of some arrays.

◆ max_epstab

integer, parameter const_gbl::max_epstab = 152

Dimensioning parameter in DQELG which controls the maximum length of the vector that can be extrapolated. The original SLATEC routine used 52.

◆ mib

integer, parameter const_gbl::mib = 1048576

Number of bytes in one Mib.

◆ mmax

integer, parameter const_gbl::mmax = 24

Default value of mmax for which is used in cgto_hgp to precalculate Fm(T) on the grid of T and m values. The minimum value of mmax is 4*L, where L is the maximum GTO L used in the calculations. The rest of the Boys function parameters have been set-up to give fully accurate results for m up to mmax=100. The value 24 below is just a sensible value that will be enough for most calculations. If the integral algorithm fails complaining that mmax is too small then set it here to the value at least 4*L where is the maximum L in your GTO atomic basis.

◆ nam_scattering_centre

character(len=nuc_nam_len) const_gbl::nam_scattering_centre = 'sc'

Name of the scattering centre.

◆ no_header

character(len=19), parameter const_gbl::no_header = 'No header specified'

Default character string used for the int_input_outputname variable.

◆ nuc_nam_len

integer, parameter const_gbl::nuc_nam_len = 2

Length of the character variable 'name' in the nucleus_type object.

◆ nuc_rep_att_ints

character(len=*), parameter const_gbl::nuc_rep_att_ints = 'Nuclear attraction integrals'

◆ number_of_types_el_ints

integer, parameter const_gbl::number_of_types_el_ints = 5

How many different types of 1-electron integrals have been implemented so far: Overlap (1), Kinetic energy (2), Kinetic energy+Nuclear attraction (3), Nuclear attraction (4), property integrals (5).

◆ one_electron_ints

character(len=*), parameter const_gbl::one_electron_ints = 'One electron integrals'

Character parameters used by the method l_mol_basismolecular_integrals to identify the type of molecular integral requested. These headers are used by the user to request the specific integral.

◆ one_elham

character(len=*), parameter const_gbl::one_elham = 'One electron Hamiltonian integrals'

◆ one_p_sym_ints

character(len=*), parameter const_gbl::one_p_sym_ints = 'One-particle AO->MO transformed integrals for symmetric operators'

◆ orbs_line

integer, parameter const_gbl::orbs_line = 6

Number of molecular orbitals for which orbital coefficients will be printed on one line by the routine molecular_orbital_basis_objprint_orbitals.

◆ overlap_ints

character(len=*), parameter const_gbl::overlap_ints = 'Overlap integrals'

◆ pg_irr_names

character(len=sym_op_nam_len), dimension(8,8), parameter const_gbl::pg_irr_names = RESHAPE ((/'A ',' ',' ',' ',' ',' ',' ',' ', "A' ",'A" ',' ',' ',' ',' ',' ',' ', 'A ','B ',' ',' ',' ',' ',' ',' ', 'Ag ','Au ',' ',' ',' ',' ',' ',' ', 'A1 ','B1 ','B2 ','A2 ',' ',' ',' ',' ', 'Ag ','Au ','Bu ','Bg ',' ',' ',' ',' ', 'A ','B3 ','B2 ','B1 ',' ',' ',' ',' ', 'Ag ','B3u','B2u','B1g','B1u','B2g','B3g','Au '/), (/8,8/))

Combined irreducible representation labels used in Psi4 Molden files.

◆ property_ints

character(len=*), parameter const_gbl::property_ints = 'Property integrals'

◆ stderr

integer, parameter const_gbl::stderr = error_unit

unit for standard error

◆ stdin

integer, parameter const_gbl::stdin = input_unit

unit for standard input

◆ stdout

integer, protected const_gbl::stdout = output_unit

unit for standard output. This value is set by the routine redirect_stdout_to_file but it can be done only once during the course of the program. Typically, the value of stdout will be chosen different for each rank (in case of MPI runs) so that each rank will use its own unit for std. output. This also ensures that no modifications in the code are needed concerning standard output.

◆ stdout_file_name

character(len=*), parameter const_gbl::stdout_file_name = "log_file."

Leading part of the file-name that will be used to collect standard output for each process. Each process gets its own file for standard output, e.g. "log_file.2" will contain standard output from MPI process with rank = 2.

◆ stdout_unit_base

integer, parameter const_gbl::stdout_unit_base = 10000

The base number that will stand for the unit number for each processes's standard output (see mpi_mod). This number must be chosen large enough so it cannot collide with any other file units that may be open during the program run.

◆ stdout_unit_step

integer, parameter const_gbl::stdout_unit_step = 1

Parameter used in redirect_stdout_to_file to compute the unit number for each processes's standard output: stdout = stdout_unit_base+myrank*stdout_unit_step This parameter allows to choose e.g. master's stdout to be unit 6 while the other processes will not interfere with any other fort units, e.g.: stdout_unit_base=6, stdout_unit_step=1000.

◆ sym_op_nam_len

integer, parameter const_gbl::sym_op_nam_len = 3

Length of the character variable sym_op specifying the symmetry operation in the geometry_obj object.

◆ taylor_k_ep

integer, parameter const_gbl::taylor_k_ep = 10

Order of the Taylor expansion for Taylor-expansion-based evaluation of the Boys function Fm(T) for quad precision calculations.

◆ taylor_k_wp

integer, parameter const_gbl::taylor_k_wp = 8

Order of the Taylor expansion for Taylor-expansion-based evaluation of the Boys function Fm(T) for double precision calculations.

◆ thrs_cf_sym_ortho_trans

real(kind=cfp), parameter const_gbl::thrs_cf_sym_ortho_trans = 10**(-(precision(cfp_dummy)-5.0_cfp))

Coefficients in the transformation matrix for the symmetric orthogonalization smaller than this value will be neglected.

◆ thrs_gs_ortho

real(kind=cfp), parameter const_gbl::thrs_gs_ortho = 10**(-(precision(cfp_dummy)-5.0_cfp))

Maximum allowed value for a cross overlap of two orbitals after the Gramm-Schmidt orthogonalization has been performed.

◆ thrs_lin_dep_gs_ortho

real(kind=cfp), parameter const_gbl::thrs_lin_dep_gs_ortho = 10**(-(precision(cfp_dummy)-8.0_cfp))

Threshold value for self-overlap of an orbital (before normalization) for the Gramm-Schmidt orthogonalization. If an orbital has a self-overlap (before normalization) smaller than this value then we assume that linear dependency in the orbital basis is present. In other words we decide that self-overlaps smaller than this value would lead to numerical problems.

◆ thrs_orb_cf

real(kind=cfp), parameter const_gbl::thrs_orb_cf = 10**(-(precision(cfp_dummy)-3.0_cfp))

Threshold for magnitude of final orbital coefficients as used by molecular_orbital_basis_objdelete_small_coefficients.

◆ thrs_symm_ortho

real(kind=cfp), parameter const_gbl::thrs_symm_ortho = 10**(-(precision(cfp_dummy)-9.0_cfp))

Maximum allowed value for a cross overlap of two orbitals after the symmetric orthogonalization has been performed.

Warning
This value should be actually equal to the threshold value for the integrals.

◆ tile

integer, parameter const_gbl::tile = 64

Tile size used for cache blocking. This is derived from cacheline size.

Todo:
define l1_cache_size and use it to manage sizes of the buffers: check if the source-target buffers fit inside the cache. If not then reallocate them (if possible) to their minimal sizes.

◆ two_el_ints

character(len=*), parameter const_gbl::two_el_ints = 'Two-electron integrals'

◆ two_el_ints_prefetched

character(len=*), parameter const_gbl::two_el_ints_prefetched = 'Prefetched two-electron integrals'

◆ two_p_sym_ints

character(len=*), parameter const_gbl::two_p_sym_ints = 'Two-particle AO->MO transformed integrals for symmetric operators'

◆ ufmat

character(len_ufmat), parameter const_gbl::ufmat = 'unformatted'

◆ wp_fit_terms

real(kind=cfp), dimension(1:fit_order+1), parameter const_gbl::wp_fit_terms = (/-18.206_cfp,0.42734_cfp,0.00075795_cfp/)

◆ wp_fit_terms_up

real(kind=cfp), dimension(1:fit_order+1), parameter const_gbl::wp_fit_terms_up = (/-18.561_cfp,0.3748_cfp,0.00088246_cfp/)

The same as above but this time for the UPWARD recursion to full double precision accuracy.

◆ y_lm_size_threshold

real(kind=cfp), parameter const_gbl::y_lm_size_threshold = 100.0_cfp

If it is requested that the Y_lm functions (evaluated by eval_BTO_CGTO_Y_lm) are saved to disk then all Y_lm with size (in MiB) greater than this value will be saved to disk. If you want to save all Y_lm regardless of their size then set this parameter to 0.0_cfp.