GBTOlib: library for evaluation of molecular integrals in mixed Gaussian / B-spline basis 111
Loading...
Searching...
No Matches
ukrmol_interface_gbl Module Reference

Functions/Subroutines

subroutine, public START_MPI
 This is used to start MPI in SWINTERF and HAMDIAG. In the serial integral code this is just a dummy routine that does nothing.
subroutine, public FINALIZE_MPI
 This is used to start MPI in SWINTERF and HAMDIAG. In the serial integral code this is just a dummy routine that does nothing.
subroutine, public DESTROY_UKRMOLP ()
subroutine, public READ_UKRMOLP_BASIS (nftint)
 Reads-in the AO and MO basis sets from the fort file with a given number.
subroutine, public GET_GEOM (nnuc, cname, xnuc, ynuc, znuc, charge)
 Returns the nuclear data for the problem as contained in the molecular_orbital_basis object.
subroutine, public GET_NAME_SYM (name, no_sym, nob, nlmq)
 Returns the number of symmetries and some other basic parameters as obtained from the initialized object molecular_orbital_basis and property_integrals.
subroutine, public READ_UKRMOLP_INTS (nfte, nfti, lembf, nint1e, nint2e, nocsf, nfta, isymtp, nsym1, nob1, iposit, scalem, name, nalm, qmoln)
 This routine reads-in all one electron and two electron integrals needed for SCATCI calculations. It is an analogoue of the INTSIN routine from SCATCI.
subroutine, public READ_UKRMOLP_PROPERTY_INTS (nftint, iwrite)
 Reads-in the basis sets and the property integrals.
real(kind=wp) function, public GET_INTEGRAL (a, b, c, d, positron_flag)
 Returns a one electron/positron integral (K+-V) or two-electron integral given the basis function indices and the positron flag. If c,d == 0 then the one electron/positron integral is returned. If positron_flag > 0 then the positronic one-particle integral (K-V) is returned.
real(kind=wp) function, public GET_ORB_ENERGY (a)
real(kind=wp) function, public GET_KINETIC_ENERGY_INTEGRAL (a, b)
 Function to return the kinetic energy integral as these integrals are used in the calculation of BEB cross sections.
real(wp) function, public GET_PROPERTY_INTEGRAL (iprop, a, b, iposit)
 Retrieve value of a one-particle property integral.
subroutine, public TMG_UKPLUS (mdel, block_data, no_blocks, nob, den, lden, prlmq, lmq, l, m, q)
 Multiplies the density matrix elements with the set of property integrals with (l,m) given by the linear index lmq = l*l+l+m+1. This routine may only be called following a call to READ_UKRMOLP_PROPERTY_INTS. We assume that this routine is called with the first call having lmq=1. This is needed for the index function initialization. The input values are: the wavefunction symmetry difference (mdel), the number of elements in the density matrix (lden), the property l,m,q sequence number (lmq) and the block data constructed in DENPROP/TMGP. The output are: the l,m,q values corresponding to lmq and the value for the electronic contribution to the property (prlmq).
subroutine, public construct_pintegrals (nftprop, iwrite, no_of_integrals, no_of_properties, nilmq, lp, mp, qp, property_name, nob, mob, mpob)
 This routine interfaces with CDENPROP. It transfers the necessary orbital data into the argument variables. We use READ_UKRMOLP_PROPERTY_INTS to read the integrals from disk and they are ready for retrieval using GET_PROPERTY_INTEGRAL.
subroutine, public write_molden_dyson_orbitals (dyson_orbitals, idyson_orbital_irrep, itarget_irrep, itarget_spin, dyson_orbital_norms, molden_format_opt)
 This routine takes the Dyson orbitals produced by CDENPROP and writes them out in the Molden file format along with the GTO basis set. The CDENPROP Dyson's are expressed as a linear combination of the target and continuum MOs. Additionally to producing the Molden file this routine always saves the complete set of Dyson orbitals (along with the AO basis set) into an object of orbital_data_obj type and then writes this to disk. This ensures that we always retain the full information about the Dyson orbitals (in case the AO basis contains high L that Molden cannot handle). The value of the R-matrix radius stored in 'options' is used to normalize the continuum functions to the R-matrix radius. These normalization coefficients are multiplied in with the orbital coefficients to produce the Molden file. The orbital coefficients saved in the UKRmol+ format are modified, too. Therefore when calculating the radial densities of the Dysons from the UKRmol+ file using the routine orbital_basis_data_objradial_charge_density the input value of rmat_radius for this routine must be set to a value .le. 0.0_cfp. See also orbital_basis_data_objradial_charge_density.
