GBTOlib: library for evaluation of molecular integrals in mixed Gaussian / B-spline basis
111
|
Public Member Functions | |
procedure | check => check_integral_options_obj |
Checks that the precision values are within a reasonable range. More... | |
procedure | read => read_integral_options_obj |
Reads the integral_options data structure from the disk given the starting position of the record. More... | |
procedure | write => write_integral_options_obj |
Writes the integral_options data structure to the disk given the starting position of the record. More... | |
procedure | print => print_integral_options_obj |
Prints to screen the contents of this object. More... | |
Public Attributes | |
real(kind=cfp) | prec = int_rel_prec |
Rel. precision for the resulting integrals. Error will be triggered if there is a suspicion that the integral fails to meet the precision specified. This is implemented for the GTO/BTO integrals. More... | |
real(kind=cfp) | tol = int_del_thr |
Threshold magnitude for the integrals: integrals with absolute values smaller than this will be zeroed out. More... | |
real(kind=cfp) | a = 0.0_cfp |
R-matrix radius. This value must be always consistent with the extent of the BTO basis. The last endpoint of the radial B-spline grid must be the value of the R-matrix radius. If no BTOs are present in the continuum basis then this value can be set to an arbitrary value. More... | |
real(kind=cfp) | max_ijrs_size = -1.0_cfp |
Maximum allowed size (in Mib) of the ijrs as allocated by molecular_orbital_basis_objtwo_electron_integrals. More... | |
logical | two_p_continuum = .false. |
If set to .true. then 2p integrals of the type [continuum,continuum|continuum,continuum] will be calculated. If set to .false. (default) then 2p integrals with only one particle in the continuum will be calculated. More... | |
logical | use_spherical_cgto_alg = .false. |
Set to .true. to use the CGTO-only integral algorithm based on the spherical CGTOs. This is numerically stable even for high L. More... | |
integer | mixed_ints_method = -1 |
Method to calculate the mixed BTO/CGTO integrals: 1 = Legendre expansion, 2 = Lebedev quadrature. todo broadcast it too in the parallel calculation! More... | |
character(len=line_len) | scratch_directory = '' |
If this variable is set then all intermediate integrals Y_lm needed for calculation of the BTO/GTO integrals will be saved to disk in the directory specified. More... | |
integer | max_l_legendre_1el = -1 |
Maximum value of L to be used in the Legendre expansion when calculating the 1-electron, resp. 2-electron mixed integrals. todo broadcast it too in the parallel calculation! More... | |
integer | max_l_legendre_2el = -1 |
real(kind=cfp) | delta_r1 = 0.25_cfp |
Grid step for the r1 coordinate used in the calculation of the mixed BTO/GTO integrals. More... | |
integer | max_property_l = -1 |
The maximum value of the angular momentum of the property operator. The default is 2, i.e. property integrals up to quadrupoles will be calculated. More... | |
logical | calculate_overlap_ints = .false. |
Set to .true. if calculation of overlap and kinetic energy integrals is required. More... | |
logical | calculate_kinetic_energy_ints = .false. |
Set to .true. if calculation of kinetic energy integrals is required. More... | |
logical | calculate_nuclear_attraction_ints = .false. |
Set to .true. if calculation of nuclear attraction integrals is required. More... | |
logical | calculate_one_el_hamiltonian_ints = .false. |
Set to .true. if calculation of one electron Hamiltonian integrals are required (sum of kinetic energy and nuclear attraction integrals). More... | |
logical | calculate_property_ints = .false. |
Set to .true. if calculation of property integrals (with max_property_l) is required. More... | |
logical | calculate_two_el_ints = .false. |
Set to .true. if calculation of two electron integrals is required. More... | |
logical | print_integrals = .false. |
Request to print the evaluated integrals at the end of the calculation. This parameter is not saved to/read from the disk. More... | |
real(kind=cfp) | dipole_damp_factor = 0.0_cfp |
If non-zero the radial part of the dipole operator will be replaced from 'r' to 'r * exp(-dipole_damp_factor * r)'. More... | |
procedure integral_storage_gbl::integral_options_obj::check |
Checks that the precision values are within a reasonable range.
procedure integral_storage_gbl::integral_options_obj::print |
Prints to screen the contents of this object.
procedure integral_storage_gbl::integral_options_obj::read |
Reads the integral_options data structure from the disk given the starting position of the record.
procedure integral_storage_gbl::integral_options_obj::write |
Writes the integral_options data structure to the disk given the starting position of the record.
real(kind=cfp) integral_storage_gbl::integral_options_obj::a = 0.0_cfp |
R-matrix radius. This value must be always consistent with the extent of the BTO basis. The last endpoint of the radial B-spline grid must be the value of the R-matrix radius. If no BTOs are present in the continuum basis then this value can be set to an arbitrary value.
logical integral_storage_gbl::integral_options_obj::calculate_kinetic_energy_ints = .false. |
Set to .true. if calculation of kinetic energy integrals is required.
logical integral_storage_gbl::integral_options_obj::calculate_nuclear_attraction_ints = .false. |
Set to .true. if calculation of nuclear attraction integrals is required.
logical integral_storage_gbl::integral_options_obj::calculate_one_el_hamiltonian_ints = .false. |
Set to .true. if calculation of one electron Hamiltonian integrals are required (sum of kinetic energy and nuclear attraction integrals).
logical integral_storage_gbl::integral_options_obj::calculate_overlap_ints = .false. |
Set to .true. if calculation of overlap and kinetic energy integrals is required.
logical integral_storage_gbl::integral_options_obj::calculate_property_ints = .false. |
Set to .true. if calculation of property integrals (with max_property_l) is required.
logical integral_storage_gbl::integral_options_obj::calculate_two_el_ints = .false. |
Set to .true. if calculation of two electron integrals is required.
real(kind=cfp) integral_storage_gbl::integral_options_obj::delta_r1 = 0.25_cfp |
Grid step for the r1 coordinate used in the calculation of the mixed BTO/GTO integrals.
real(kind=cfp) integral_storage_gbl::integral_options_obj::dipole_damp_factor = 0.0_cfp |
If non-zero the radial part of the dipole operator will be replaced from 'r' to 'r * exp(-dipole_damp_factor * r)'.
real(kind=cfp) integral_storage_gbl::integral_options_obj::max_ijrs_size = -1.0_cfp |
Maximum allowed size (in Mib) of the ijrs as allocated by molecular_orbital_basis_objtwo_electron_integrals.
integer integral_storage_gbl::integral_options_obj::max_l_legendre_1el = -1 |
Maximum value of L to be used in the Legendre expansion when calculating the 1-electron, resp. 2-electron mixed integrals. todo broadcast it too in the parallel calculation!
integer integral_storage_gbl::integral_options_obj::max_l_legendre_2el = -1 |
integer integral_storage_gbl::integral_options_obj::max_property_l = -1 |
The maximum value of the angular momentum of the property operator. The default is 2, i.e. property integrals up to quadrupoles will be calculated.
integer integral_storage_gbl::integral_options_obj::mixed_ints_method = -1 |
Method to calculate the mixed BTO/CGTO integrals: 1 = Legendre expansion, 2 = Lebedev quadrature. todo broadcast it too in the parallel calculation!
real(kind=cfp) integral_storage_gbl::integral_options_obj::prec = int_rel_prec |
Rel. precision for the resulting integrals. Error will be triggered if there is a suspicion that the integral fails to meet the precision specified. This is implemented for the GTO/BTO integrals.
logical integral_storage_gbl::integral_options_obj::print_integrals = .false. |
Request to print the evaluated integrals at the end of the calculation. This parameter is not saved to/read from the disk.
character(len=line_len) integral_storage_gbl::integral_options_obj::scratch_directory = '' |
If this variable is set then all intermediate integrals Y_lm needed for calculation of the BTO/GTO integrals will be saved to disk in the directory specified.
real(kind=cfp) integral_storage_gbl::integral_options_obj::tol = int_del_thr |
Threshold magnitude for the integrals: integrals with absolute values smaller than this will be zeroed out.
logical integral_storage_gbl::integral_options_obj::two_p_continuum = .false. |
If set to .true. then 2p integrals of the type [continuum,continuum|continuum,continuum] will be calculated. If set to .false. (default) then 2p integrals with only one particle in the continuum will be calculated.
logical integral_storage_gbl::integral_options_obj::use_spherical_cgto_alg = .false. |
Set to .true. to use the CGTO-only integral algorithm based on the spherical CGTOs. This is numerically stable even for high L.