(git:374b731)
Loading...
Searching...
No Matches
qs_integral_utils.F
Go to the documentation of this file.
1!--------------------------------------------------------------------------------------------------!
2! CP2K: A general program to perform molecular dynamics simulations !
3! Copyright 2000-2024 CP2K developers group <https://cp2k.org> !
4! !
5! SPDX-License-Identifier: GPL-2.0-or-later !
6!--------------------------------------------------------------------------------------------------!
7
8! **************************************************************************************************
9!> \brief Some utility functions for the calculation of integrals
10!> \par History
11!> JGH: initial version
12!> \author JGH (10.07.2014)
13! **************************************************************************************************
15
19 USE qs_kind_types, ONLY: get_qs_kind,&
22#include "./base/base_uses.f90"
23
24 IMPLICIT NONE
25
26 PRIVATE
27
28! *** Global parameters ***
29
30 CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'qs_integral_utils'
31
32! *** Interfaces ***
33
35 MODULE PROCEDURE get_memory_usage_a, get_memory_usage_ab, &
36 get_memory_usage_abc, get_memory_usage_abcd
37 END INTERFACE
38
39! *** Public subroutines ***
40
42
43CONTAINS
44
45! **************************************************************************************************
46!> \brief Return the maximum memory usage in integral calculations
47!> \param qs_kind_set The info for all atomic kinds
48!> \param basis_type_a Type of basis
49!> \return Result
50! **************************************************************************************************
51 FUNCTION get_memory_usage_a(qs_kind_set, basis_type_a) RESULT(ldmem)
52
53 TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
54 CHARACTER(LEN=*), INTENT(IN) :: basis_type_a
55 INTEGER :: ldmem
56
57 INTEGER :: maxc, maxl, maxs
58
59 CALL get_qs_kind_set(qs_kind_set=qs_kind_set, &
60 maxco=maxc, maxlgto=maxl, maxsgf=maxs, &
61 basis_type=basis_type_a)
62 ldmem = max(maxc, maxs)
63
64 CALL init_orbital_pointers(maxl + 2)
65
66 END FUNCTION get_memory_usage_a
67
68! **************************************************************************************************
69!> \brief Return the maximum memory usage in integral calculations
70!> \param qs_kind_set The info for all atomic kinds
71!> \param basis_type_a Type of basis
72!> \param basis_type_b Type of basis
73!> \return Result
74! **************************************************************************************************
75 FUNCTION get_memory_usage_ab(qs_kind_set, basis_type_a, basis_type_b) RESULT(ldmem)
76
77 TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
78 CHARACTER(LEN=*), INTENT(IN) :: basis_type_a, basis_type_b
79 INTEGER :: ldmem
80
81 INTEGER :: lda, ldb
82
83 lda = get_memory_usage_a(qs_kind_set, basis_type_a)
84 ldb = get_memory_usage_a(qs_kind_set, basis_type_b)
85 ldmem = max(lda, ldb)
86
87 END FUNCTION get_memory_usage_ab
88
89! **************************************************************************************************
90!> \brief Return the maximum memory usage in integral calculations
91!> \param qs_kind_set The info for all atomic kinds
92!> \param basis_type_a Type of basis
93!> \param basis_type_b Type of basis
94!> \param basis_type_c Type of basis
95!> \return Result
96! **************************************************************************************************
97 FUNCTION get_memory_usage_abc(qs_kind_set, basis_type_a, &
98 basis_type_b, basis_type_c) RESULT(ldmem)
99
100 TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
101 CHARACTER(LEN=*), INTENT(IN) :: basis_type_a, basis_type_b, basis_type_c
102 INTEGER :: ldmem
103
104 INTEGER :: lda, ldb, ldc
105
106 lda = get_memory_usage_a(qs_kind_set, basis_type_a)
107 ldb = get_memory_usage_a(qs_kind_set, basis_type_b)
108 ldc = get_memory_usage_a(qs_kind_set, basis_type_c)
109 ldmem = max(lda, ldb, ldc)
110
111 END FUNCTION get_memory_usage_abc
112
113! **************************************************************************************************
114!> \brief Return the maximum memory usage in integral calculations
115!> \param qs_kind_set The info for all atomic kinds
116!> \param basis_type_a Type of basis
117!> \param basis_type_b Type of basis
118!> \param basis_type_c Type of basis
119!> \param basis_type_d Type of basis
120!> \return Result
121! **************************************************************************************************
122 FUNCTION get_memory_usage_abcd(qs_kind_set, basis_type_a, &
123 basis_type_b, basis_type_c, basis_type_d) RESULT(ldmem)
124
125 TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
126 CHARACTER(LEN=*), INTENT(IN) :: basis_type_a, basis_type_b, &
127 basis_type_c, basis_type_d
128 INTEGER :: ldmem
129
130 INTEGER :: lda, ldb, ldc, ldd
131
132 lda = get_memory_usage_a(qs_kind_set, basis_type_a)
133 ldb = get_memory_usage_a(qs_kind_set, basis_type_b)
134 ldc = get_memory_usage_a(qs_kind_set, basis_type_c)
135 ldd = get_memory_usage_a(qs_kind_set, basis_type_d)
136 ldmem = max(lda, ldb, ldc, ldd)
137
138 END FUNCTION get_memory_usage_abcd
139
140! **************************************************************************************************
141
142! **************************************************************************************************
143!> \brief Set up an easy accessible list of the basis sets for all kinds
144!> \param basis_set_list The basis set list
145!> \param basis_type ...
146!> \param qs_kind_set Kind information, the basis is used
147! **************************************************************************************************
148 SUBROUTINE basis_set_list_setup(basis_set_list, basis_type, qs_kind_set)
149
150 TYPE(gto_basis_set_p_type), DIMENSION(:) :: basis_set_list
151 CHARACTER(len=*), INTENT(IN) :: basis_type
152 TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
153
154 INTEGER :: ikind
155 TYPE(gto_basis_set_type), POINTER :: basis_set
156 TYPE(qs_kind_type), POINTER :: qs_kind
157
158 ! set up basis sets
159 DO ikind = 1, SIZE(qs_kind_set)
160 qs_kind => qs_kind_set(ikind)
161 CALL get_qs_kind(qs_kind=qs_kind, basis_set=basis_set, &
162 basis_type=basis_type)
163 NULLIFY (basis_set_list(ikind)%gto_basis_set)
164 IF (ASSOCIATED(basis_set)) basis_set_list(ikind)%gto_basis_set => basis_set
165 END DO
166
167 END SUBROUTINE basis_set_list_setup
168
169! **************************************************************************************************
170
171END MODULE qs_integral_utils
172
Provides Cartesian and spherical orbital pointers and indices.
subroutine, public init_orbital_pointers(maxl)
Initialize or update the orbital pointers.
Some utility functions for the calculation of integrals.
subroutine, public basis_set_list_setup(basis_set_list, basis_type, qs_kind_set)
Set up an easy accessible list of the basis sets for all kinds.
Define the quickstep kind type and their sub types.
subroutine, public get_qs_kind(qs_kind, basis_set, basis_type, ncgf, nsgf, all_potential, tnadd_potential, gth_potential, sgp_potential, upf_potential, se_parameter, dftb_parameter, xtb_parameter, dftb3_param, zeff, elec_conf, mao, lmax_dftb, alpha_core_charge, ccore_charge, core_charge, core_charge_radius, paw_proj_set, paw_atom, hard_radius, hard0_radius, max_rad_local, covalent_radius, vdw_radius, gpw_r3d_rs_type_forced, harmonics, max_iso_not0, max_s_harm, grid_atom, ngrid_ang, ngrid_rad, lmax_rho0, dft_plus_u_atom, l_of_dft_plus_u, n_of_dft_plus_u, u_minus_j, u_of_dft_plus_u, j_of_dft_plus_u, alpha_of_dft_plus_u, beta_of_dft_plus_u, j0_of_dft_plus_u, occupation_of_dft_plus_u, dispersion, bs_occupation, magnetization, no_optimize, addel, laddel, naddel, orbitals, max_scf, eps_scf, smear, u_ramping, u_minus_j_target, eps_u_ramping, init_u_ramping_each_scf, reltmat, ghost, floating, name, element_symbol, pao_basis_size, pao_potentials, pao_descriptors, nelec)
Get attributes of an atomic kind.
subroutine, public get_qs_kind_set(qs_kind_set, all_potential_present, tnadd_potential_present, gth_potential_present, sgp_potential_present, paw_atom_present, dft_plus_u_atom_present, maxcgf, maxsgf, maxco, maxco_proj, maxgtops, maxlgto, maxlprj, maxnset, maxsgf_set, ncgf, npgf, nset, nsgf, nshell, maxpol, maxlppl, maxlppnl, maxppnl, nelectron, maxder, max_ngrid_rad, max_sph_harm, maxg_iso_not0, lmax_rho0, basis_rcut, basis_type, total_zeff_corr)
Get attributes of an atomic kind set.
Provides all information about a quickstep kind.