subroutine, public EVAL_AMPLITUDES (a, normalize_to_a)
 Calculates the orbital amplitudes and determines the channel information. The result is stored in the module private variables amplitudes,continuum_channels_m_l_irr for later use by UKP_PREAMP, UKP_READAMP. Call to this routine must be preceeded by calling e.g. READ_UKRMOLP_BASIS in order to initialize the molecular_orbital_basis object which performs the amplitude evaluation.
subroutine, public UKP_PREAMP (irrcont, ix, lchl, mchl, qchl, nch, maxnmo)
 Transfers the channel angular and IRR numbers to the channel variables from SWINTERF. It is the equivalent of the PREAMP routine from SWINTERF. The call to this routine must be preceeded by call to EVAL_AMPLITUDES.
subroutine, public UKP_READAMP (wamps, nchan, irrchl, lchl, mchl, qchl, ncontcsf, mcont, iprint)
 This routine transfers into WAMPS the R-matrix amplitudes for the orbitals. It is equivalent of the READAMP routine. Note that the orbital amplitudes are divided by sqrt(2). See the routine outerio/READRM routine and the equation for the R-matrix there. It implies that the 1/2 factor from the R-matrix formula is absorbed into the product of the two amplitudes in the summation over the R-matrix poles. Hence the factor 1/sqrt(2) by which the amplitudes must be multiplied. The call to this routine must be preceeded by call to UKP_PREAMP.
subroutine, public GET_ECP_CORE_CHARGES (core_charges)
subroutine, public READ_ADC_INTS (nfti, nfta, nint1e, nint2e, nsym)
 This routine reads-in all two electron integrals needed for ADC codes. It is derived from the READ_UKRMOLP_INTS – parameter interface is significantly simplified and redundant/UKRmol-specific functionality is removed. no positrons are expected, corresponding switch parameter (iposit) removed.
subroutine, public GET_ENERSYMOCC (pgroup, nsym, nbas, nocc, ener, sym, occ, iscont, mt)
 Returns molecular orbital data needed for the construction of ADC Hamiltonian.
real(kind=wp) function, public LPQRS (p, q, r, s)
 Returns double-precision two-electron integral given the four INTEGER(2) indices p,q,r,s – required by the ADC codes; one-electron integrals or positron flag not supported.

Variables

real(kind=cfp), dimension(:,:), allocatable, public fock
 These two variables are used in SCATCI_ROUTINES to identify the format of the integral file on input.
type(molecular_orbital_basis_obj), public molecular_orbital_basis
 Basis set of molecular orbitals. Read-in by READ_UKRMOLP_BASIS.

Function/Subroutine Documentation

◆ construct_pintegrals()

subroutine, public ukrmol_interface_gbl::construct_pintegrals ( integer, intent(in) nftprop,
integer, intent(in) iwrite,
integer, intent(out) no_of_integrals,
integer, intent(out) no_of_properties,
integer, dimension(:), allocatable nilmq,
integer, dimension(:), allocatable lp,
integer, dimension(:), allocatable mp,
integer, dimension(:), allocatable qp,
character(len=8), dimension(:), allocatable property_name,
integer, dimension(:), allocatable nob,
integer, dimension(:), allocatable mob,
integer, dimension(:), allocatable mpob )

This routine interfaces with CDENPROP. It transfers the necessary orbital data into the argument variables. We use READ_UKRMOLP_PROPERTY_INTS to read the integrals from disk and they are ready for retrieval using GET_PROPERTY_INTEGRAL.

◆ DESTROY_UKRMOLP()

subroutine, public ukrmol_interface_gbl::DESTROY_UKRMOLP

◆ EVAL_AMPLITUDES()

subroutine, public ukrmol_interface_gbl::EVAL_AMPLITUDES ( real(kind=wp), intent(in) a,
logical, intent(in) normalize_to_a )

Calculates the orbital amplitudes and determines the channel information. The result is stored in the module private variables amplitudes,continuum_channels_m_l_irr for later use by UKP_PREAMP, UKP_READAMP. Call to this routine must be preceeded by calling e.g. READ_UKRMOLP_BASIS in order to initialize the molecular_orbital_basis object which performs the amplitude evaluation.

Warning
This routine ALWAYS requires double precision value of the R-matrix radius on input due to the need for compatibility with UKRmol-out which is always compiled using double precision.

◆ FINALIZE_MPI()

subroutine, public ukrmol_interface_gbl::FINALIZE_MPI

This is used to start MPI in SWINTERF and HAMDIAG. In the serial integral code this is just a dummy routine that does nothing.

Here is the call graph for this function:

◆ GET_ECP_CORE_CHARGES()

subroutine, public ukrmol_interface_gbl::GET_ECP_CORE_CHARGES ( integer, dimension(:), allocatable core_charges)

◆ GET_ENERSYMOCC()

subroutine, public ukrmol_interface_gbl::GET_ENERSYMOCC ( integer(kind=shortint), intent(out) pgroup,
integer(kind=shortint), intent(in) nsym,
integer(kind=shortint), intent(out) nbas,
integer(kind=shortint), intent(out) nocc,
real(kind=wp), dimension(:), intent(out), allocatable ener,
integer(kind=shortint), dimension(:), intent(out), allocatable sym,
real(kind=wp), dimension(:), intent(out), allocatable occ,
logical, dimension(:), intent(out), allocatable iscont,
integer(kind=shortint), dimension(8,8), intent(out) mt )

Returns molecular orbital data needed for the construction of ADC Hamiltonian.

◆ GET_GEOM()

subroutine, public ukrmol_interface_gbl::GET_GEOM ( integer, intent(out) nnuc,
character(len=8), dimension(:), intent(out), allocatable cname,
real(kind=wp), dimension(:), intent(out), allocatable xnuc,
real(kind=wp), dimension(:), intent(out), allocatable ynuc,
real(kind=wp), dimension(:), intent(out), allocatable znuc,
real(kind=wp), dimension(:), intent(out), allocatable charge )

Returns the nuclear data for the problem as contained in the molecular_orbital_basis object.

Warning
This routine ALWAYS returns double precision values due to the need for compatibility with UKRmol-in/out which is always compiled using double precision.

◆ GET_INTEGRAL()

real(kind=wp) function, public ukrmol_interface_gbl::GET_INTEGRAL ( integer, intent(in) a,
integer, intent(in) b,
integer, intent(in) c,
integer, intent(in) d,
integer, intent(in) positron_flag )

Returns a one electron/positron integral (K+-V) or two-electron integral given the basis function indices and the positron flag. If c,d == 0 then the one electron/positron integral is returned. If positron_flag > 0 then the positronic one-particle integral (K-V) is returned.

Warning
This routine ALWAYS returns double precision value due to the need for compatibility with UKRmol-in/out which is always compiled using double precision. NEW FEATURE: WORKS WITH SHARED MEMORY OF THE MPI-3.0 STANDARD

◆ GET_KINETIC_ENERGY_INTEGRAL()

real(kind=wp) function, public ukrmol_interface_gbl::GET_KINETIC_ENERGY_INTEGRAL ( integer, intent(in) a,
integer, intent(in) b )

Function to return the kinetic energy integral as these integrals are used in the calculation of BEB cross sections.

Authors
Bridgette Cooper
Date
October 2019
Warning
This routine ALWAYS returns double precision value due to the need for compatibility with UKRmol-in/out which is always compiled using double precision.

◆ GET_NAME_SYM()

subroutine, public ukrmol_interface_gbl::GET_NAME_SYM ( character(len=132), intent(out) name,
integer, intent(out) no_sym,
integer, dimension(:), intent(out) nob,
integer, intent(out) nlmq )

Returns the number of symmetries and some other basic parameters as obtained from the initialized object molecular_orbital_basis and property_integrals.

◆ GET_ORB_ENERGY()

real(kind=wp) function, public ukrmol_interface_gbl::GET_ORB_ENERGY ( integer, intent(in) a)

◆ GET_PROPERTY_INTEGRAL()

real(wp) function, public ukrmol_interface_gbl::GET_PROPERTY_INTEGRAL ( integer, intent(in) iprop,
integer, intent(in) a,
integer, intent(in) b,
integer, optional iposit )

Retrieve value of a one-particle property integral.

Authors
J Benda
Date
2023

Assuming that the subroutine READ_UKRMOLP_PROPERTY_INTS has been called before, retrieve the one-electron property integral. This function is currently used in CDENPROP.

Parameters
[in]ipropIndex of the property (1 = s0, 2 = py, 3 = pz, 4 = px, etc)
[in]aFirst molecular orbital index (continuous index in concatenated symmetries)
[in]bSecond molecular orbital index (continuous index in concatenated symmetries)
[in]ipositReturn electronic property (if absent or zero) or positronic one (otherwise)
Returns
Value of the integal as a double precision real number.

◆ LPQRS()

real(kind=wp) function, public ukrmol_interface_gbl::LPQRS ( integer(kind=shortint2), intent(in) p,
integer(kind=shortint2), intent(in) q,
integer(kind=shortint2), intent(in) r,
integer(kind=shortint2), intent(in) s )

Returns double-precision two-electron integral given the four INTEGER(2) indices p,q,r,s – required by the ADC codes; one-electron integrals or positron flag not supported.

◆ READ_ADC_INTS()

subroutine, public ukrmol_interface_gbl::READ_ADC_INTS ( integer, intent(in) nfti,
integer, intent(in) nfta,
integer, intent(inout) nint1e,
integer, intent(inout) nint2e,
integer(kind=shortint), intent(inout) nsym )

This routine reads-in all two electron integrals needed for ADC codes. It is derived from the READ_UKRMOLP_INTS – parameter interface is significantly simplified and redundant/UKRmol-specific functionality is removed. no positrons are expected, corresponding switch parameter (iposit) removed.


  • ROUTINES PROVIDING SIMPLIFIED INTERFACE FOR ADC CODES *

  • Warning
    This routine ALWAYS requires double precision value of scalem on input due to the need for compatibility with UKRmol-in/out which is always compiled using double precision. Note that this results in reducing the values in ke_integralsa into double precision. NEW FEATURE: WORKS WITH SHARED MEMORY OF THE MPI-3.0 STANDARD. This means that all integrals are read-in the shared mode so in the MPI regime there is only one copy of the integral arrays PER NODE. If we do not have the shared memory capability then the code resorts to the local mode, i.e. each MPI task keeps its own copy of all integral arrays.
Here is the call graph for this function:

◆ READ_UKRMOLP_BASIS()

subroutine, public ukrmol_interface_gbl::READ_UKRMOLP_BASIS ( integer, intent(in) nftint)

Reads-in the AO and MO basis sets from the fort file with a given number.

◆ READ_UKRMOLP_INTS()

subroutine, public ukrmol_interface_gbl::READ_UKRMOLP_INTS ( integer, intent(in) nfte,
integer, intent(in) nfti,
integer, intent(in) lembf,
integer, intent(inout) nint1e,
integer, intent(inout) nint2e,
integer, intent(in) nocsf,
integer, intent(in) nfta,
integer, intent(in) isymtp,
integer, intent(in) nsym1,
integer, dimension(nsym1), intent(in) nob1,
integer, intent(in) iposit,
real(kind=wp), intent(in) scalem,
character(len=120), intent(in) name,
integer, intent(inout) nalm,
logical, intent(in) qmoln )

This routine reads-in all one electron and two electron integrals needed for SCATCI calculations. It is an analogoue of the INTSIN routine from SCATCI.

Warning
This routine ALWAYS requires double precision value of scalem on input due to the need for compatibility with UKRmol-in/out which is always compiled using double precision. Note that this results in reducing the values in ke_integralsa into double precision. NEW FEATURE: WORKS WITH SHARED MEMORY OF THE MPI-3.0 STANDARD. This means that all integrals are read-in the shared mode so in the MPI regime there is only one copy of the integral arrays PER NODE. If we do not have the shared memory capability then the code resorts to the local mode, i.e. each MPI task keeps its own copy of all integral arrays.

◆ READ_UKRMOLP_PROPERTY_INTS()

subroutine, public ukrmol_interface_gbl::READ_UKRMOLP_PROPERTY_INTS ( integer, intent(in) nftint,
integer, intent(in) iwrite )

Reads-in the basis sets and the property integrals.

◆ START_MPI()

subroutine, public ukrmol_interface_gbl::START_MPI

This is used to start MPI in SWINTERF and HAMDIAG. In the serial integral code this is just a dummy routine that does nothing.

Here is the call graph for this function:

◆ TMG_UKPLUS()

subroutine, public ukrmol_interface_gbl::TMG_UKPLUS ( integer, intent(in) mdel,
integer, dimension(:,:), intent(in) block_data,
integer, intent(in) no_blocks,
integer, dimension(:), intent(in) nob,
real(kind=wp), dimension(lden) den,
integer, intent(in) lden,
real(kind=wp) prlmq,
integer, intent(in) lmq,
integer, intent(out) l,
integer, intent(out) m,
integer, intent(out) q )

Multiplies the density matrix elements with the set of property integrals with (l,m) given by the linear index lmq = l*l+l+m+1. This routine may only be called following a call to READ_UKRMOLP_PROPERTY_INTS. We assume that this routine is called with the first call having lmq=1. This is needed for the index function initialization. The input values are: the wavefunction symmetry difference (mdel), the number of elements in the density matrix (lden), the property l,m,q sequence number (lmq) and the block data constructed in DENPROP/TMGP. The output are: the l,m,q values corresponding to lmq and the value for the electronic contribution to the property (prlmq).

Warning
This routine ALWAYS uses double precision values of the property integrals due to the need for compatibility with UKRmol-in/out which is always compiled using double precision.

◆ UKP_PREAMP()

subroutine, public ukrmol_interface_gbl::UKP_PREAMP ( integer, intent(in) irrcont,
integer, intent(in) ix,
integer, dimension(*), intent(out) lchl,
integer, dimension(*), intent(out) mchl,
integer, dimension(*), intent(out) qchl,
integer, intent(out) nch,
integer, intent(inout) maxnmo )

Transfers the channel angular and IRR numbers to the channel variables from SWINTERF. It is the equivalent of the PREAMP routine from SWINTERF. The call to this routine must be preceeded by call to EVAL_AMPLITUDES.

◆ UKP_READAMP()

subroutine, public ukrmol_interface_gbl::UKP_READAMP ( real(kind=wp), dimension(nchan,*) wamps,
integer nchan,
integer, dimension(*) irrchl,
integer, dimension(*) lchl,
integer, dimension(*) mchl,
integer, dimension(*) qchl,
integer, dimension(8) ncontcsf,
integer, dimension(:) mcont,
integer iprint )

This routine transfers into WAMPS the R-matrix amplitudes for the orbitals. It is equivalent of the READAMP routine. Note that the orbital amplitudes are divided by sqrt(2). See the routine outerio/READRM routine and the equation for the R-matrix there. It implies that the 1/2 factor from the R-matrix formula is absorbed into the product of the two amplitudes in the summation over the R-matrix poles. Hence the factor 1/sqrt(2) by which the amplitudes must be multiplied. The call to this routine must be preceeded by call to UKP_PREAMP.

Warning
This routine ALWAYS outputs double precision values of the amplitudes in WAMPS due to the need for compatibility with UKRmol-in/out which is always compiled using double precision.

◆ write_molden_dyson_orbitals()

subroutine, public ukrmol_interface_gbl::write_molden_dyson_orbitals ( real(kind=wp), dimension(:,:,:), intent(in) dyson_orbitals,
integer, dimension(:), intent(in) idyson_orbital_irrep,
integer, dimension(:), intent(in) itarget_irrep,
integer, dimension(:), intent(in) itarget_spin,
real(kind=wp), dimension(:,:), intent(in) dyson_orbital_norms,
integer, intent(in), optional molden_format_opt )

This routine takes the Dyson orbitals produced by CDENPROP and writes them out in the Molden file format along with the GTO basis set. The CDENPROP Dyson's are expressed as a linear combination of the target and continuum MOs. Additionally to producing the Molden file this routine always saves the complete set of Dyson orbitals (along with the AO basis set) into an object of orbital_data_obj type and then writes this to disk. This ensures that we always retain the full information about the Dyson orbitals (in case the AO basis contains high L that Molden cannot handle). The value of the R-matrix radius stored in 'options' is used to normalize the continuum functions to the R-matrix radius. These normalization coefficients are multiplied in with the orbital coefficients to produce the Molden file. The orbital coefficients saved in the UKRmol+ format are modified, too. Therefore when calculating the radial densities of the Dysons from the UKRmol+ file using the routine orbital_basis_data_objradial_charge_density the input value of rmat_radius for this routine must be set to a value .le. 0.0_cfp. See also orbital_basis_data_objradial_charge_density.

By default, the Molden file will use Cartesian basis. This can be overriden by setting the optional parameter molden_format_opt to 2, which will instead produce a spherical Molden file.

Warning
This routine multiplies in the spin-related CG coefficient (see commit 879 to cdenprop_io) but it is hard coded to the value 1/sqrt(2). This must be changed for your particular case.
This routine ALWAYS REQUIRES double precision values of the coefficients on input due to the need for compatibility with CDENPROP which is always compiled using double precision.
Here is the call graph for this function:

Variable Documentation

◆ fock

real(kind=cfp), dimension(:,:), allocatable, public ukrmol_interface_gbl::fock

These two variables are used in SCATCI_ROUTINES to identify the format of the integral file on input.

The following routines are used by UKRmol-in an UKRmol-out programs (SCATCI, DENPROP, CDENPROP) and routines (SWINTERF). ADC-RELATED ROUTINES AND VARIABLES – START Routines to read basis information and two-electron integrals as needed in ADC codes Returns two-electron integrals – simplifies GET_INTEGRAL, compatible with integer(2) indices in ADC codes Fock matrix

◆ molecular_orbital_basis

type(molecular_orbital_basis_obj), public ukrmol_interface_gbl::molecular_orbital_basis

Basis set of molecular orbitals. Read-in by READ_UKRMOLP_BASIS.