Line data Source code
1 : !--------------------------------------------------------------------------------------------------!
2 : ! CP2K: A general program to perform molecular dynamics simulations !
3 : ! Copyright 2000-2025 CP2K developers group <https://cp2k.org> !
4 : ! !
5 : ! SPDX-License-Identifier: GPL-2.0-or-later !
6 : !--------------------------------------------------------------------------------------------------!
7 :
8 : ! **************************************************************************************************
9 : !> \par History
10 : !> - 02.2004 flexible normalization of basis sets [jgh]
11 : !> - 07.2014 Add a set of contraction coefficient that only work on active
12 : !> functions
13 : !> \author Matthias Krack (04.07.2000)
14 : ! **************************************************************************************************
15 : MODULE basis_set_types
16 :
17 : USE ai_coulomb, ONLY: coulomb2
18 : USE bibliography, ONLY: VandeVondele2007,&
19 : cite_reference
20 : USE cp_linked_list_input, ONLY: cp_sll_val_next,&
21 : cp_sll_val_type
22 : USE cp_parser_methods, ONLY: parser_get_object,&
23 : parser_search_string
24 : USE cp_parser_types, ONLY: cp_parser_type,&
25 : parser_create,&
26 : parser_release
27 : USE input_section_types, ONLY: section_vals_list_get,&
28 : section_vals_type,&
29 : section_vals_val_get
30 : USE input_val_types, ONLY: val_get,&
31 : val_type
32 : USE kinds, ONLY: default_path_length,&
33 : default_string_length,&
34 : dp
35 : USE mathconstants, ONLY: dfac,&
36 : pi
37 : USE memory_utilities, ONLY: reallocate
38 : USE message_passing, ONLY: mp_para_env_type
39 : USE orbital_pointers, ONLY: coset,&
40 : indco,&
41 : init_orbital_pointers,&
42 : nco,&
43 : ncoset,&
44 : nso,&
45 : nsoset
46 : USE orbital_symbols, ONLY: cgf_symbol,&
47 : sgf_symbol
48 : USE orbital_transformation_matrices, ONLY: init_spherical_harmonics,&
49 : orbtramat
50 : USE sto_ng, ONLY: get_sto_ng
51 : USE string_utilities, ONLY: integer_to_string,&
52 : remove_word,&
53 : uppercase
54 : USE util, ONLY: sort
55 : #include "../base/base_uses.f90"
56 :
57 : IMPLICIT NONE
58 :
59 : PRIVATE
60 :
61 : ! Global parameters (only in this module)
62 :
63 : CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'basis_set_types'
64 :
65 : ! basis set sort criteria
66 : INTEGER, PARAMETER, PUBLIC :: basis_sort_default = 0, &
67 : basis_sort_zet = 1
68 :
69 : ! **************************************************************************************************
70 : ! Define the Gaussian-type orbital basis set type
71 :
72 : TYPE gto_basis_set_type
73 : !MK PRIVATE
74 : CHARACTER(LEN=default_string_length) :: name = ""
75 : CHARACTER(LEN=default_string_length) :: aliases = ""
76 : REAL(KIND=dp) :: kind_radius = 0.0_dp
77 : REAL(KIND=dp) :: short_kind_radius = 0.0_dp
78 : INTEGER :: norm_type = -1
79 : INTEGER :: ncgf = -1, nset = -1, nsgf = -1
80 : CHARACTER(LEN=12), DIMENSION(:), POINTER :: cgf_symbol => NULL()
81 : CHARACTER(LEN=6), DIMENSION(:), POINTER :: sgf_symbol => NULL()
82 : REAL(KIND=dp), DIMENSION(:), POINTER :: norm_cgf => NULL(), set_radius => NULL()
83 : INTEGER, DIMENSION(:), POINTER :: lmax => NULL(), lmin => NULL(), &
84 : lx => NULL(), ly => NULL(), lz => NULL(), &
85 : m => NULL(), ncgf_set => NULL(), &
86 : npgf => NULL(), nsgf_set => NULL(), nshell => NULL()
87 : REAL(KIND=dp), DIMENSION(:, :), POINTER :: cphi => NULL(), pgf_radius => NULL(), sphi => NULL(), &
88 : scon => NULL(), zet => NULL()
89 : INTEGER, DIMENSION(:, :), POINTER :: first_cgf => NULL(), first_sgf => NULL(), l => NULL(), &
90 : last_cgf => NULL(), last_sgf => NULL(), n => NULL()
91 : REAL(KIND=dp), DIMENSION(:, :, :), POINTER :: gcc => NULL()
92 : END TYPE gto_basis_set_type
93 :
94 : TYPE gto_basis_set_p_type
95 : TYPE(gto_basis_set_type), POINTER :: gto_basis_set => NULL()
96 : END TYPE gto_basis_set_p_type
97 :
98 : ! **************************************************************************************************
99 : ! Define the Slater-type orbital basis set type
100 :
101 : TYPE sto_basis_set_type
102 : PRIVATE
103 : CHARACTER(LEN=default_string_length) :: name = ""
104 : INTEGER :: nshell = -1
105 : CHARACTER(LEN=6), DIMENSION(:), POINTER :: symbol => NULL()
106 : INTEGER, DIMENSION(:), POINTER :: nq => NULL(), lq => NULL()
107 : REAL(KIND=dp), DIMENSION(:), POINTER :: zet => NULL()
108 : END TYPE sto_basis_set_type
109 :
110 : ! **************************************************************************************************
111 : INTERFACE read_gto_basis_set
112 : MODULE PROCEDURE read_gto_basis_set1, read_gto_basis_set2
113 : END INTERFACE
114 : ! **************************************************************************************************
115 :
116 : ! Public subroutines
117 : PUBLIC :: allocate_gto_basis_set, &
118 : deallocate_gto_basis_set, &
119 : get_gto_basis_set, &
120 : init_aux_basis_set, &
121 : init_cphi_and_sphi, &
122 : init_orb_basis_set, &
123 : read_gto_basis_set, &
124 : copy_gto_basis_set, &
125 : create_primitive_basis_set, &
126 : combine_basis_sets, &
127 : set_gto_basis_set, &
128 : sort_gto_basis_set, &
129 : write_gto_basis_set, &
130 : write_orb_basis_set
131 :
132 : PUBLIC :: allocate_sto_basis_set, &
133 : read_sto_basis_set, &
134 : create_gto_from_sto_basis, &
135 : deallocate_sto_basis_set, &
136 : set_sto_basis_set, &
137 : srules, process_gto_basis
138 :
139 : ! Public data types
140 : PUBLIC :: gto_basis_set_p_type, &
141 : gto_basis_set_type, &
142 : sto_basis_set_type
143 :
144 : CONTAINS
145 :
146 : ! **************************************************************************************************
147 : !> \brief ...
148 : !> \param gto_basis_set ...
149 : ! **************************************************************************************************
150 21158 : SUBROUTINE allocate_gto_basis_set(gto_basis_set)
151 :
152 : ! Allocate a Gaussian-type orbital (GTO) basis set data set.
153 :
154 : ! - Creation (26.10.2000,MK)
155 :
156 : TYPE(gto_basis_set_type), POINTER :: gto_basis_set
157 :
158 21158 : CALL deallocate_gto_basis_set(gto_basis_set)
159 :
160 21158 : ALLOCATE (gto_basis_set)
161 :
162 21158 : END SUBROUTINE allocate_gto_basis_set
163 :
164 : ! **************************************************************************************************
165 : !> \brief ...
166 : !> \param gto_basis_set ...
167 : ! **************************************************************************************************
168 42624 : SUBROUTINE deallocate_gto_basis_set(gto_basis_set)
169 :
170 : ! Deallocate a Gaussian-type orbital (GTO) basis set data set.
171 :
172 : ! - Creation (03.11.2000,MK)
173 :
174 : TYPE(gto_basis_set_type), POINTER :: gto_basis_set
175 :
176 42624 : IF (ASSOCIATED(gto_basis_set)) THEN
177 21466 : IF (ASSOCIATED(gto_basis_set%cgf_symbol)) DEALLOCATE (gto_basis_set%cgf_symbol)
178 21466 : IF (ASSOCIATED(gto_basis_set%sgf_symbol)) DEALLOCATE (gto_basis_set%sgf_symbol)
179 21466 : IF (ASSOCIATED(gto_basis_set%norm_cgf)) DEALLOCATE (gto_basis_set%norm_cgf)
180 21466 : IF (ASSOCIATED(gto_basis_set%set_radius)) DEALLOCATE (gto_basis_set%set_radius)
181 21466 : IF (ASSOCIATED(gto_basis_set%lmax)) DEALLOCATE (gto_basis_set%lmax)
182 21466 : IF (ASSOCIATED(gto_basis_set%lmin)) DEALLOCATE (gto_basis_set%lmin)
183 21466 : IF (ASSOCIATED(gto_basis_set%lx)) DEALLOCATE (gto_basis_set%lx)
184 21466 : IF (ASSOCIATED(gto_basis_set%ly)) DEALLOCATE (gto_basis_set%ly)
185 21466 : IF (ASSOCIATED(gto_basis_set%lz)) DEALLOCATE (gto_basis_set%lz)
186 21466 : IF (ASSOCIATED(gto_basis_set%m)) DEALLOCATE (gto_basis_set%m)
187 21466 : IF (ASSOCIATED(gto_basis_set%ncgf_set)) DEALLOCATE (gto_basis_set%ncgf_set)
188 21466 : IF (ASSOCIATED(gto_basis_set%npgf)) DEALLOCATE (gto_basis_set%npgf)
189 21466 : IF (ASSOCIATED(gto_basis_set%nsgf_set)) DEALLOCATE (gto_basis_set%nsgf_set)
190 21466 : IF (ASSOCIATED(gto_basis_set%nshell)) DEALLOCATE (gto_basis_set%nshell)
191 21466 : IF (ASSOCIATED(gto_basis_set%cphi)) DEALLOCATE (gto_basis_set%cphi)
192 21466 : IF (ASSOCIATED(gto_basis_set%pgf_radius)) DEALLOCATE (gto_basis_set%pgf_radius)
193 21466 : IF (ASSOCIATED(gto_basis_set%sphi)) DEALLOCATE (gto_basis_set%sphi)
194 21466 : IF (ASSOCIATED(gto_basis_set%scon)) DEALLOCATE (gto_basis_set%scon)
195 21466 : IF (ASSOCIATED(gto_basis_set%zet)) DEALLOCATE (gto_basis_set%zet)
196 21466 : IF (ASSOCIATED(gto_basis_set%first_cgf)) DEALLOCATE (gto_basis_set%first_cgf)
197 21466 : IF (ASSOCIATED(gto_basis_set%first_sgf)) DEALLOCATE (gto_basis_set%first_sgf)
198 21466 : IF (ASSOCIATED(gto_basis_set%l)) DEALLOCATE (gto_basis_set%l)
199 21466 : IF (ASSOCIATED(gto_basis_set%last_cgf)) DEALLOCATE (gto_basis_set%last_cgf)
200 21466 : IF (ASSOCIATED(gto_basis_set%last_sgf)) DEALLOCATE (gto_basis_set%last_sgf)
201 21466 : IF (ASSOCIATED(gto_basis_set%n)) DEALLOCATE (gto_basis_set%n)
202 21466 : IF (ASSOCIATED(gto_basis_set%gcc)) DEALLOCATE (gto_basis_set%gcc)
203 21466 : DEALLOCATE (gto_basis_set)
204 : END IF
205 42624 : END SUBROUTINE deallocate_gto_basis_set
206 :
207 : ! **************************************************************************************************
208 : !> \brief ...
209 : !> \param basis_set_in ...
210 : !> \param basis_set_out ...
211 : ! **************************************************************************************************
212 2984 : SUBROUTINE copy_gto_basis_set(basis_set_in, basis_set_out)
213 :
214 : ! Copy a Gaussian-type orbital (GTO) basis set data set.
215 :
216 : TYPE(gto_basis_set_type), INTENT(IN) :: basis_set_in
217 : TYPE(gto_basis_set_type), POINTER :: basis_set_out
218 :
219 : INTEGER :: maxco, maxpgf, maxshell, ncgf, nset, nsgf
220 :
221 2984 : CALL allocate_gto_basis_set(basis_set_out)
222 :
223 2984 : basis_set_out%name = basis_set_in%name
224 2984 : basis_set_out%aliases = basis_set_in%aliases
225 2984 : basis_set_out%kind_radius = basis_set_in%kind_radius
226 2984 : basis_set_out%norm_type = basis_set_in%norm_type
227 2984 : basis_set_out%nset = basis_set_in%nset
228 2984 : basis_set_out%ncgf = basis_set_in%ncgf
229 2984 : basis_set_out%nsgf = basis_set_in%nsgf
230 2984 : nset = basis_set_in%nset
231 2984 : ncgf = basis_set_in%ncgf
232 2984 : nsgf = basis_set_in%nsgf
233 8952 : ALLOCATE (basis_set_out%cgf_symbol(ncgf))
234 8952 : ALLOCATE (basis_set_out%sgf_symbol(nsgf))
235 47852 : basis_set_out%cgf_symbol = basis_set_in%cgf_symbol
236 42702 : basis_set_out%sgf_symbol = basis_set_in%sgf_symbol
237 8952 : ALLOCATE (basis_set_out%norm_cgf(ncgf))
238 47852 : basis_set_out%norm_cgf = basis_set_in%norm_cgf
239 8952 : ALLOCATE (basis_set_out%set_radius(nset))
240 11926 : basis_set_out%set_radius = basis_set_in%set_radius
241 26856 : ALLOCATE (basis_set_out%lmax(nset), basis_set_out%lmin(nset), basis_set_out%npgf(nset), basis_set_out%nshell(nset))
242 11926 : basis_set_out%lmax = basis_set_in%lmax
243 11926 : basis_set_out%lmin = basis_set_in%lmin
244 11926 : basis_set_out%npgf = basis_set_in%npgf
245 11926 : basis_set_out%nshell = basis_set_in%nshell
246 26856 : ALLOCATE (basis_set_out%lx(ncgf), basis_set_out%ly(ncgf), basis_set_out%lz(ncgf), basis_set_out%m(nsgf))
247 47852 : basis_set_out%lx = basis_set_in%lx
248 47852 : basis_set_out%ly = basis_set_in%ly
249 47852 : basis_set_out%lz = basis_set_in%lz
250 42702 : basis_set_out%m = basis_set_in%m
251 14920 : ALLOCATE (basis_set_out%ncgf_set(nset), basis_set_out%nsgf_set(nset))
252 11926 : basis_set_out%ncgf_set = basis_set_in%ncgf_set
253 11926 : basis_set_out%nsgf_set = basis_set_in%nsgf_set
254 2984 : maxco = SIZE(basis_set_in%cphi, 1)
255 29840 : ALLOCATE (basis_set_out%cphi(maxco, ncgf), basis_set_out%sphi(maxco, nsgf), basis_set_out%scon(maxco, nsgf))
256 1243232 : basis_set_out%cphi = basis_set_in%cphi
257 1046462 : basis_set_out%sphi = basis_set_in%sphi
258 1046462 : basis_set_out%scon = basis_set_in%scon
259 11926 : maxpgf = MAXVAL(basis_set_in%npgf)
260 20888 : ALLOCATE (basis_set_out%pgf_radius(maxpgf, nset), basis_set_out%zet(maxpgf, nset))
261 45042 : basis_set_out%pgf_radius = basis_set_in%pgf_radius
262 45042 : basis_set_out%zet = basis_set_in%zet
263 11926 : maxshell = MAXVAL(basis_set_in%nshell)
264 20888 : ALLOCATE (basis_set_out%first_cgf(maxshell, nset), basis_set_out%first_sgf(maxshell, nset))
265 20888 : ALLOCATE (basis_set_out%last_cgf(maxshell, nset), basis_set_out%last_sgf(maxshell, nset))
266 32332 : basis_set_out%first_cgf = basis_set_in%first_cgf
267 32332 : basis_set_out%first_sgf = basis_set_in%first_sgf
268 32332 : basis_set_out%last_cgf = basis_set_in%last_cgf
269 32332 : basis_set_out%last_sgf = basis_set_in%last_sgf
270 20888 : ALLOCATE (basis_set_out%n(maxshell, nset), basis_set_out%l(maxshell, nset))
271 32332 : basis_set_out%n = basis_set_in%n
272 32332 : basis_set_out%l = basis_set_in%l
273 14920 : ALLOCATE (basis_set_out%gcc(maxpgf, maxshell, nset))
274 115906 : basis_set_out%gcc = basis_set_in%gcc
275 :
276 2984 : END SUBROUTINE copy_gto_basis_set
277 :
278 : ! **************************************************************************************************
279 : !> \brief ...
280 : !> \param basis_set ...
281 : !> \param pbasis ...
282 : !> \param lmax ...
283 : ! **************************************************************************************************
284 70 : SUBROUTINE create_primitive_basis_set(basis_set, pbasis, lmax)
285 :
286 : ! Create a primitives only basis set
287 :
288 : TYPE(gto_basis_set_type), INTENT(IN) :: basis_set
289 : TYPE(gto_basis_set_type), POINTER :: pbasis
290 : INTEGER, INTENT(IN), OPTIONAL :: lmax
291 :
292 : INTEGER :: i, ico, ip, ipgf, iset, ishell, l, lm, &
293 : lshell, m, maxco, mpgf, nc, ncgf, ns, &
294 : nset, nsgf
295 70 : INTEGER, ALLOCATABLE, DIMENSION(:) :: nindex, nprim
296 : REAL(KIND=dp) :: zet0
297 70 : REAL(KIND=dp), ALLOCATABLE, DIMENSION(:, :) :: zet, zeta
298 :
299 240 : mpgf = SUM(basis_set%npgf)
300 240 : lm = MAXVAL(basis_set%lmax)
301 770 : ALLOCATE (zet(mpgf, 0:lm), zeta(mpgf, lm + 1), nindex(mpgf), nprim(0:lm))
302 2320 : zet = 0.0_dp
303 2320 : zeta = 0.0_dp
304 270 : DO l = 0, lm
305 200 : ip = 0
306 734 : DO iset = 1, basis_set%nset
307 734 : IF (basis_set%lmin(iset) <= l .AND. basis_set%lmax(iset) >= l) THEN
308 964 : DO ipgf = 1, basis_set%npgf(iset)
309 764 : ip = ip + 1
310 964 : zet(ip, l) = basis_set%zet(ipgf, iset)
311 : END DO
312 : END IF
313 : END DO
314 270 : nprim(l) = ip
315 : END DO
316 :
317 : ! sort exponents
318 270 : DO l = 0, lm
319 964 : zet(1:nprim(l), l) = -zet(1:nprim(l), l)
320 200 : CALL sort(zet(1:nprim(l), l), nprim(l), nindex)
321 : ! remove duplicates
322 200 : ip = 0
323 200 : zet0 = 0.0_dp
324 964 : DO i = 1, nprim(l)
325 964 : IF (ABS(zet0 - zet(i, l)) > 1.e-6_dp) THEN
326 764 : ip = ip + 1
327 764 : zeta(ip, l + 1) = zet(i, l)
328 : END IF
329 : END DO
330 200 : nprim(l) = ip
331 : !
332 1034 : zeta(1:ip, l + 1) = -zeta(1:ip, l + 1)
333 : END DO
334 :
335 70 : CALL allocate_gto_basis_set(pbasis)
336 :
337 70 : IF (PRESENT(lmax)) THEN
338 : ! if requested, reduce max l val
339 16 : lm = MIN(lm, lmax)
340 : END IF
341 :
342 70 : pbasis%name = basis_set%name//"_primitive"
343 70 : pbasis%kind_radius = basis_set%kind_radius
344 70 : pbasis%short_kind_radius = basis_set%short_kind_radius
345 70 : pbasis%norm_type = basis_set%norm_type
346 70 : nset = lm + 1
347 70 : pbasis%nset = nset
348 420 : ALLOCATE (pbasis%lmax(nset), pbasis%lmin(nset), pbasis%npgf(nset), pbasis%nshell(nset))
349 228 : DO iset = 1, nset
350 158 : pbasis%lmax(iset) = iset - 1
351 158 : pbasis%lmin(iset) = iset - 1
352 158 : pbasis%npgf(iset) = nprim(iset - 1)
353 228 : pbasis%nshell(iset) = nprim(iset - 1)
354 : END DO
355 70 : pbasis%ncgf = 0
356 70 : pbasis%nsgf = 0
357 228 : DO l = 0, lm
358 158 : pbasis%ncgf = pbasis%ncgf + nprim(l)*((l + 1)*(l + 2))/2
359 228 : pbasis%nsgf = pbasis%nsgf + nprim(l)*(2*l + 1)
360 : END DO
361 270 : mpgf = MAXVAL(nprim)
362 280 : ALLOCATE (pbasis%zet(mpgf, nset))
363 1034 : pbasis%zet(1:mpgf, 1:nset) = zeta(1:mpgf, 1:nset)
364 :
365 420 : ALLOCATE (pbasis%l(mpgf, nset), pbasis%n(mpgf, nset))
366 228 : DO iset = 1, nset
367 826 : DO ip = 1, nprim(iset - 1)
368 598 : pbasis%l(ip, iset) = iset - 1
369 756 : pbasis%n(ip, iset) = iset + ip - 1
370 : END DO
371 : END DO
372 :
373 210 : ALLOCATE (pbasis%cgf_symbol(pbasis%ncgf))
374 210 : ALLOCATE (pbasis%lx(pbasis%ncgf))
375 210 : ALLOCATE (pbasis%ly(pbasis%ncgf))
376 210 : ALLOCATE (pbasis%lz(pbasis%ncgf))
377 210 : ALLOCATE (pbasis%m(pbasis%nsgf))
378 210 : ALLOCATE (pbasis%sgf_symbol(pbasis%nsgf))
379 210 : ALLOCATE (pbasis%ncgf_set(nset), pbasis%nsgf_set(nset))
380 :
381 70 : ncgf = 0
382 70 : nsgf = 0
383 228 : DO iset = 1, nset
384 158 : l = iset - 1
385 158 : pbasis%ncgf_set(iset) = nprim(l)*((l + 1)*(l + 2))/2
386 158 : pbasis%nsgf_set(iset) = nprim(l)*(2*l + 1)
387 826 : DO ishell = 1, pbasis%nshell(iset)
388 598 : lshell = pbasis%l(ishell, iset)
389 2670 : DO ico = ncoset(lshell - 1) + 1, ncoset(lshell)
390 2072 : ncgf = ncgf + 1
391 2072 : pbasis%lx(ncgf) = indco(1, ico)
392 2072 : pbasis%ly(ncgf) = indco(2, ico)
393 2072 : pbasis%lz(ncgf) = indco(3, ico)
394 : pbasis%cgf_symbol(ncgf) = &
395 8886 : cgf_symbol(pbasis%n(ishell, iset), (/pbasis%lx(ncgf), pbasis%ly(ncgf), pbasis%lz(ncgf)/))
396 : END DO
397 2510 : DO m = -lshell, lshell
398 1754 : nsgf = nsgf + 1
399 1754 : pbasis%m(nsgf) = m
400 2352 : pbasis%sgf_symbol(nsgf) = sgf_symbol(pbasis%n(ishell, iset), lshell, m)
401 : END DO
402 : END DO
403 : END DO
404 70 : CPASSERT(ncgf == pbasis%ncgf)
405 70 : CPASSERT(nsgf == pbasis%nsgf)
406 :
407 350 : ALLOCATE (pbasis%gcc(mpgf, mpgf, nset))
408 6152 : pbasis%gcc = 0.0_dp
409 228 : DO iset = 1, nset
410 1034 : DO i = 1, mpgf
411 964 : pbasis%gcc(i, i, iset) = 1.0_dp
412 : END DO
413 : END DO
414 :
415 210 : ALLOCATE (pbasis%first_cgf(mpgf, nset))
416 210 : ALLOCATE (pbasis%first_sgf(mpgf, nset))
417 210 : ALLOCATE (pbasis%last_cgf(mpgf, nset))
418 210 : ALLOCATE (pbasis%last_sgf(mpgf, nset))
419 70 : nc = 0
420 70 : ns = 0
421 70 : maxco = 0
422 228 : DO iset = 1, nset
423 756 : DO ishell = 1, pbasis%nshell(iset)
424 598 : lshell = pbasis%l(ishell, iset)
425 598 : pbasis%first_cgf(ishell, iset) = nc + 1
426 598 : nc = nc + nco(lshell)
427 598 : pbasis%last_cgf(ishell, iset) = nc
428 598 : pbasis%first_sgf(ishell, iset) = ns + 1
429 598 : ns = ns + nso(lshell)
430 756 : pbasis%last_sgf(ishell, iset) = ns
431 : END DO
432 228 : maxco = MAX(maxco, pbasis%npgf(iset)*ncoset(pbasis%lmax(iset)))
433 : END DO
434 :
435 210 : ALLOCATE (pbasis%norm_cgf(ncgf))
436 280 : ALLOCATE (pbasis%cphi(maxco, ncgf))
437 172300 : pbasis%cphi = 0.0_dp
438 280 : ALLOCATE (pbasis%sphi(maxco, nsgf))
439 137098 : pbasis%sphi = 0.0_dp
440 210 : ALLOCATE (pbasis%scon(maxco, ncgf))
441 172300 : pbasis%scon = 0.0_dp
442 210 : ALLOCATE (pbasis%set_radius(nset))
443 210 : ALLOCATE (pbasis%pgf_radius(mpgf, nset))
444 1034 : pbasis%pgf_radius = 0.0_dp
445 :
446 70 : CALL init_orb_basis_set(pbasis)
447 :
448 70 : DEALLOCATE (zet, zeta, nindex, nprim)
449 :
450 70 : END SUBROUTINE create_primitive_basis_set
451 :
452 : ! **************************************************************************************************
453 : !> \brief ...
454 : !> \param basis_set ...
455 : !> \param basis_set_add ...
456 : ! **************************************************************************************************
457 42 : SUBROUTINE combine_basis_sets(basis_set, basis_set_add)
458 :
459 : ! Combine two Gaussian-type orbital (GTO) basis sets.
460 :
461 : TYPE(gto_basis_set_type), INTENT(INOUT) :: basis_set
462 : TYPE(gto_basis_set_type), INTENT(IN) :: basis_set_add
463 :
464 42 : CHARACTER(LEN=12), DIMENSION(:), POINTER :: cgf_symbol
465 42 : CHARACTER(LEN=6), DIMENSION(:), POINTER :: sgf_symbol
466 : INTEGER :: iset, ishell, lshell, maxco, maxpgf, &
467 : maxshell, nc, ncgf, ncgfn, ncgfo, ns, &
468 : nset, nsetn, nseto, nsgf, nsgfn, nsgfo
469 :
470 84 : basis_set%name = basis_set%name//basis_set_add%name
471 42 : basis_set%nset = basis_set%nset + basis_set_add%nset
472 42 : basis_set%ncgf = basis_set%ncgf + basis_set_add%ncgf
473 42 : basis_set%nsgf = basis_set%nsgf + basis_set_add%nsgf
474 42 : nset = basis_set%nset
475 42 : ncgf = basis_set%ncgf
476 42 : nsgf = basis_set%nsgf
477 :
478 42 : nsetn = basis_set_add%nset
479 42 : nseto = nset - nsetn
480 42 : CALL reallocate(basis_set%set_radius, 1, nset) ! to be defined later
481 42 : CALL reallocate(basis_set%lmax, 1, nset)
482 42 : CALL reallocate(basis_set%lmin, 1, nset)
483 42 : CALL reallocate(basis_set%npgf, 1, nset)
484 42 : CALL reallocate(basis_set%nshell, 1, nset)
485 160 : basis_set%lmax(nseto + 1:nset) = basis_set_add%lmax(1:nsetn)
486 160 : basis_set%lmin(nseto + 1:nset) = basis_set_add%lmin(1:nsetn)
487 160 : basis_set%npgf(nseto + 1:nset) = basis_set_add%npgf(1:nsetn)
488 160 : basis_set%nshell(nseto + 1:nset) = basis_set_add%nshell(1:nsetn)
489 42 : CALL reallocate(basis_set%ncgf_set, 1, nset)
490 42 : CALL reallocate(basis_set%nsgf_set, 1, nset)
491 160 : basis_set%ncgf_set(nseto + 1:nset) = basis_set_add%ncgf_set(1:nsetn)
492 160 : basis_set%nsgf_set(nseto + 1:nset) = basis_set_add%nsgf_set(1:nsetn)
493 :
494 42 : nsgfn = basis_set_add%nsgf
495 42 : nsgfo = nsgf - nsgfn
496 42 : ncgfn = basis_set_add%ncgf
497 42 : ncgfo = ncgf - ncgfn
498 :
499 210 : ALLOCATE (cgf_symbol(ncgf), sgf_symbol(nsgf))
500 562 : cgf_symbol(1:ncgfo) = basis_set%cgf_symbol(1:ncgfo)
501 1868 : cgf_symbol(ncgfo + 1:ncgf) = basis_set_add%cgf_symbol(1:ncgfn)
502 530 : sgf_symbol(1:nsgfo) = basis_set%sgf_symbol(1:nsgfo)
503 1554 : sgf_symbol(nsgfo + 1:nsgf) = basis_set_add%sgf_symbol(1:nsgfn)
504 42 : DEALLOCATE (basis_set%cgf_symbol, basis_set%sgf_symbol)
505 126 : ALLOCATE (basis_set%cgf_symbol(ncgf), basis_set%sgf_symbol(nsgf))
506 2388 : basis_set%cgf_symbol = cgf_symbol
507 2042 : basis_set%sgf_symbol = sgf_symbol
508 42 : DEALLOCATE (cgf_symbol, sgf_symbol)
509 :
510 42 : CALL reallocate(basis_set%lx, 1, ncgf)
511 42 : CALL reallocate(basis_set%ly, 1, ncgf)
512 42 : CALL reallocate(basis_set%lz, 1, ncgf)
513 42 : CALL reallocate(basis_set%m, 1, nsgf)
514 1868 : basis_set%lx(ncgfo + 1:ncgf) = basis_set_add%lx(1:ncgfn)
515 1868 : basis_set%ly(ncgfo + 1:ncgf) = basis_set_add%ly(1:ncgfn)
516 1868 : basis_set%lz(ncgfo + 1:ncgf) = basis_set_add%lz(1:ncgfn)
517 1554 : basis_set%m(nsgfo + 1:nsgf) = basis_set_add%m(1:nsgfn)
518 :
519 238 : maxpgf = MAXVAL(basis_set%npgf)
520 42 : CALL reallocate(basis_set%zet, 1, maxpgf, 1, nset)
521 42 : nc = SIZE(basis_set_add%zet, 1)
522 160 : DO iset = 1, nsetn
523 770 : basis_set%zet(1:nc, nseto + iset) = basis_set_add%zet(1:nc, iset)
524 : END DO
525 :
526 238 : maxshell = MAXVAL(basis_set%nshell)
527 42 : CALL reallocate(basis_set%l, 1, maxshell, 1, nset)
528 42 : CALL reallocate(basis_set%n, 1, maxshell, 1, nset)
529 42 : nc = SIZE(basis_set_add%l, 1)
530 160 : DO iset = 1, nsetn
531 728 : basis_set%l(1:nc, nseto + iset) = basis_set_add%l(1:nc, iset)
532 770 : basis_set%n(1:nc, nseto + iset) = basis_set_add%n(1:nc, iset)
533 : END DO
534 :
535 42 : CALL reallocate(basis_set%first_cgf, 1, maxshell, 1, nset)
536 42 : CALL reallocate(basis_set%first_sgf, 1, maxshell, 1, nset)
537 42 : CALL reallocate(basis_set%last_cgf, 1, maxshell, 1, nset)
538 42 : CALL reallocate(basis_set%last_sgf, 1, maxshell, 1, nset)
539 42 : nc = 0
540 42 : ns = 0
541 238 : DO iset = 1, nset
542 866 : DO ishell = 1, basis_set%nshell(iset)
543 628 : lshell = basis_set%l(ishell, iset)
544 628 : basis_set%first_cgf(ishell, iset) = nc + 1
545 628 : nc = nc + nco(lshell)
546 628 : basis_set%last_cgf(ishell, iset) = nc
547 628 : basis_set%first_sgf(ishell, iset) = ns + 1
548 628 : ns = ns + nso(lshell)
549 824 : basis_set%last_sgf(ishell, iset) = ns
550 : END DO
551 : END DO
552 :
553 42 : CALL reallocate(basis_set%gcc, 1, maxpgf, 1, maxshell, 1, nset)
554 42 : nc = SIZE(basis_set_add%gcc, 1)
555 42 : ns = SIZE(basis_set_add%gcc, 2)
556 160 : DO iset = 1, nsetn
557 4840 : basis_set%gcc(1:nc, 1:ns, nseto + iset) = basis_set_add%gcc(1:nc, 1:ns, iset)
558 : END DO
559 :
560 : ! these arrays are determined later using initialization calls
561 42 : CALL reallocate(basis_set%norm_cgf, 1, ncgf)
562 42 : maxco = MAX(SIZE(basis_set%cphi, 1), SIZE(basis_set_add%cphi, 1))
563 42 : CALL reallocate(basis_set%cphi, 1, maxco, 1, ncgf)
564 42 : CALL reallocate(basis_set%sphi, 1, maxco, 1, nsgf)
565 42 : CALL reallocate(basis_set%scon, 1, maxco, 1, nsgf)
566 42 : CALL reallocate(basis_set%pgf_radius, 1, maxpgf, 1, nset)
567 :
568 42 : END SUBROUTINE combine_basis_sets
569 :
570 : ! **************************************************************************************************
571 : !> \brief ...
572 : !> \param gto_basis_set ...
573 : !> \param name ...
574 : !> \param aliases ...
575 : !> \param norm_type ...
576 : !> \param kind_radius ...
577 : !> \param ncgf ...
578 : !> \param nset ...
579 : !> \param nsgf ...
580 : !> \param cgf_symbol ...
581 : !> \param sgf_symbol ...
582 : !> \param norm_cgf ...
583 : !> \param set_radius ...
584 : !> \param lmax ...
585 : !> \param lmin ...
586 : !> \param lx ...
587 : !> \param ly ...
588 : !> \param lz ...
589 : !> \param m ...
590 : !> \param ncgf_set ...
591 : !> \param npgf ...
592 : !> \param nsgf_set ...
593 : !> \param nshell ...
594 : !> \param cphi ...
595 : !> \param pgf_radius ...
596 : !> \param sphi ...
597 : !> \param scon ...
598 : !> \param zet ...
599 : !> \param first_cgf ...
600 : !> \param first_sgf ...
601 : !> \param l ...
602 : !> \param last_cgf ...
603 : !> \param last_sgf ...
604 : !> \param n ...
605 : !> \param gcc ...
606 : !> \param maxco ...
607 : !> \param maxl ...
608 : !> \param maxpgf ...
609 : !> \param maxsgf_set ...
610 : !> \param maxshell ...
611 : !> \param maxso ...
612 : !> \param nco_sum ...
613 : !> \param npgf_sum ...
614 : !> \param nshell_sum ...
615 : !> \param maxder ...
616 : !> \param short_kind_radius ...
617 : !> \param npgf_seg_sum number of primitives in "segmented contraction format"
618 : ! **************************************************************************************************
619 32004036 : SUBROUTINE get_gto_basis_set(gto_basis_set, name, aliases, norm_type, kind_radius, ncgf, &
620 : nset, nsgf, cgf_symbol, sgf_symbol, norm_cgf, set_radius, lmax, lmin, lx, ly, lz, &
621 : m, ncgf_set, npgf, nsgf_set, nshell, cphi, pgf_radius, sphi, scon, zet, first_cgf, first_sgf, l, &
622 : last_cgf, last_sgf, n, gcc, maxco, maxl, maxpgf, maxsgf_set, maxshell, maxso, nco_sum, &
623 : npgf_sum, nshell_sum, maxder, short_kind_radius, npgf_seg_sum)
624 :
625 : ! Get informations about a Gaussian-type orbital (GTO) basis set.
626 :
627 : ! - Creation (10.01.2002,MK)
628 :
629 : TYPE(gto_basis_set_type), INTENT(IN) :: gto_basis_set
630 : CHARACTER(LEN=default_string_length), &
631 : INTENT(OUT), OPTIONAL :: name, aliases
632 : INTEGER, INTENT(OUT), OPTIONAL :: norm_type
633 : REAL(KIND=dp), INTENT(OUT), OPTIONAL :: kind_radius
634 : INTEGER, INTENT(OUT), OPTIONAL :: ncgf, nset, nsgf
635 : CHARACTER(LEN=12), DIMENSION(:), OPTIONAL, POINTER :: cgf_symbol
636 : CHARACTER(LEN=6), DIMENSION(:), OPTIONAL, POINTER :: sgf_symbol
637 : REAL(KIND=dp), DIMENSION(:), OPTIONAL, POINTER :: norm_cgf, set_radius
638 : INTEGER, DIMENSION(:), OPTIONAL, POINTER :: lmax, lmin, lx, ly, lz, m, ncgf_set, &
639 : npgf, nsgf_set, nshell
640 : REAL(KIND=dp), DIMENSION(:, :), OPTIONAL, POINTER :: cphi, pgf_radius, sphi, scon, zet
641 : INTEGER, DIMENSION(:, :), OPTIONAL, POINTER :: first_cgf, first_sgf, l, last_cgf, &
642 : last_sgf, n
643 : REAL(KIND=dp), DIMENSION(:, :, :), OPTIONAL, &
644 : POINTER :: gcc
645 : INTEGER, INTENT(OUT), OPTIONAL :: maxco, maxl, maxpgf, maxsgf_set, &
646 : maxshell, maxso, nco_sum, npgf_sum, &
647 : nshell_sum
648 : INTEGER, INTENT(IN), OPTIONAL :: maxder
649 : REAL(KIND=dp), INTENT(OUT), OPTIONAL :: short_kind_radius
650 : INTEGER, INTENT(OUT), OPTIONAL :: npgf_seg_sum
651 :
652 : INTEGER :: iset, nder
653 :
654 32004036 : IF (PRESENT(name)) name = gto_basis_set%name
655 32004036 : IF (PRESENT(aliases)) aliases = gto_basis_set%aliases
656 32004036 : IF (PRESENT(norm_type)) norm_type = gto_basis_set%norm_type
657 32004036 : IF (PRESENT(kind_radius)) kind_radius = gto_basis_set%kind_radius
658 32004036 : IF (PRESENT(short_kind_radius)) short_kind_radius = gto_basis_set%short_kind_radius
659 32004036 : IF (PRESENT(ncgf)) ncgf = gto_basis_set%ncgf
660 32004036 : IF (PRESENT(nset)) nset = gto_basis_set%nset
661 32004036 : IF (PRESENT(nsgf)) nsgf = gto_basis_set%nsgf
662 32004036 : IF (PRESENT(cgf_symbol)) cgf_symbol => gto_basis_set%cgf_symbol
663 32004036 : IF (PRESENT(sgf_symbol)) sgf_symbol => gto_basis_set%sgf_symbol
664 32004036 : IF (PRESENT(norm_cgf)) norm_cgf => gto_basis_set%norm_cgf
665 32004036 : IF (PRESENT(set_radius)) set_radius => gto_basis_set%set_radius
666 32004036 : IF (PRESENT(lmax)) lmax => gto_basis_set%lmax
667 32004036 : IF (PRESENT(lmin)) lmin => gto_basis_set%lmin
668 32004036 : IF (PRESENT(lx)) lx => gto_basis_set%lx
669 32004036 : IF (PRESENT(ly)) ly => gto_basis_set%ly
670 32004036 : IF (PRESENT(lz)) lz => gto_basis_set%lz
671 32004036 : IF (PRESENT(m)) m => gto_basis_set%m
672 32004036 : IF (PRESENT(ncgf_set)) ncgf_set => gto_basis_set%ncgf_set
673 32004036 : IF (PRESENT(npgf)) npgf => gto_basis_set%npgf
674 32004036 : IF (PRESENT(nsgf_set)) nsgf_set => gto_basis_set%nsgf_set
675 32004036 : IF (PRESENT(nshell)) nshell => gto_basis_set%nshell
676 32004036 : IF (PRESENT(cphi)) cphi => gto_basis_set%cphi
677 32004036 : IF (PRESENT(pgf_radius)) pgf_radius => gto_basis_set%pgf_radius
678 32004036 : IF (PRESENT(sphi)) sphi => gto_basis_set%sphi
679 32004036 : IF (PRESENT(scon)) scon => gto_basis_set%scon
680 32004036 : IF (PRESENT(zet)) zet => gto_basis_set%zet
681 32004036 : IF (PRESENT(first_cgf)) first_cgf => gto_basis_set%first_cgf
682 32004036 : IF (PRESENT(first_sgf)) first_sgf => gto_basis_set%first_sgf
683 32004036 : IF (PRESENT(l)) l => gto_basis_set%l
684 32004036 : IF (PRESENT(last_cgf)) last_cgf => gto_basis_set%last_cgf
685 32004036 : IF (PRESENT(last_sgf)) last_sgf => gto_basis_set%last_sgf
686 32004036 : IF (PRESENT(n)) n => gto_basis_set%n
687 32004036 : IF (PRESENT(gcc)) gcc => gto_basis_set%gcc
688 32004036 : IF (PRESENT(maxco)) THEN
689 7167627 : maxco = 0
690 7167627 : IF (PRESENT(maxder)) THEN
691 0 : nder = maxder
692 : ELSE
693 : nder = 0
694 : END IF
695 23514226 : DO iset = 1, gto_basis_set%nset
696 : maxco = MAX(maxco, gto_basis_set%npgf(iset)* &
697 23514226 : ncoset(gto_basis_set%lmax(iset) + nder))
698 : END DO
699 : END IF
700 32004036 : IF (PRESENT(maxl)) THEN
701 6920053 : maxl = -1
702 21934582 : DO iset = 1, gto_basis_set%nset
703 21934582 : maxl = MAX(maxl, gto_basis_set%lmax(iset))
704 : END DO
705 : END IF
706 32004036 : IF (PRESENT(maxpgf)) THEN
707 4248 : maxpgf = 0
708 18708 : DO iset = 1, gto_basis_set%nset
709 18708 : maxpgf = MAX(maxpgf, gto_basis_set%npgf(iset))
710 : END DO
711 : END IF
712 32004036 : IF (PRESENT(maxsgf_set)) THEN
713 443059 : maxsgf_set = 0
714 2231718 : DO iset = 1, gto_basis_set%nset
715 2231718 : maxsgf_set = MAX(maxsgf_set, gto_basis_set%nsgf_set(iset))
716 : END DO
717 : END IF
718 32004036 : IF (PRESENT(maxshell)) THEN ! MAXVAL on structure component avoided
719 2624 : maxshell = 0
720 11958 : DO iset = 1, gto_basis_set%nset
721 11958 : maxshell = MAX(maxshell, gto_basis_set%nshell(iset))
722 : END DO
723 : END IF
724 32004036 : IF (PRESENT(maxso)) THEN
725 1843072 : maxso = 0
726 6773228 : DO iset = 1, gto_basis_set%nset
727 : maxso = MAX(maxso, gto_basis_set%npgf(iset)* &
728 6773228 : nsoset(gto_basis_set%lmax(iset)))
729 : END DO
730 : END IF
731 :
732 32004036 : IF (PRESENT(nco_sum)) THEN
733 72940 : nco_sum = 0
734 257166 : DO iset = 1, gto_basis_set%nset
735 : nco_sum = nco_sum + gto_basis_set%npgf(iset)* &
736 257166 : ncoset(gto_basis_set%lmax(iset))
737 : END DO
738 : END IF
739 32020024 : IF (PRESENT(npgf_sum)) npgf_sum = SUM(gto_basis_set%npgf)
740 32020082 : IF (PRESENT(nshell_sum)) nshell_sum = SUM(gto_basis_set%nshell)
741 32004036 : IF (PRESENT(npgf_seg_sum)) THEN
742 10 : npgf_seg_sum = 0
743 68 : DO iset = 1, gto_basis_set%nset
744 68 : npgf_seg_sum = npgf_seg_sum + gto_basis_set%npgf(iset)*gto_basis_set%nshell(iset)
745 : END DO
746 : END IF
747 :
748 32004036 : END SUBROUTINE get_gto_basis_set
749 :
750 : ! **************************************************************************************************
751 : !> \brief ...
752 : !> \param gto_basis_set ...
753 : ! **************************************************************************************************
754 12 : SUBROUTINE init_aux_basis_set(gto_basis_set)
755 :
756 : ! Initialise a Gaussian-type orbital (GTO) basis set data set.
757 :
758 : ! - Creation (06.12.2000,MK)
759 :
760 : TYPE(gto_basis_set_type), POINTER :: gto_basis_set
761 :
762 : CHARACTER(len=*), PARAMETER :: routineN = 'init_aux_basis_set'
763 :
764 : INTEGER :: handle
765 :
766 : ! -------------------------------------------------------------------------
767 :
768 6 : IF (.NOT. ASSOCIATED(gto_basis_set)) RETURN
769 :
770 6 : CALL timeset(routineN, handle)
771 :
772 12 : SELECT CASE (gto_basis_set%norm_type)
773 : CASE (0)
774 : ! No normalisation requested
775 : CASE (1)
776 6 : CALL init_norm_cgf_aux_2(gto_basis_set)
777 : CASE (2)
778 : ! WARNING this was never tested
779 0 : CALL init_norm_cgf_aux(gto_basis_set)
780 : CASE DEFAULT
781 6 : CPABORT("Normalization method not specified")
782 : END SELECT
783 :
784 : ! Initialise the transformation matrices "pgf" -> "cgf"
785 6 : CALL init_cphi_and_sphi(gto_basis_set)
786 :
787 6 : CALL timestop(handle)
788 :
789 : END SUBROUTINE init_aux_basis_set
790 :
791 : ! **************************************************************************************************
792 : !> \brief ...
793 : !> \param gto_basis_set ...
794 : ! **************************************************************************************************
795 18491 : SUBROUTINE init_cphi_and_sphi(gto_basis_set)
796 :
797 : ! Initialise the matrices for the transformation of primitive Cartesian
798 : ! Gaussian-type functions to contracted Cartesian (cphi) and spherical
799 : ! (sphi) Gaussian-type functions.
800 :
801 : ! - Creation (20.09.2000,MK)
802 :
803 : TYPE(gto_basis_set_type), INTENT(INOUT) :: gto_basis_set
804 :
805 : CHARACTER(len=*), PARAMETER :: routineN = 'init_cphi_and_sphi'
806 :
807 : INTEGER :: first_cgf, first_sgf, handle, icgf, ico, &
808 : ipgf, iset, ishell, l, last_sgf, lmax, &
809 : lmin, n, n1, n2, ncgf, nn, nn1, nn2, &
810 : npgf, nsgf
811 :
812 : ! -------------------------------------------------------------------------
813 : ! Build the Cartesian transformation matrix "cphi"
814 :
815 18491 : CALL timeset(routineN, handle)
816 :
817 6882298 : gto_basis_set%cphi = 0.0_dp
818 68596 : DO iset = 1, gto_basis_set%nset
819 50105 : n = ncoset(gto_basis_set%lmax(iset))
820 148421 : DO ishell = 1, gto_basis_set%nshell(iset)
821 295389 : DO icgf = gto_basis_set%first_cgf(ishell, iset), &
822 129930 : gto_basis_set%last_cgf(ishell, iset)
823 : ico = coset(gto_basis_set%lx(icgf), &
824 : gto_basis_set%ly(icgf), &
825 215564 : gto_basis_set%lz(icgf))
826 988813 : DO ipgf = 1, gto_basis_set%npgf(iset)
827 : gto_basis_set%cphi(ico, icgf) = gto_basis_set%norm_cgf(icgf)* &
828 693424 : gto_basis_set%gcc(ipgf, ishell, iset)
829 908988 : ico = ico + n
830 : END DO
831 : END DO
832 : END DO
833 : END DO
834 :
835 : ! Build the spherical transformation matrix "sphi"
836 :
837 18491 : n = SIZE(gto_basis_set%cphi, 1)
838 :
839 6033843 : gto_basis_set%sphi = 0.0_dp
840 18491 : IF (n .GT. 0) THEN
841 18481 : lmax = -1
842 : ! Ensure proper setup of orbtramat
843 68580 : DO iset = 1, gto_basis_set%nset
844 148399 : DO ishell = 1, gto_basis_set%nshell(iset)
845 129918 : lmax = MAX(lmax, gto_basis_set%l(ishell, iset))
846 : END DO
847 : END DO
848 18481 : CALL init_spherical_harmonics(lmax, -1)
849 :
850 68580 : DO iset = 1, gto_basis_set%nset
851 148399 : DO ishell = 1, gto_basis_set%nshell(iset)
852 79819 : l = gto_basis_set%l(ishell, iset)
853 79819 : first_cgf = gto_basis_set%first_cgf(ishell, iset)
854 79819 : first_sgf = gto_basis_set%first_sgf(ishell, iset)
855 79819 : ncgf = nco(l)
856 79819 : nsgf = nso(l)
857 : CALL dgemm("N", "T", n, nsgf, ncgf, &
858 : 1.0_dp, gto_basis_set%cphi(1, first_cgf), n, &
859 : orbtramat(l)%c2s(1, 1), nsgf, &
860 129918 : 0.0_dp, gto_basis_set%sphi(1, first_sgf), n)
861 : END DO
862 : END DO
863 : END IF
864 :
865 : ! Build the reduced transformation matrix "scon"
866 : ! This matrix transforms from cartesian primitifs to spherical contracted functions
867 : ! "scon" only includes primitifs (lmin -> lmax), whereas sphi is (0 -> lmax)
868 :
869 18491 : n = SIZE(gto_basis_set%scon, 1)
870 :
871 6069045 : gto_basis_set%scon = 0.0_dp
872 18491 : IF (n .GT. 0) THEN
873 68580 : DO iset = 1, gto_basis_set%nset
874 50099 : lmin = gto_basis_set%lmin(iset)
875 50099 : lmax = gto_basis_set%lmax(iset)
876 50099 : npgf = gto_basis_set%npgf(iset)
877 50099 : nn = ncoset(lmax) - ncoset(lmin - 1)
878 148399 : DO ishell = 1, gto_basis_set%nshell(iset)
879 79819 : first_sgf = gto_basis_set%first_sgf(ishell, iset)
880 79819 : last_sgf = gto_basis_set%last_sgf(ishell, iset)
881 421381 : DO ipgf = 1, npgf
882 291463 : nn1 = (ipgf - 1)*ncoset(lmax) + ncoset(lmin - 1) + 1
883 291463 : nn2 = ipgf*ncoset(lmax)
884 291463 : n1 = (ipgf - 1)*nn + 1
885 291463 : n2 = ipgf*nn
886 5040551 : gto_basis_set%scon(n1:n2, first_sgf:last_sgf) = gto_basis_set%sphi(nn1:nn2, first_sgf:last_sgf)
887 : END DO
888 : END DO
889 : END DO
890 : END IF
891 :
892 18491 : CALL timestop(handle)
893 :
894 18491 : END SUBROUTINE init_cphi_and_sphi
895 :
896 : ! **************************************************************************************************
897 : !> \brief ...
898 : !> \param gto_basis_set ...
899 : ! **************************************************************************************************
900 0 : SUBROUTINE init_norm_cgf_aux(gto_basis_set)
901 :
902 : ! Initialise the normalization factors of the contracted Cartesian Gaussian
903 : ! functions, if the Gaussian functions represent charge distributions.
904 :
905 : ! - Creation (07.12.2000,MK)
906 :
907 : TYPE(gto_basis_set_type), INTENT(INOUT) :: gto_basis_set
908 :
909 : INTEGER :: icgf, ico, ipgf, iset, ishell, jco, &
910 : jpgf, ll, lmax, lmin, lx, ly, lz, n, &
911 : npgfa
912 : REAL(KIND=dp) :: fnorm, gcca, gccb
913 0 : REAL(KIND=dp), ALLOCATABLE, DIMENSION(:) :: ff
914 : REAL(KIND=dp), ALLOCATABLE, DIMENSION(:, :) :: gaa
915 : REAL(KIND=dp), ALLOCATABLE, DIMENSION(:, :, :) :: vv
916 0 : REAL(KIND=dp), DIMENSION(:), POINTER :: rpgfa, zeta
917 :
918 : ! -------------------------------------------------------------------------
919 :
920 0 : n = 0
921 0 : ll = 0
922 0 : DO iset = 1, gto_basis_set%nset
923 0 : n = MAX(n, gto_basis_set%npgf(iset)*ncoset(gto_basis_set%lmax(iset)))
924 0 : ll = MAX(ll, gto_basis_set%lmax(iset))
925 : END DO
926 :
927 0 : ALLOCATE (gaa(n, n))
928 0 : ALLOCATE (vv(ncoset(ll), ncoset(ll), ll + ll + 1))
929 0 : ALLOCATE (ff(0:ll + ll))
930 :
931 0 : DO iset = 1, gto_basis_set%nset
932 0 : lmax = gto_basis_set%lmax(iset)
933 0 : lmin = gto_basis_set%lmin(iset)
934 0 : n = ncoset(lmax)
935 0 : npgfa = gto_basis_set%npgf(iset)
936 0 : rpgfa => gto_basis_set%pgf_radius(1:npgfa, iset)
937 0 : zeta => gto_basis_set%zet(1:npgfa, iset)
938 : CALL coulomb2(lmax, npgfa, zeta, rpgfa, lmin, &
939 : lmax, npgfa, zeta, rpgfa, lmin, &
940 0 : (/0.0_dp, 0.0_dp, 0.0_dp/), 0.0_dp, gaa, vv, ff(0:))
941 0 : DO ishell = 1, gto_basis_set%nshell(iset)
942 0 : DO icgf = gto_basis_set%first_cgf(ishell, iset), &
943 0 : gto_basis_set%last_cgf(ishell, iset)
944 0 : lx = gto_basis_set%lx(icgf)
945 0 : ly = gto_basis_set%ly(icgf)
946 0 : lz = gto_basis_set%lz(icgf)
947 0 : ico = coset(lx, ly, lz)
948 0 : fnorm = 0.0_dp
949 0 : DO ipgf = 1, npgfa
950 0 : gcca = gto_basis_set%gcc(ipgf, ishell, iset)
951 0 : jco = coset(lx, ly, lz)
952 0 : DO jpgf = 1, npgfa
953 0 : gccb = gto_basis_set%gcc(jpgf, ishell, iset)
954 0 : fnorm = fnorm + gcca*gccb*gaa(ico, jco)
955 0 : jco = jco + n
956 : END DO
957 0 : ico = ico + n
958 : END DO
959 0 : gto_basis_set%norm_cgf(icgf) = 1.0_dp/SQRT(fnorm)
960 : END DO
961 : END DO
962 : END DO
963 :
964 0 : DEALLOCATE (vv, ff)
965 0 : DEALLOCATE (gaa)
966 :
967 0 : END SUBROUTINE init_norm_cgf_aux
968 :
969 : ! **************************************************************************************************
970 : !> \brief ...
971 : !> \param gto_basis_set ...
972 : ! **************************************************************************************************
973 6 : ELEMENTAL SUBROUTINE init_norm_cgf_aux_2(gto_basis_set)
974 :
975 : ! Initialise the normalization factors of the auxiliary Cartesian Gaussian
976 : ! functions (Kim-Gordon polarization basis) Norm = 1.
977 :
978 : ! - Creation (07.12.2000,GT)
979 :
980 : TYPE(gto_basis_set_type), INTENT(INOUT) :: gto_basis_set
981 :
982 : INTEGER :: icgf, iset, ishell
983 :
984 92 : DO iset = 1, gto_basis_set%nset
985 178 : DO ishell = 1, gto_basis_set%nshell(iset)
986 396 : DO icgf = gto_basis_set%first_cgf(ishell, iset), &
987 172 : gto_basis_set%last_cgf(ishell, iset)
988 396 : gto_basis_set%norm_cgf(icgf) = 1.0_dp
989 : END DO
990 : END DO
991 : END DO
992 :
993 6 : END SUBROUTINE init_norm_cgf_aux_2
994 :
995 : ! **************************************************************************************************
996 : !> \brief Initialise the normalization factors of the contracted Cartesian Gaussian functions.
997 : !> \param gto_basis_set ...
998 : !> \author MK
999 : ! **************************************************************************************************
1000 16861 : ELEMENTAL SUBROUTINE init_norm_cgf_orb(gto_basis_set)
1001 :
1002 : TYPE(gto_basis_set_type), INTENT(INOUT) :: gto_basis_set
1003 :
1004 : INTEGER :: icgf, ipgf, iset, ishell, jpgf, l, lx, &
1005 : ly, lz
1006 : REAL(KIND=dp) :: expzet, fnorm, gcca, gccb, prefac, zeta, &
1007 : zetb
1008 :
1009 61830 : DO iset = 1, gto_basis_set%nset
1010 133551 : DO ishell = 1, gto_basis_set%nshell(iset)
1011 :
1012 71721 : l = gto_basis_set%l(ishell, iset)
1013 :
1014 71721 : expzet = 0.5_dp*REAL(2*l + 3, dp)
1015 :
1016 71721 : fnorm = 0.0_dp
1017 :
1018 347482 : DO ipgf = 1, gto_basis_set%npgf(iset)
1019 275761 : gcca = gto_basis_set%gcc(ipgf, ishell, iset)
1020 275761 : zeta = gto_basis_set%zet(ipgf, iset)
1021 2003157 : DO jpgf = 1, gto_basis_set%npgf(iset)
1022 1655675 : gccb = gto_basis_set%gcc(jpgf, ishell, iset)
1023 1655675 : zetb = gto_basis_set%zet(jpgf, iset)
1024 1931436 : fnorm = fnorm + gcca*gccb/(zeta + zetb)**expzet
1025 : END DO
1026 : END DO
1027 :
1028 71721 : fnorm = 0.5_dp**l*pi**1.5_dp*fnorm
1029 :
1030 268017 : DO icgf = gto_basis_set%first_cgf(ishell, iset), &
1031 116690 : gto_basis_set%last_cgf(ishell, iset)
1032 196296 : lx = gto_basis_set%lx(icgf)
1033 196296 : ly = gto_basis_set%ly(icgf)
1034 196296 : lz = gto_basis_set%lz(icgf)
1035 196296 : prefac = dfac(2*lx - 1)*dfac(2*ly - 1)*dfac(2*lz - 1)
1036 268017 : gto_basis_set%norm_cgf(icgf) = 1.0_dp/SQRT(prefac*fnorm)
1037 : END DO
1038 :
1039 : END DO
1040 : END DO
1041 :
1042 16861 : END SUBROUTINE init_norm_cgf_orb
1043 :
1044 : ! **************************************************************************************************
1045 : !> \brief Initialise the normalization factors of the contracted Cartesian Gaussian
1046 : !> functions used for frozen density representation.
1047 : !> \param gto_basis_set ...
1048 : !> \author GT
1049 : ! **************************************************************************************************
1050 0 : ELEMENTAL SUBROUTINE init_norm_cgf_orb_den(gto_basis_set)
1051 :
1052 : TYPE(gto_basis_set_type), INTENT(INOUT) :: gto_basis_set
1053 :
1054 : INTEGER :: icgf, ipgf, iset, ishell, l
1055 : REAL(KIND=dp) :: expzet, gcca, prefac, zeta
1056 :
1057 0 : DO iset = 1, gto_basis_set%nset
1058 0 : DO ishell = 1, gto_basis_set%nshell(iset)
1059 0 : l = gto_basis_set%l(ishell, iset)
1060 0 : expzet = 0.5_dp*REAL(2*l + 3, dp)
1061 0 : prefac = (1.0_dp/pi)**1.5_dp
1062 0 : DO ipgf = 1, gto_basis_set%npgf(iset)
1063 0 : gcca = gto_basis_set%gcc(ipgf, ishell, iset)
1064 0 : zeta = gto_basis_set%zet(ipgf, iset)
1065 0 : gto_basis_set%gcc(ipgf, ishell, iset) = prefac*zeta**expzet*gcca
1066 : END DO
1067 0 : DO icgf = gto_basis_set%first_cgf(ishell, iset), &
1068 0 : gto_basis_set%last_cgf(ishell, iset)
1069 0 : gto_basis_set%norm_cgf(icgf) = 1.0_dp
1070 : END DO
1071 : END DO
1072 : END DO
1073 :
1074 0 : END SUBROUTINE init_norm_cgf_orb_den
1075 :
1076 : ! **************************************************************************************************
1077 : !> \brief Initialise a Gaussian-type orbital (GTO) basis set data set.
1078 : !> \param gto_basis_set ...
1079 : !> \author MK
1080 : ! **************************************************************************************************
1081 33722 : SUBROUTINE init_orb_basis_set(gto_basis_set)
1082 :
1083 : TYPE(gto_basis_set_type), POINTER :: gto_basis_set
1084 :
1085 : CHARACTER(len=*), PARAMETER :: routineN = 'init_orb_basis_set'
1086 :
1087 : INTEGER :: handle
1088 :
1089 : ! -------------------------------------------------------------------------
1090 :
1091 16861 : IF (.NOT. ASSOCIATED(gto_basis_set)) RETURN
1092 :
1093 16861 : CALL timeset(routineN, handle)
1094 :
1095 16861 : SELECT CASE (gto_basis_set%norm_type)
1096 : CASE (0)
1097 : ! No normalisation requested
1098 : CASE (1)
1099 0 : CALL init_norm_cgf_orb_den(gto_basis_set)
1100 : CASE (2)
1101 : ! Normalise the primitive Gaussian functions
1102 16749 : CALL normalise_gcc_orb(gto_basis_set)
1103 : ! Compute the normalization factors of the contracted Gaussian-type
1104 : ! functions
1105 16749 : CALL init_norm_cgf_orb(gto_basis_set)
1106 : CASE (3)
1107 112 : CALL init_norm_cgf_orb(gto_basis_set)
1108 : CASE DEFAULT
1109 16861 : CPABORT("Normalization method not specified")
1110 : END SELECT
1111 :
1112 : ! Initialise the transformation matrices "pgf" -> "cgf"
1113 :
1114 16861 : CALL init_cphi_and_sphi(gto_basis_set)
1115 :
1116 16861 : CALL timestop(handle)
1117 :
1118 : END SUBROUTINE init_orb_basis_set
1119 :
1120 : ! **************************************************************************************************
1121 : !> \brief Normalise the primitive Cartesian Gaussian functions. The normalization
1122 : !> factor is included in the Gaussian contraction coefficients.
1123 : !> \param gto_basis_set ...
1124 : !> \author MK
1125 : ! **************************************************************************************************
1126 16749 : SUBROUTINE normalise_gcc_orb(gto_basis_set)
1127 :
1128 : TYPE(gto_basis_set_type), POINTER :: gto_basis_set
1129 :
1130 : INTEGER :: ipgf, iset, ishell, l
1131 : REAL(KIND=dp) :: expzet, gcca, prefac, zeta
1132 :
1133 61468 : DO iset = 1, gto_basis_set%nset
1134 132939 : DO ishell = 1, gto_basis_set%nshell(iset)
1135 71471 : l = gto_basis_set%l(ishell, iset)
1136 71471 : expzet = 0.25_dp*REAL(2*l + 3, dp)
1137 71471 : prefac = 2.0_dp**l*(2.0_dp/pi)**0.75_dp
1138 390555 : DO ipgf = 1, gto_basis_set%npgf(iset)
1139 274365 : gcca = gto_basis_set%gcc(ipgf, ishell, iset)
1140 274365 : zeta = gto_basis_set%zet(ipgf, iset)
1141 345836 : gto_basis_set%gcc(ipgf, ishell, iset) = prefac*zeta**expzet*gcca
1142 : END DO
1143 : END DO
1144 : END DO
1145 :
1146 16749 : END SUBROUTINE normalise_gcc_orb
1147 :
1148 : ! **************************************************************************************************
1149 : !> \brief Read a Gaussian-type orbital (GTO) basis set from the database file.
1150 : !> \param element_symbol ...
1151 : !> \param basis_set_name ...
1152 : !> \param gto_basis_set ...
1153 : !> \param para_env ...
1154 : !> \param dft_section ...
1155 : !> \author MK
1156 : ! **************************************************************************************************
1157 11284 : SUBROUTINE read_gto_basis_set1(element_symbol, basis_set_name, gto_basis_set, &
1158 : para_env, dft_section)
1159 :
1160 : CHARACTER(LEN=*), INTENT(IN) :: element_symbol, basis_set_name
1161 : TYPE(gto_basis_set_type), INTENT(INOUT) :: gto_basis_set
1162 : TYPE(mp_para_env_type), POINTER :: para_env
1163 : TYPE(section_vals_type), POINTER :: dft_section
1164 :
1165 : CHARACTER(LEN=240) :: line
1166 : CHARACTER(LEN=242) :: line2
1167 : CHARACTER(len=default_path_length) :: basis_set_file_name, tmp
1168 : CHARACTER(LEN=default_path_length), DIMENSION(:), &
1169 11284 : POINTER :: cbasis
1170 11284 : CHARACTER(LEN=LEN(basis_set_name)) :: bsname
1171 11284 : CHARACTER(LEN=LEN(basis_set_name)+2) :: bsname2
1172 11284 : CHARACTER(LEN=LEN(element_symbol)) :: symbol
1173 11284 : CHARACTER(LEN=LEN(element_symbol)+2) :: symbol2
1174 : INTEGER :: i, ibasis, ico, ipgf, irep, iset, ishell, lshell, m, maxco, maxl, maxpgf, &
1175 : maxshell, nbasis, ncgf, nmin, nset, nsgf, sort_method, strlen1, strlen2
1176 11284 : INTEGER, DIMENSION(:), POINTER :: lmax, lmin, npgf, nshell
1177 11284 : INTEGER, DIMENSION(:, :), POINTER :: l, n
1178 : LOGICAL :: basis_found, found, match
1179 11284 : REAL(KIND=dp), DIMENSION(:, :), POINTER :: zet
1180 11284 : REAL(KIND=dp), DIMENSION(:, :, :), POINTER :: gcc
1181 : TYPE(cp_parser_type) :: parser
1182 :
1183 11284 : line = ""
1184 11284 : line2 = ""
1185 11284 : symbol = ""
1186 11284 : symbol2 = ""
1187 11284 : bsname = ""
1188 11284 : bsname2 = ""
1189 :
1190 11284 : nbasis = 1
1191 :
1192 11284 : gto_basis_set%name = basis_set_name
1193 11284 : gto_basis_set%aliases = basis_set_name
1194 : CALL section_vals_val_get(dft_section, "BASIS_SET_FILE_NAME", &
1195 11284 : n_rep_val=nbasis)
1196 33852 : ALLOCATE (cbasis(nbasis))
1197 24440 : DO ibasis = 1, nbasis
1198 : CALL section_vals_val_get(dft_section, "BASIS_SET_FILE_NAME", &
1199 13156 : i_rep_val=ibasis, c_val=cbasis(ibasis))
1200 13156 : basis_set_file_name = cbasis(ibasis)
1201 13156 : tmp = basis_set_file_name
1202 13156 : CALL uppercase(tmp)
1203 24440 : IF (INDEX(tmp, "MOLOPT") .NE. 0) CALL cite_reference(VandeVondele2007)
1204 : END DO
1205 :
1206 : ! Search for the requested basis set in the basis set file
1207 : ! until the basis set is found or the end of file is reached
1208 :
1209 11284 : basis_found = .FALSE.
1210 23424 : basis_loop: DO ibasis = 1, nbasis
1211 13090 : IF (basis_found) EXIT basis_loop
1212 12140 : basis_set_file_name = cbasis(ibasis)
1213 12140 : CALL parser_create(parser, basis_set_file_name, para_env=para_env)
1214 :
1215 12140 : bsname = basis_set_name
1216 12140 : symbol = element_symbol
1217 12140 : irep = 0
1218 :
1219 12140 : tmp = basis_set_name
1220 12140 : CALL uppercase(tmp)
1221 12140 : IF (INDEX(tmp, "MOLOPT") .NE. 0) CALL cite_reference(VandeVondele2007)
1222 :
1223 12140 : nset = 0
1224 12140 : maxshell = 0
1225 12140 : maxpgf = 0
1226 12140 : maxco = 0
1227 12140 : ncgf = 0
1228 12140 : nsgf = 0
1229 12140 : gto_basis_set%nset = nset
1230 12140 : gto_basis_set%ncgf = ncgf
1231 12140 : gto_basis_set%nsgf = nsgf
1232 12140 : CALL reallocate(gto_basis_set%lmax, 1, nset)
1233 12140 : CALL reallocate(gto_basis_set%lmin, 1, nset)
1234 12140 : CALL reallocate(gto_basis_set%npgf, 1, nset)
1235 12140 : CALL reallocate(gto_basis_set%nshell, 1, nset)
1236 12140 : CALL reallocate(gto_basis_set%n, 1, maxshell, 1, nset)
1237 12140 : CALL reallocate(gto_basis_set%l, 1, maxshell, 1, nset)
1238 12140 : CALL reallocate(gto_basis_set%zet, 1, maxpgf, 1, nset)
1239 12140 : CALL reallocate(gto_basis_set%gcc, 1, maxpgf, 1, maxshell, 1, nset)
1240 12140 : CALL reallocate(gto_basis_set%set_radius, 1, nset)
1241 12140 : CALL reallocate(gto_basis_set%pgf_radius, 1, maxpgf, 1, nset)
1242 12140 : CALL reallocate(gto_basis_set%first_cgf, 1, maxshell, 1, nset)
1243 12140 : CALL reallocate(gto_basis_set%first_sgf, 1, maxshell, 1, nset)
1244 12140 : CALL reallocate(gto_basis_set%last_cgf, 1, maxshell, 1, nset)
1245 12140 : CALL reallocate(gto_basis_set%last_sgf, 1, maxshell, 1, nset)
1246 12140 : CALL reallocate(gto_basis_set%ncgf_set, 1, nset)
1247 12140 : CALL reallocate(gto_basis_set%nsgf_set, 1, nset)
1248 12140 : CALL reallocate(gto_basis_set%cphi, 1, maxco, 1, ncgf)
1249 12140 : CALL reallocate(gto_basis_set%sphi, 1, maxco, 1, nsgf)
1250 12140 : CALL reallocate(gto_basis_set%scon, 1, maxco, 1, nsgf)
1251 12140 : CALL reallocate(gto_basis_set%lx, 1, ncgf)
1252 12140 : CALL reallocate(gto_basis_set%ly, 1, ncgf)
1253 12140 : CALL reallocate(gto_basis_set%lz, 1, ncgf)
1254 12140 : CALL reallocate(gto_basis_set%m, 1, nsgf)
1255 12140 : CALL reallocate(gto_basis_set%norm_cgf, 1, ncgf)
1256 :
1257 12140 : IF (tmp .NE. "NONE") THEN
1258 : search_loop: DO
1259 :
1260 60129 : CALL parser_search_string(parser, TRIM(bsname), .TRUE., found, line)
1261 60129 : IF (found) THEN
1262 59273 : CALL uppercase(symbol)
1263 59273 : CALL uppercase(bsname)
1264 :
1265 59273 : match = .FALSE.
1266 59273 : CALL uppercase(line)
1267 : ! Check both the element symbol and the basis set name
1268 59273 : line2 = " "//line//" "
1269 59273 : symbol2 = " "//TRIM(symbol)//" "
1270 59273 : bsname2 = " "//TRIM(bsname)//" "
1271 59273 : strlen1 = LEN_TRIM(symbol2) + 1
1272 59273 : strlen2 = LEN_TRIM(bsname2) + 1
1273 :
1274 59273 : IF ((INDEX(line2, symbol2(:strlen1)) > 0) .AND. &
1275 11276 : (INDEX(line2, bsname2(:strlen2)) > 0)) match = .TRUE.
1276 59273 : IF (match) THEN
1277 : ! copy all names into aliases field
1278 11276 : i = INDEX(line2, symbol2(:strlen1))
1279 11276 : i = i + 1 + INDEX(line2(i + 1:), " ")
1280 11276 : gto_basis_set%aliases = line2(i:)
1281 :
1282 11276 : NULLIFY (gcc, l, lmax, lmin, n, npgf, nshell, zet)
1283 : ! Read the basis set information
1284 11276 : CALL parser_get_object(parser, nset, newline=.TRUE.)
1285 :
1286 11276 : CALL reallocate(npgf, 1, nset)
1287 11276 : CALL reallocate(nshell, 1, nset)
1288 11276 : CALL reallocate(lmax, 1, nset)
1289 11276 : CALL reallocate(lmin, 1, nset)
1290 11276 : CALL reallocate(n, 1, 1, 1, nset)
1291 :
1292 11276 : maxl = 0
1293 11276 : maxpgf = 0
1294 11276 : maxshell = 0
1295 :
1296 44248 : DO iset = 1, nset
1297 32972 : CALL parser_get_object(parser, n(1, iset), newline=.TRUE.)
1298 32972 : CALL parser_get_object(parser, lmin(iset))
1299 32972 : CALL parser_get_object(parser, lmax(iset))
1300 32972 : CALL parser_get_object(parser, npgf(iset))
1301 32972 : maxl = MAX(maxl, lmax(iset))
1302 32972 : IF (npgf(iset) > maxpgf) THEN
1303 11460 : maxpgf = npgf(iset)
1304 11460 : CALL reallocate(zet, 1, maxpgf, 1, nset)
1305 11460 : CALL reallocate(gcc, 1, maxpgf, 1, maxshell, 1, nset)
1306 : END IF
1307 32972 : nshell(iset) = 0
1308 74890 : DO lshell = lmin(iset), lmax(iset)
1309 41918 : nmin = n(1, iset) + lshell - lmin(iset)
1310 41918 : CALL parser_get_object(parser, ishell)
1311 41918 : nshell(iset) = nshell(iset) + ishell
1312 41918 : IF (nshell(iset) > maxshell) THEN
1313 17604 : maxshell = nshell(iset)
1314 17604 : CALL reallocate(n, 1, maxshell, 1, nset)
1315 17604 : CALL reallocate(l, 1, maxshell, 1, nset)
1316 17604 : CALL reallocate(gcc, 1, maxpgf, 1, maxshell, 1, nset)
1317 : END IF
1318 169249 : DO i = 1, ishell
1319 52441 : n(nshell(iset) - ishell + i, iset) = nmin + i - 1
1320 94359 : l(nshell(iset) - ishell + i, iset) = lshell
1321 : END DO
1322 : END DO
1323 113026 : DO ipgf = 1, npgf(iset)
1324 68778 : CALL parser_get_object(parser, zet(ipgf, iset), newline=.TRUE.)
1325 250721 : DO ishell = 1, nshell(iset)
1326 217749 : CALL parser_get_object(parser, gcc(ipgf, ishell, iset))
1327 : END DO
1328 : END DO
1329 : END DO
1330 :
1331 : ! Maximum angular momentum quantum number of the atomic kind
1332 :
1333 11276 : CALL init_orbital_pointers(maxl)
1334 :
1335 : ! Allocate the global variables
1336 :
1337 11276 : gto_basis_set%nset = nset
1338 11276 : CALL reallocate(gto_basis_set%lmax, 1, nset)
1339 11276 : CALL reallocate(gto_basis_set%lmin, 1, nset)
1340 11276 : CALL reallocate(gto_basis_set%npgf, 1, nset)
1341 11276 : CALL reallocate(gto_basis_set%nshell, 1, nset)
1342 11276 : CALL reallocate(gto_basis_set%n, 1, maxshell, 1, nset)
1343 11276 : CALL reallocate(gto_basis_set%l, 1, maxshell, 1, nset)
1344 11276 : CALL reallocate(gto_basis_set%zet, 1, maxpgf, 1, nset)
1345 11276 : CALL reallocate(gto_basis_set%gcc, 1, maxpgf, 1, maxshell, 1, nset)
1346 :
1347 : ! Copy the basis set information into the data structure
1348 :
1349 44248 : DO iset = 1, nset
1350 32972 : gto_basis_set%lmax(iset) = lmax(iset)
1351 32972 : gto_basis_set%lmin(iset) = lmin(iset)
1352 32972 : gto_basis_set%npgf(iset) = npgf(iset)
1353 32972 : gto_basis_set%nshell(iset) = nshell(iset)
1354 85413 : DO ishell = 1, nshell(iset)
1355 52441 : gto_basis_set%n(ishell, iset) = n(ishell, iset)
1356 52441 : gto_basis_set%l(ishell, iset) = l(ishell, iset)
1357 234384 : DO ipgf = 1, npgf(iset)
1358 201412 : gto_basis_set%gcc(ipgf, ishell, iset) = gcc(ipgf, ishell, iset)
1359 : END DO
1360 : END DO
1361 113026 : DO ipgf = 1, npgf(iset)
1362 101750 : gto_basis_set%zet(ipgf, iset) = zet(ipgf, iset)
1363 : END DO
1364 : END DO
1365 :
1366 : ! Initialise the depending atomic kind information
1367 :
1368 11276 : CALL reallocate(gto_basis_set%set_radius, 1, nset)
1369 11276 : CALL reallocate(gto_basis_set%pgf_radius, 1, maxpgf, 1, nset)
1370 11276 : CALL reallocate(gto_basis_set%first_cgf, 1, maxshell, 1, nset)
1371 11276 : CALL reallocate(gto_basis_set%first_sgf, 1, maxshell, 1, nset)
1372 11276 : CALL reallocate(gto_basis_set%last_cgf, 1, maxshell, 1, nset)
1373 11276 : CALL reallocate(gto_basis_set%last_sgf, 1, maxshell, 1, nset)
1374 11276 : CALL reallocate(gto_basis_set%ncgf_set, 1, nset)
1375 11276 : CALL reallocate(gto_basis_set%nsgf_set, 1, nset)
1376 :
1377 11276 : maxco = 0
1378 11276 : ncgf = 0
1379 11276 : nsgf = 0
1380 :
1381 44248 : DO iset = 1, nset
1382 32972 : gto_basis_set%ncgf_set(iset) = 0
1383 32972 : gto_basis_set%nsgf_set(iset) = 0
1384 85413 : DO ishell = 1, nshell(iset)
1385 52441 : lshell = gto_basis_set%l(ishell, iset)
1386 52441 : gto_basis_set%first_cgf(ishell, iset) = ncgf + 1
1387 52441 : ncgf = ncgf + nco(lshell)
1388 52441 : gto_basis_set%last_cgf(ishell, iset) = ncgf
1389 : gto_basis_set%ncgf_set(iset) = &
1390 52441 : gto_basis_set%ncgf_set(iset) + nco(lshell)
1391 52441 : gto_basis_set%first_sgf(ishell, iset) = nsgf + 1
1392 52441 : nsgf = nsgf + nso(lshell)
1393 52441 : gto_basis_set%last_sgf(ishell, iset) = nsgf
1394 : gto_basis_set%nsgf_set(iset) = &
1395 85413 : gto_basis_set%nsgf_set(iset) + nso(lshell)
1396 : END DO
1397 44248 : maxco = MAX(maxco, npgf(iset)*ncoset(lmax(iset)))
1398 : END DO
1399 :
1400 11276 : gto_basis_set%ncgf = ncgf
1401 11276 : gto_basis_set%nsgf = nsgf
1402 :
1403 11276 : CALL reallocate(gto_basis_set%cphi, 1, maxco, 1, ncgf)
1404 11276 : CALL reallocate(gto_basis_set%sphi, 1, maxco, 1, nsgf)
1405 11276 : CALL reallocate(gto_basis_set%scon, 1, maxco, 1, nsgf)
1406 11276 : CALL reallocate(gto_basis_set%lx, 1, ncgf)
1407 11276 : CALL reallocate(gto_basis_set%ly, 1, ncgf)
1408 11276 : CALL reallocate(gto_basis_set%lz, 1, ncgf)
1409 11276 : CALL reallocate(gto_basis_set%m, 1, nsgf)
1410 11276 : CALL reallocate(gto_basis_set%norm_cgf, 1, ncgf)
1411 33828 : ALLOCATE (gto_basis_set%cgf_symbol(ncgf))
1412 :
1413 33828 : ALLOCATE (gto_basis_set%sgf_symbol(nsgf))
1414 :
1415 11276 : ncgf = 0
1416 11276 : nsgf = 0
1417 :
1418 44248 : DO iset = 1, nset
1419 96689 : DO ishell = 1, nshell(iset)
1420 52441 : lshell = gto_basis_set%l(ishell, iset)
1421 195925 : DO ico = ncoset(lshell - 1) + 1, ncoset(lshell)
1422 143484 : ncgf = ncgf + 1
1423 143484 : gto_basis_set%lx(ncgf) = indco(1, ico)
1424 143484 : gto_basis_set%ly(ncgf) = indco(2, ico)
1425 143484 : gto_basis_set%lz(ncgf) = indco(3, ico)
1426 : gto_basis_set%cgf_symbol(ncgf) = &
1427 : cgf_symbol(n(ishell, iset), (/gto_basis_set%lx(ncgf), &
1428 : gto_basis_set%ly(ncgf), &
1429 626377 : gto_basis_set%lz(ncgf)/))
1430 : END DO
1431 214992 : DO m = -lshell, lshell
1432 129579 : nsgf = nsgf + 1
1433 129579 : gto_basis_set%m(nsgf) = m
1434 : gto_basis_set%sgf_symbol(nsgf) = &
1435 182020 : sgf_symbol(n(ishell, iset), lshell, m)
1436 : END DO
1437 : END DO
1438 : END DO
1439 :
1440 11276 : DEALLOCATE (gcc, l, lmax, lmin, n, npgf, nshell, zet)
1441 :
1442 11276 : basis_found = .TRUE.
1443 11276 : EXIT search_loop
1444 : END IF
1445 : ELSE
1446 : EXIT search_loop
1447 : END IF
1448 : END DO search_loop
1449 : ELSE
1450 8 : match = .FALSE.
1451 16 : ALLOCATE (gto_basis_set%cgf_symbol(ncgf))
1452 16 : ALLOCATE (gto_basis_set%sgf_symbol(nsgf))
1453 : END IF
1454 :
1455 35564 : CALL parser_release(parser)
1456 :
1457 : END DO basis_loop
1458 :
1459 11284 : IF (tmp .NE. "NONE") THEN
1460 11276 : IF (.NOT. basis_found) THEN
1461 0 : basis_set_file_name = ""
1462 0 : DO ibasis = 1, nbasis
1463 0 : basis_set_file_name = TRIM(basis_set_file_name)//"<"//TRIM(cbasis(ibasis))//"> "
1464 : END DO
1465 : CALL cp_abort(__LOCATION__, &
1466 : "The requested basis set <"//TRIM(bsname)// &
1467 : "> for element <"//TRIM(symbol)//"> was not "// &
1468 : "found in the basis set files "// &
1469 0 : TRIM(basis_set_file_name))
1470 : END IF
1471 : END IF
1472 11284 : DEALLOCATE (cbasis)
1473 :
1474 11284 : CALL section_vals_val_get(dft_section, "SORT_BASIS", i_val=sort_method)
1475 11284 : CALL sort_gto_basis_set(gto_basis_set, sort_method)
1476 :
1477 56420 : END SUBROUTINE read_gto_basis_set1
1478 :
1479 : ! **************************************************************************************************
1480 : !> \brief Read a Gaussian-type orbital (GTO) basis set from the database file.
1481 : !> \param element_symbol ...
1482 : !> \param basis_type ...
1483 : !> \param gto_basis_set ...
1484 : !> \param basis_section ...
1485 : !> \param irep ...
1486 : !> \param dft_section ...
1487 : !> \author MK
1488 : ! **************************************************************************************************
1489 174 : SUBROUTINE read_gto_basis_set2(element_symbol, basis_type, gto_basis_set, &
1490 : basis_section, irep, dft_section)
1491 :
1492 : CHARACTER(LEN=*), INTENT(IN) :: element_symbol
1493 : CHARACTER(LEN=*), INTENT(INOUT) :: basis_type
1494 : TYPE(gto_basis_set_type), INTENT(INOUT) :: gto_basis_set
1495 : TYPE(section_vals_type), OPTIONAL, POINTER :: basis_section
1496 : INTEGER :: irep
1497 : TYPE(section_vals_type), OPTIONAL, POINTER :: dft_section
1498 :
1499 : CHARACTER(len=20*default_string_length) :: line_att
1500 : CHARACTER(LEN=240) :: line
1501 : CHARACTER(LEN=242) :: line2
1502 : CHARACTER(LEN=default_path_length) :: bsname, bsname2
1503 174 : CHARACTER(LEN=LEN(element_symbol)) :: symbol
1504 174 : CHARACTER(LEN=LEN(element_symbol)+2) :: symbol2
1505 : INTEGER :: i, ico, ipgf, iset, ishell, lshell, m, &
1506 : maxco, maxl, maxpgf, maxshell, ncgf, &
1507 : nmin, nset, nsgf, sort_method
1508 174 : INTEGER, DIMENSION(:), POINTER :: lmax, lmin, npgf, nshell
1509 174 : INTEGER, DIMENSION(:, :), POINTER :: l, n
1510 : LOGICAL :: is_ok
1511 174 : REAL(KIND=dp), DIMENSION(:, :), POINTER :: zet
1512 174 : REAL(KIND=dp), DIMENSION(:, :, :), POINTER :: gcc
1513 : TYPE(cp_sll_val_type), POINTER :: list
1514 : TYPE(val_type), POINTER :: val
1515 :
1516 174 : line = ""
1517 174 : line2 = ""
1518 174 : symbol = ""
1519 174 : symbol2 = ""
1520 : bsname = ""
1521 174 : bsname2 = ""
1522 :
1523 174 : gto_basis_set%name = " "
1524 174 : gto_basis_set%aliases = " "
1525 :
1526 174 : bsname = " "
1527 174 : symbol = element_symbol
1528 :
1529 174 : nset = 0
1530 174 : maxshell = 0
1531 174 : maxpgf = 0
1532 174 : maxco = 0
1533 174 : ncgf = 0
1534 174 : nsgf = 0
1535 174 : gto_basis_set%nset = nset
1536 174 : gto_basis_set%ncgf = ncgf
1537 174 : gto_basis_set%nsgf = nsgf
1538 174 : CALL reallocate(gto_basis_set%lmax, 1, nset)
1539 174 : CALL reallocate(gto_basis_set%lmin, 1, nset)
1540 174 : CALL reallocate(gto_basis_set%npgf, 1, nset)
1541 174 : CALL reallocate(gto_basis_set%nshell, 1, nset)
1542 174 : CALL reallocate(gto_basis_set%n, 1, maxshell, 1, nset)
1543 174 : CALL reallocate(gto_basis_set%l, 1, maxshell, 1, nset)
1544 174 : CALL reallocate(gto_basis_set%zet, 1, maxpgf, 1, nset)
1545 174 : CALL reallocate(gto_basis_set%gcc, 1, maxpgf, 1, maxshell, 1, nset)
1546 174 : CALL reallocate(gto_basis_set%set_radius, 1, nset)
1547 174 : CALL reallocate(gto_basis_set%pgf_radius, 1, maxpgf, 1, nset)
1548 174 : CALL reallocate(gto_basis_set%first_cgf, 1, maxshell, 1, nset)
1549 174 : CALL reallocate(gto_basis_set%first_sgf, 1, maxshell, 1, nset)
1550 174 : CALL reallocate(gto_basis_set%last_cgf, 1, maxshell, 1, nset)
1551 174 : CALL reallocate(gto_basis_set%last_sgf, 1, maxshell, 1, nset)
1552 174 : CALL reallocate(gto_basis_set%ncgf_set, 1, nset)
1553 174 : CALL reallocate(gto_basis_set%nsgf_set, 1, nset)
1554 174 : CALL reallocate(gto_basis_set%cphi, 1, maxco, 1, ncgf)
1555 174 : CALL reallocate(gto_basis_set%sphi, 1, maxco, 1, nsgf)
1556 174 : CALL reallocate(gto_basis_set%scon, 1, maxco, 1, nsgf)
1557 174 : CALL reallocate(gto_basis_set%lx, 1, ncgf)
1558 174 : CALL reallocate(gto_basis_set%ly, 1, ncgf)
1559 174 : CALL reallocate(gto_basis_set%lz, 1, ncgf)
1560 174 : CALL reallocate(gto_basis_set%m, 1, nsgf)
1561 174 : CALL reallocate(gto_basis_set%norm_cgf, 1, ncgf)
1562 :
1563 174 : basis_type = ""
1564 174 : CALL section_vals_val_get(basis_section, "_SECTION_PARAMETERS_", i_rep_section=irep, c_val=basis_type)
1565 174 : IF (basis_type == "Orbital") basis_type = "ORB"
1566 :
1567 174 : NULLIFY (list, val)
1568 174 : CALL section_vals_list_get(basis_section, "_DEFAULT_KEYWORD_", i_rep_section=irep, list=list)
1569 174 : CALL uppercase(symbol)
1570 174 : CALL uppercase(bsname)
1571 :
1572 174 : NULLIFY (gcc, l, lmax, lmin, n, npgf, nshell, zet)
1573 : ! Read the basis set information
1574 174 : is_ok = cp_sll_val_next(list, val)
1575 174 : IF (.NOT. is_ok) CPABORT("Error reading the Basis set from input file!")
1576 174 : CALL val_get(val, c_val=line_att)
1577 174 : READ (line_att, *) nset
1578 :
1579 174 : CALL reallocate(npgf, 1, nset)
1580 174 : CALL reallocate(nshell, 1, nset)
1581 174 : CALL reallocate(lmax, 1, nset)
1582 174 : CALL reallocate(lmin, 1, nset)
1583 174 : CALL reallocate(n, 1, 1, 1, nset)
1584 :
1585 174 : maxl = 0
1586 174 : maxpgf = 0
1587 174 : maxshell = 0
1588 :
1589 500 : DO iset = 1, nset
1590 326 : is_ok = cp_sll_val_next(list, val)
1591 326 : IF (.NOT. is_ok) CPABORT("Error reading the Basis set from input file!")
1592 326 : CALL val_get(val, c_val=line_att)
1593 326 : READ (line_att, *) n(1, iset)
1594 326 : CALL remove_word(line_att)
1595 326 : READ (line_att, *) lmin(iset)
1596 326 : CALL remove_word(line_att)
1597 326 : READ (line_att, *) lmax(iset)
1598 326 : CALL remove_word(line_att)
1599 326 : READ (line_att, *) npgf(iset)
1600 326 : CALL remove_word(line_att)
1601 326 : maxl = MAX(maxl, lmax(iset))
1602 326 : IF (npgf(iset) > maxpgf) THEN
1603 174 : maxpgf = npgf(iset)
1604 174 : CALL reallocate(zet, 1, maxpgf, 1, nset)
1605 174 : CALL reallocate(gcc, 1, maxpgf, 1, maxshell, 1, nset)
1606 : END IF
1607 326 : nshell(iset) = 0
1608 1082 : DO lshell = lmin(iset), lmax(iset)
1609 756 : nmin = n(1, iset) + lshell - lmin(iset)
1610 756 : READ (line_att, *) ishell
1611 756 : CALL remove_word(line_att)
1612 756 : nshell(iset) = nshell(iset) + ishell
1613 756 : IF (nshell(iset) > maxshell) THEN
1614 334 : maxshell = nshell(iset)
1615 334 : CALL reallocate(n, 1, maxshell, 1, nset)
1616 334 : CALL reallocate(l, 1, maxshell, 1, nset)
1617 334 : CALL reallocate(gcc, 1, maxpgf, 1, maxshell, 1, nset)
1618 : END IF
1619 2126 : DO i = 1, ishell
1620 1044 : n(nshell(iset) - ishell + i, iset) = nmin + i - 1
1621 1800 : l(nshell(iset) - ishell + i, iset) = lshell
1622 : END DO
1623 : END DO
1624 326 : IF (LEN_TRIM(line_att) /= 0) &
1625 0 : CPABORT("Error reading the Basis from input file!")
1626 1308 : DO ipgf = 1, npgf(iset)
1627 808 : is_ok = cp_sll_val_next(list, val)
1628 808 : IF (.NOT. is_ok) CPABORT("Error reading the Basis set from input file!")
1629 808 : CALL val_get(val, c_val=line_att)
1630 1134 : READ (line_att, *) zet(ipgf, iset), (gcc(ipgf, ishell, iset), ishell=1, nshell(iset))
1631 : END DO
1632 : END DO
1633 :
1634 : ! Maximum angular momentum quantum number of the atomic kind
1635 :
1636 174 : CALL init_orbital_pointers(maxl)
1637 :
1638 : ! Allocate the global variables
1639 :
1640 174 : gto_basis_set%nset = nset
1641 174 : CALL reallocate(gto_basis_set%lmax, 1, nset)
1642 174 : CALL reallocate(gto_basis_set%lmin, 1, nset)
1643 174 : CALL reallocate(gto_basis_set%npgf, 1, nset)
1644 174 : CALL reallocate(gto_basis_set%nshell, 1, nset)
1645 174 : CALL reallocate(gto_basis_set%n, 1, maxshell, 1, nset)
1646 174 : CALL reallocate(gto_basis_set%l, 1, maxshell, 1, nset)
1647 174 : CALL reallocate(gto_basis_set%zet, 1, maxpgf, 1, nset)
1648 174 : CALL reallocate(gto_basis_set%gcc, 1, maxpgf, 1, maxshell, 1, nset)
1649 :
1650 : ! Copy the basis set information into the data structure
1651 :
1652 500 : DO iset = 1, nset
1653 326 : gto_basis_set%lmax(iset) = lmax(iset)
1654 326 : gto_basis_set%lmin(iset) = lmin(iset)
1655 326 : gto_basis_set%npgf(iset) = npgf(iset)
1656 326 : gto_basis_set%nshell(iset) = nshell(iset)
1657 1370 : DO ishell = 1, nshell(iset)
1658 1044 : gto_basis_set%n(ishell, iset) = n(ishell, iset)
1659 1044 : gto_basis_set%l(ishell, iset) = l(ishell, iset)
1660 5656 : DO ipgf = 1, npgf(iset)
1661 5330 : gto_basis_set%gcc(ipgf, ishell, iset) = gcc(ipgf, ishell, iset)
1662 : END DO
1663 : END DO
1664 1308 : DO ipgf = 1, npgf(iset)
1665 1134 : gto_basis_set%zet(ipgf, iset) = zet(ipgf, iset)
1666 : END DO
1667 : END DO
1668 :
1669 : ! Initialise the depending atomic kind information
1670 :
1671 174 : CALL reallocate(gto_basis_set%set_radius, 1, nset)
1672 174 : CALL reallocate(gto_basis_set%pgf_radius, 1, maxpgf, 1, nset)
1673 174 : CALL reallocate(gto_basis_set%first_cgf, 1, maxshell, 1, nset)
1674 174 : CALL reallocate(gto_basis_set%first_sgf, 1, maxshell, 1, nset)
1675 174 : CALL reallocate(gto_basis_set%last_cgf, 1, maxshell, 1, nset)
1676 174 : CALL reallocate(gto_basis_set%last_sgf, 1, maxshell, 1, nset)
1677 174 : CALL reallocate(gto_basis_set%ncgf_set, 1, nset)
1678 174 : CALL reallocate(gto_basis_set%nsgf_set, 1, nset)
1679 :
1680 174 : maxco = 0
1681 174 : ncgf = 0
1682 174 : nsgf = 0
1683 :
1684 500 : DO iset = 1, nset
1685 326 : gto_basis_set%ncgf_set(iset) = 0
1686 326 : gto_basis_set%nsgf_set(iset) = 0
1687 1370 : DO ishell = 1, nshell(iset)
1688 1044 : lshell = gto_basis_set%l(ishell, iset)
1689 1044 : gto_basis_set%first_cgf(ishell, iset) = ncgf + 1
1690 1044 : ncgf = ncgf + nco(lshell)
1691 1044 : gto_basis_set%last_cgf(ishell, iset) = ncgf
1692 : gto_basis_set%ncgf_set(iset) = &
1693 1044 : gto_basis_set%ncgf_set(iset) + nco(lshell)
1694 1044 : gto_basis_set%first_sgf(ishell, iset) = nsgf + 1
1695 1044 : nsgf = nsgf + nso(lshell)
1696 1044 : gto_basis_set%last_sgf(ishell, iset) = nsgf
1697 : gto_basis_set%nsgf_set(iset) = &
1698 1370 : gto_basis_set%nsgf_set(iset) + nso(lshell)
1699 : END DO
1700 500 : maxco = MAX(maxco, npgf(iset)*ncoset(lmax(iset)))
1701 : END DO
1702 :
1703 174 : gto_basis_set%ncgf = ncgf
1704 174 : gto_basis_set%nsgf = nsgf
1705 :
1706 174 : CALL reallocate(gto_basis_set%cphi, 1, maxco, 1, ncgf)
1707 174 : CALL reallocate(gto_basis_set%sphi, 1, maxco, 1, nsgf)
1708 174 : CALL reallocate(gto_basis_set%scon, 1, maxco, 1, nsgf)
1709 174 : CALL reallocate(gto_basis_set%lx, 1, ncgf)
1710 174 : CALL reallocate(gto_basis_set%ly, 1, ncgf)
1711 174 : CALL reallocate(gto_basis_set%lz, 1, ncgf)
1712 174 : CALL reallocate(gto_basis_set%m, 1, nsgf)
1713 174 : CALL reallocate(gto_basis_set%norm_cgf, 1, ncgf)
1714 522 : ALLOCATE (gto_basis_set%cgf_symbol(ncgf))
1715 :
1716 522 : ALLOCATE (gto_basis_set%sgf_symbol(nsgf))
1717 :
1718 174 : ncgf = 0
1719 174 : nsgf = 0
1720 :
1721 500 : DO iset = 1, nset
1722 1544 : DO ishell = 1, nshell(iset)
1723 1044 : lshell = gto_basis_set%l(ishell, iset)
1724 5040 : DO ico = ncoset(lshell - 1) + 1, ncoset(lshell)
1725 3996 : ncgf = ncgf + 1
1726 3996 : gto_basis_set%lx(ncgf) = indco(1, ico)
1727 3996 : gto_basis_set%ly(ncgf) = indco(2, ico)
1728 3996 : gto_basis_set%lz(ncgf) = indco(3, ico)
1729 : gto_basis_set%cgf_symbol(ncgf) = &
1730 : cgf_symbol(n(ishell, iset), (/gto_basis_set%lx(ncgf), &
1731 : gto_basis_set%ly(ncgf), &
1732 17028 : gto_basis_set%lz(ncgf)/))
1733 : END DO
1734 4606 : DO m = -lshell, lshell
1735 3236 : nsgf = nsgf + 1
1736 3236 : gto_basis_set%m(nsgf) = m
1737 : gto_basis_set%sgf_symbol(nsgf) = &
1738 4280 : sgf_symbol(n(ishell, iset), lshell, m)
1739 : END DO
1740 : END DO
1741 : END DO
1742 :
1743 174 : DEALLOCATE (gcc, l, lmax, lmin, n, npgf, nshell, zet)
1744 :
1745 174 : IF (PRESENT(dft_section)) THEN
1746 162 : CALL section_vals_val_get(dft_section, "SORT_BASIS", i_val=sort_method)
1747 162 : CALL sort_gto_basis_set(gto_basis_set, sort_method)
1748 : END IF
1749 :
1750 174 : END SUBROUTINE read_gto_basis_set2
1751 :
1752 : ! **************************************************************************************************
1753 : !> \brief Set the components of Gaussian-type orbital (GTO) basis set data set.
1754 : !> \param gto_basis_set ...
1755 : !> \param name ...
1756 : !> \param aliases ...
1757 : !> \param norm_type ...
1758 : !> \param kind_radius ...
1759 : !> \param ncgf ...
1760 : !> \param nset ...
1761 : !> \param nsgf ...
1762 : !> \param cgf_symbol ...
1763 : !> \param sgf_symbol ...
1764 : !> \param norm_cgf ...
1765 : !> \param set_radius ...
1766 : !> \param lmax ...
1767 : !> \param lmin ...
1768 : !> \param lx ...
1769 : !> \param ly ...
1770 : !> \param lz ...
1771 : !> \param m ...
1772 : !> \param ncgf_set ...
1773 : !> \param npgf ...
1774 : !> \param nsgf_set ...
1775 : !> \param nshell ...
1776 : !> \param cphi ...
1777 : !> \param pgf_radius ...
1778 : !> \param sphi ...
1779 : !> \param scon ...
1780 : !> \param zet ...
1781 : !> \param first_cgf ...
1782 : !> \param first_sgf ...
1783 : !> \param l ...
1784 : !> \param last_cgf ...
1785 : !> \param last_sgf ...
1786 : !> \param n ...
1787 : !> \param gcc ...
1788 : !> \param short_kind_radius ...
1789 : !> \author MK
1790 : ! **************************************************************************************************
1791 34538 : SUBROUTINE set_gto_basis_set(gto_basis_set, name, aliases, norm_type, kind_radius, ncgf, &
1792 : nset, nsgf, cgf_symbol, sgf_symbol, norm_cgf, set_radius, lmax, &
1793 : lmin, lx, ly, lz, m, ncgf_set, npgf, nsgf_set, nshell, &
1794 : cphi, pgf_radius, sphi, scon, zet, first_cgf, first_sgf, l, &
1795 : last_cgf, last_sgf, n, gcc, short_kind_radius)
1796 :
1797 : TYPE(gto_basis_set_type), INTENT(INOUT) :: gto_basis_set
1798 : CHARACTER(LEN=default_string_length), INTENT(IN), &
1799 : OPTIONAL :: name, aliases
1800 : INTEGER, INTENT(IN), OPTIONAL :: norm_type
1801 : REAL(KIND=dp), INTENT(IN), OPTIONAL :: kind_radius
1802 : INTEGER, INTENT(IN), OPTIONAL :: ncgf, nset, nsgf
1803 : CHARACTER(LEN=12), DIMENSION(:), OPTIONAL, POINTER :: cgf_symbol
1804 : CHARACTER(LEN=6), DIMENSION(:), OPTIONAL, POINTER :: sgf_symbol
1805 : REAL(KIND=dp), DIMENSION(:), OPTIONAL, POINTER :: norm_cgf, set_radius
1806 : INTEGER, DIMENSION(:), OPTIONAL, POINTER :: lmax, lmin, lx, ly, lz, m, ncgf_set, &
1807 : npgf, nsgf_set, nshell
1808 : REAL(KIND=dp), DIMENSION(:, :), OPTIONAL, POINTER :: cphi, pgf_radius, sphi, scon, zet
1809 : INTEGER, DIMENSION(:, :), OPTIONAL, POINTER :: first_cgf, first_sgf, l, last_cgf, &
1810 : last_sgf, n
1811 : REAL(KIND=dp), DIMENSION(:, :, :), OPTIONAL, &
1812 : POINTER :: gcc
1813 : REAL(KIND=dp), INTENT(IN), OPTIONAL :: short_kind_radius
1814 :
1815 34538 : IF (PRESENT(name)) gto_basis_set%name = name
1816 34538 : IF (PRESENT(aliases)) gto_basis_set%aliases = aliases
1817 34538 : IF (PRESENT(norm_type)) gto_basis_set%norm_type = norm_type
1818 34538 : IF (PRESENT(kind_radius)) gto_basis_set%kind_radius = kind_radius
1819 34538 : IF (PRESENT(short_kind_radius)) gto_basis_set%short_kind_radius = short_kind_radius
1820 34538 : IF (PRESENT(ncgf)) gto_basis_set%ncgf = ncgf
1821 34538 : IF (PRESENT(nset)) gto_basis_set%nset = nset
1822 34538 : IF (PRESENT(nsgf)) gto_basis_set%nsgf = nsgf
1823 34538 : IF (PRESENT(cgf_symbol)) gto_basis_set%cgf_symbol(:) = cgf_symbol(:)
1824 34538 : IF (PRESENT(sgf_symbol)) gto_basis_set%sgf_symbol(:) = sgf_symbol(:)
1825 34538 : IF (PRESENT(norm_cgf)) gto_basis_set%norm_cgf(:) = norm_cgf(:)
1826 154646 : IF (PRESENT(set_radius)) gto_basis_set%set_radius(:) = set_radius(:)
1827 34538 : IF (PRESENT(lmax)) gto_basis_set%lmax(:) = lmax(:)
1828 34538 : IF (PRESENT(lmin)) gto_basis_set%lmin(:) = lmin(:)
1829 34538 : IF (PRESENT(lx)) gto_basis_set%lx(:) = lx(:)
1830 34538 : IF (PRESENT(ly)) gto_basis_set%ly(:) = ly(:)
1831 34538 : IF (PRESENT(lz)) gto_basis_set%lz(:) = lz(:)
1832 34538 : IF (PRESENT(m)) gto_basis_set%m(:) = m(:)
1833 34538 : IF (PRESENT(ncgf_set)) gto_basis_set%ncgf_set(:) = ncgf_set(:)
1834 34538 : IF (PRESENT(npgf)) gto_basis_set%npgf(:) = npgf(:)
1835 34538 : IF (PRESENT(nsgf_set)) gto_basis_set%nsgf_set(:) = nsgf_set(:)
1836 34538 : IF (PRESENT(nshell)) gto_basis_set%nshell(:) = nshell(:)
1837 34538 : IF (PRESENT(cphi)) gto_basis_set%cphi(:, :) = cphi(:, :)
1838 477280 : IF (PRESENT(pgf_radius)) gto_basis_set%pgf_radius(:, :) = pgf_radius(:, :)
1839 34538 : IF (PRESENT(sphi)) gto_basis_set%sphi(:, :) = sphi(:, :)
1840 34538 : IF (PRESENT(scon)) gto_basis_set%scon(:, :) = scon(:, :)
1841 34538 : IF (PRESENT(zet)) gto_basis_set%zet(:, :) = zet(:, :)
1842 34538 : IF (PRESENT(first_cgf)) gto_basis_set%first_cgf(:, :) = first_cgf(:, :)
1843 34538 : IF (PRESENT(first_sgf)) gto_basis_set%first_sgf(:, :) = first_sgf(:, :)
1844 34538 : IF (PRESENT(l)) l(:, :) = gto_basis_set%l(:, :)
1845 34538 : IF (PRESENT(last_cgf)) gto_basis_set%last_cgf(:, :) = last_cgf(:, :)
1846 34538 : IF (PRESENT(last_sgf)) gto_basis_set%last_sgf(:, :) = last_sgf(:, :)
1847 34538 : IF (PRESENT(n)) gto_basis_set%n(:, :) = n(:, :)
1848 34538 : IF (PRESENT(gcc)) gto_basis_set%gcc(:, :, :) = gcc(:, :, :)
1849 :
1850 34538 : END SUBROUTINE set_gto_basis_set
1851 :
1852 : ! **************************************************************************************************
1853 : !> \brief Write a Gaussian-type orbital (GTO) basis set data set to the output unit.
1854 : !> \param gto_basis_set ...
1855 : !> \param output_unit ...
1856 : !> \param header ...
1857 : !> \author MK
1858 : ! **************************************************************************************************
1859 139 : SUBROUTINE write_gto_basis_set(gto_basis_set, output_unit, header)
1860 :
1861 : TYPE(gto_basis_set_type), INTENT(IN) :: gto_basis_set
1862 : INTEGER, INTENT(in) :: output_unit
1863 : CHARACTER(len=*), OPTIONAL :: header
1864 :
1865 : INTEGER :: ipgf, iset, ishell
1866 :
1867 139 : IF (output_unit > 0) THEN
1868 :
1869 139 : IF (PRESENT(header)) THEN
1870 : WRITE (UNIT=output_unit, FMT="(/,T6,A,T41,A40)") &
1871 139 : TRIM(header), TRIM(gto_basis_set%name)
1872 : END IF
1873 :
1874 : WRITE (UNIT=output_unit, FMT="(/,(T8,A,T71,I10))") &
1875 139 : "Number of orbital shell sets: ", &
1876 139 : gto_basis_set%nset, &
1877 139 : "Number of orbital shells: ", &
1878 489 : SUM(gto_basis_set%nshell(:)), &
1879 139 : "Number of primitive Cartesian functions: ", &
1880 489 : SUM(gto_basis_set%npgf(:)), &
1881 139 : "Number of Cartesian basis functions: ", &
1882 139 : gto_basis_set%ncgf, &
1883 139 : "Number of spherical basis functions: ", &
1884 139 : gto_basis_set%nsgf, &
1885 139 : "Norm type: ", &
1886 278 : gto_basis_set%norm_type
1887 :
1888 : WRITE (UNIT=output_unit, FMT="(/,T6,A,T41,A40,/,/,T25,A)") &
1889 139 : "GTO basis set information for", TRIM(gto_basis_set%name), &
1890 278 : "Set Shell n l Exponent Coefficient"
1891 :
1892 489 : DO iset = 1, gto_basis_set%nset
1893 350 : WRITE (UNIT=output_unit, FMT="(A)") ""
1894 976 : DO ishell = 1, gto_basis_set%nshell(iset)
1895 : WRITE (UNIT=output_unit, &
1896 : FMT="(T25,I3,4X,I4,4X,I2,2X,I2,(T51,2F15.6))") &
1897 487 : iset, ishell, &
1898 487 : gto_basis_set%n(ishell, iset), &
1899 487 : gto_basis_set%l(ishell, iset), &
1900 1839 : (gto_basis_set%zet(ipgf, iset), &
1901 2326 : gto_basis_set%gcc(ipgf, ishell, iset), &
1902 3650 : ipgf=1, gto_basis_set%npgf(iset))
1903 : END DO
1904 : END DO
1905 :
1906 : END IF
1907 :
1908 139 : END SUBROUTINE write_gto_basis_set
1909 :
1910 : ! **************************************************************************************************
1911 :
1912 : ! **************************************************************************************************
1913 : !> \brief Write a Gaussian-type orbital (GTO) basis set data set to the output unit.
1914 : !> \param orb_basis_set ...
1915 : !> \param output_unit ...
1916 : !> \param header ...
1917 : !> \author MK
1918 : ! **************************************************************************************************
1919 4550 : SUBROUTINE write_orb_basis_set(orb_basis_set, output_unit, header)
1920 :
1921 : TYPE(gto_basis_set_type), INTENT(IN) :: orb_basis_set
1922 : INTEGER, INTENT(in) :: output_unit
1923 : CHARACTER(len=*), OPTIONAL :: header
1924 :
1925 : INTEGER :: icgf, ico, ipgf, iset, ishell
1926 :
1927 4550 : IF (output_unit > 0) THEN
1928 4550 : IF (PRESENT(header)) THEN
1929 : WRITE (UNIT=output_unit, FMT="(/,T6,A,T41,A40)") &
1930 4550 : TRIM(header), TRIM(orb_basis_set%name)
1931 : END IF
1932 :
1933 : WRITE (UNIT=output_unit, FMT="(/,(T8,A,T71,I10))") &
1934 4550 : "Number of orbital shell sets: ", &
1935 4550 : orb_basis_set%nset, &
1936 4550 : "Number of orbital shells: ", &
1937 17663 : SUM(orb_basis_set%nshell(:)), &
1938 4550 : "Number of primitive Cartesian functions: ", &
1939 17663 : SUM(orb_basis_set%npgf(:)), &
1940 4550 : "Number of Cartesian basis functions: ", &
1941 4550 : orb_basis_set%ncgf, &
1942 4550 : "Number of spherical basis functions: ", &
1943 4550 : orb_basis_set%nsgf, &
1944 4550 : "Norm type: ", &
1945 9100 : orb_basis_set%norm_type
1946 :
1947 : WRITE (UNIT=output_unit, FMT="(/,T8,A,/,/,T25,A)") &
1948 4550 : "Normalised Cartesian orbitals:", &
1949 9100 : "Set Shell Orbital Exponent Coefficient"
1950 :
1951 4550 : icgf = 0
1952 :
1953 17663 : DO iset = 1, orb_basis_set%nset
1954 37067 : DO ishell = 1, orb_basis_set%nshell(iset)
1955 19404 : WRITE (UNIT=output_unit, FMT="(A)") ""
1956 85369 : DO ico = 1, nco(orb_basis_set%l(ishell, iset))
1957 52852 : icgf = icgf + 1
1958 : WRITE (UNIT=output_unit, &
1959 : FMT="(T25,I3,4X,I4,3X,A12,(T51,2F15.6))") &
1960 52852 : iset, ishell, orb_basis_set%cgf_symbol(icgf), &
1961 161945 : (orb_basis_set%zet(ipgf, iset), &
1962 : orb_basis_set%norm_cgf(icgf)* &
1963 214797 : orb_basis_set%gcc(ipgf, ishell, iset), &
1964 339905 : ipgf=1, orb_basis_set%npgf(iset))
1965 : END DO
1966 : END DO
1967 : END DO
1968 : END IF
1969 :
1970 4550 : END SUBROUTINE write_orb_basis_set
1971 :
1972 : ! **************************************************************************************************
1973 : !> \brief ...
1974 : !> \param sto_basis_set ...
1975 : ! **************************************************************************************************
1976 4546 : SUBROUTINE allocate_sto_basis_set(sto_basis_set)
1977 :
1978 : TYPE(sto_basis_set_type), POINTER :: sto_basis_set
1979 :
1980 4546 : CALL deallocate_sto_basis_set(sto_basis_set)
1981 :
1982 4546 : ALLOCATE (sto_basis_set)
1983 :
1984 4546 : END SUBROUTINE allocate_sto_basis_set
1985 :
1986 : ! **************************************************************************************************
1987 : !> \brief ...
1988 : !> \param sto_basis_set ...
1989 : ! **************************************************************************************************
1990 10816 : SUBROUTINE deallocate_sto_basis_set(sto_basis_set)
1991 :
1992 : TYPE(sto_basis_set_type), POINTER :: sto_basis_set
1993 :
1994 : ! -------------------------------------------------------------------------
1995 :
1996 10816 : IF (ASSOCIATED(sto_basis_set)) THEN
1997 4546 : IF (ASSOCIATED(sto_basis_set%symbol)) THEN
1998 4406 : DEALLOCATE (sto_basis_set%symbol)
1999 : END IF
2000 4546 : IF (ASSOCIATED(sto_basis_set%nq)) THEN
2001 4546 : DEALLOCATE (sto_basis_set%nq)
2002 : END IF
2003 4546 : IF (ASSOCIATED(sto_basis_set%lq)) THEN
2004 4546 : DEALLOCATE (sto_basis_set%lq)
2005 : END IF
2006 4546 : IF (ASSOCIATED(sto_basis_set%zet)) THEN
2007 4546 : DEALLOCATE (sto_basis_set%zet)
2008 : END IF
2009 :
2010 4546 : DEALLOCATE (sto_basis_set)
2011 : END IF
2012 10816 : END SUBROUTINE deallocate_sto_basis_set
2013 :
2014 : ! **************************************************************************************************
2015 : !> \brief ...
2016 : !> \param sto_basis_set ...
2017 : !> \param name ...
2018 : !> \param nshell ...
2019 : !> \param symbol ...
2020 : !> \param nq ...
2021 : !> \param lq ...
2022 : !> \param zet ...
2023 : !> \param maxlq ...
2024 : !> \param numsto ...
2025 : ! **************************************************************************************************
2026 4546 : SUBROUTINE get_sto_basis_set(sto_basis_set, name, nshell, symbol, nq, lq, zet, maxlq, numsto)
2027 :
2028 : TYPE(sto_basis_set_type), INTENT(IN) :: sto_basis_set
2029 : CHARACTER(LEN=default_string_length), &
2030 : INTENT(OUT), OPTIONAL :: name
2031 : INTEGER, INTENT(OUT), OPTIONAL :: nshell
2032 : CHARACTER(LEN=6), DIMENSION(:), OPTIONAL, POINTER :: symbol
2033 : INTEGER, DIMENSION(:), OPTIONAL, POINTER :: nq, lq
2034 : REAL(KIND=dp), DIMENSION(:), OPTIONAL, POINTER :: zet
2035 : INTEGER, INTENT(OUT), OPTIONAL :: maxlq, numsto
2036 :
2037 : INTEGER :: iset
2038 :
2039 4546 : IF (PRESENT(name)) name = sto_basis_set%name
2040 4546 : IF (PRESENT(nshell)) nshell = sto_basis_set%nshell
2041 4546 : IF (PRESENT(symbol)) symbol => sto_basis_set%symbol
2042 4546 : IF (PRESENT(nq)) nq => sto_basis_set%nq
2043 4546 : IF (PRESENT(lq)) lq => sto_basis_set%lq
2044 4546 : IF (PRESENT(zet)) zet => sto_basis_set%zet
2045 4546 : IF (PRESENT(maxlq)) THEN
2046 0 : maxlq = MAXVAL(sto_basis_set%lq(1:sto_basis_set%nshell))
2047 : END IF
2048 4546 : IF (PRESENT(numsto)) THEN
2049 0 : numsto = 0
2050 0 : DO iset = 1, sto_basis_set%nshell
2051 0 : numsto = numsto + 2*sto_basis_set%lq(iset) + 1
2052 : END DO
2053 : END IF
2054 :
2055 4546 : END SUBROUTINE get_sto_basis_set
2056 :
2057 : ! **************************************************************************************************
2058 : !> \brief ...
2059 : !> \param sto_basis_set ...
2060 : !> \param name ...
2061 : !> \param nshell ...
2062 : !> \param symbol ...
2063 : !> \param nq ...
2064 : !> \param lq ...
2065 : !> \param zet ...
2066 : ! **************************************************************************************************
2067 4542 : SUBROUTINE set_sto_basis_set(sto_basis_set, name, nshell, symbol, nq, lq, zet)
2068 :
2069 : TYPE(sto_basis_set_type), INTENT(INOUT) :: sto_basis_set
2070 : CHARACTER(LEN=default_string_length), INTENT(IN), &
2071 : OPTIONAL :: name
2072 : INTEGER, INTENT(IN), OPTIONAL :: nshell
2073 : CHARACTER(LEN=6), DIMENSION(:), OPTIONAL, POINTER :: symbol
2074 : INTEGER, DIMENSION(:), OPTIONAL, POINTER :: nq, lq
2075 : REAL(KIND=dp), DIMENSION(:), OPTIONAL, POINTER :: zet
2076 :
2077 : INTEGER :: ns
2078 :
2079 4542 : IF (PRESENT(name)) sto_basis_set%name = name
2080 4542 : IF (PRESENT(nshell)) sto_basis_set%nshell = nshell
2081 4542 : IF (PRESENT(symbol)) THEN
2082 4402 : ns = SIZE(symbol)
2083 4402 : IF (ASSOCIATED(sto_basis_set%symbol)) DEALLOCATE (sto_basis_set%symbol)
2084 13206 : ALLOCATE (sto_basis_set%symbol(1:ns))
2085 18260 : sto_basis_set%symbol(:) = symbol(:)
2086 : END IF
2087 4542 : IF (PRESENT(nq)) THEN
2088 4542 : ns = SIZE(nq)
2089 4542 : CALL reallocate(sto_basis_set%nq, 1, ns)
2090 18540 : sto_basis_set%nq = nq(:)
2091 : END IF
2092 4542 : IF (PRESENT(lq)) THEN
2093 4542 : ns = SIZE(lq)
2094 4542 : CALL reallocate(sto_basis_set%lq, 1, ns)
2095 18540 : sto_basis_set%lq = lq(:)
2096 : END IF
2097 4542 : IF (PRESENT(zet)) THEN
2098 4542 : ns = SIZE(zet)
2099 4542 : CALL reallocate(sto_basis_set%zet, 1, ns)
2100 18540 : sto_basis_set%zet = zet(:)
2101 : END IF
2102 :
2103 4542 : END SUBROUTINE set_sto_basis_set
2104 :
2105 : ! **************************************************************************************************
2106 : !> \brief ...
2107 : !> \param element_symbol ...
2108 : !> \param basis_set_name ...
2109 : !> \param sto_basis_set ...
2110 : !> \param para_env ...
2111 : !> \param dft_section ...
2112 : ! **************************************************************************************************
2113 4 : SUBROUTINE read_sto_basis_set(element_symbol, basis_set_name, sto_basis_set, para_env, dft_section)
2114 :
2115 : ! Read a Slater-type orbital (STO) basis set from the database file.
2116 :
2117 : CHARACTER(LEN=*), INTENT(IN) :: element_symbol, basis_set_name
2118 : TYPE(sto_basis_set_type), INTENT(INOUT) :: sto_basis_set
2119 : TYPE(mp_para_env_type), POINTER :: para_env
2120 : TYPE(section_vals_type), POINTER :: dft_section
2121 :
2122 : CHARACTER(LEN=10) :: nlsym
2123 : CHARACTER(LEN=2) :: lsym
2124 : CHARACTER(LEN=240) :: line
2125 : CHARACTER(LEN=242) :: line2
2126 : CHARACTER(len=default_path_length) :: basis_set_file_name, tmp
2127 : CHARACTER(LEN=default_path_length), DIMENSION(:), &
2128 4 : POINTER :: cbasis
2129 4 : CHARACTER(LEN=LEN(basis_set_name)) :: bsname
2130 4 : CHARACTER(LEN=LEN(basis_set_name)+2) :: bsname2
2131 4 : CHARACTER(LEN=LEN(element_symbol)) :: symbol
2132 4 : CHARACTER(LEN=LEN(element_symbol)+2) :: symbol2
2133 : INTEGER :: ibasis, irep, iset, nbasis, nq, nset, &
2134 : strlen1, strlen2
2135 : LOGICAL :: basis_found, found, match
2136 : REAL(KIND=dp) :: zet
2137 : TYPE(cp_parser_type) :: parser
2138 :
2139 4 : line = ""
2140 4 : line2 = ""
2141 4 : symbol = ""
2142 4 : symbol2 = ""
2143 4 : bsname = ""
2144 4 : bsname2 = ""
2145 :
2146 4 : nbasis = 1
2147 :
2148 4 : sto_basis_set%name = basis_set_name
2149 : CALL section_vals_val_get(dft_section, "BASIS_SET_FILE_NAME", &
2150 4 : n_rep_val=nbasis)
2151 12 : ALLOCATE (cbasis(nbasis))
2152 8 : DO ibasis = 1, nbasis
2153 : CALL section_vals_val_get(dft_section, "BASIS_SET_FILE_NAME", &
2154 4 : i_rep_val=ibasis, c_val=cbasis(ibasis))
2155 4 : basis_set_file_name = cbasis(ibasis)
2156 4 : tmp = basis_set_file_name
2157 8 : CALL uppercase(tmp)
2158 : END DO
2159 :
2160 : ! Search for the requested basis set in the basis set file
2161 : ! until the basis set is found or the end of file is reached
2162 :
2163 4 : basis_found = .FALSE.
2164 8 : basis_loop: DO ibasis = 1, nbasis
2165 4 : IF (basis_found) EXIT basis_loop
2166 4 : basis_set_file_name = cbasis(ibasis)
2167 4 : CALL parser_create(parser, basis_set_file_name, para_env=para_env)
2168 :
2169 4 : bsname = basis_set_name
2170 4 : symbol = element_symbol
2171 4 : irep = 0
2172 :
2173 4 : tmp = basis_set_name
2174 4 : CALL uppercase(tmp)
2175 :
2176 4 : IF (tmp .NE. "NONE") THEN
2177 : search_loop: DO
2178 :
2179 6 : CALL parser_search_string(parser, TRIM(bsname), .TRUE., found, line)
2180 6 : IF (found) THEN
2181 6 : CALL uppercase(symbol)
2182 6 : CALL uppercase(bsname)
2183 :
2184 6 : match = .FALSE.
2185 6 : CALL uppercase(line)
2186 : ! Check both the element symbol and the basis set name
2187 6 : line2 = " "//line//" "
2188 6 : symbol2 = " "//TRIM(symbol)//" "
2189 6 : bsname2 = " "//TRIM(bsname)//" "
2190 6 : strlen1 = LEN_TRIM(symbol2) + 1
2191 6 : strlen2 = LEN_TRIM(bsname2) + 1
2192 :
2193 6 : IF ((INDEX(line2, symbol2(:strlen1)) > 0) .AND. &
2194 4 : (INDEX(line2, bsname2(:strlen2)) > 0)) match = .TRUE.
2195 6 : IF (match) THEN
2196 : ! Read the basis set information
2197 4 : CALL parser_get_object(parser, nset, newline=.TRUE.)
2198 4 : sto_basis_set%nshell = nset
2199 :
2200 4 : CALL reallocate(sto_basis_set%nq, 1, nset)
2201 4 : CALL reallocate(sto_basis_set%lq, 1, nset)
2202 4 : CALL reallocate(sto_basis_set%zet, 1, nset)
2203 12 : ALLOCATE (sto_basis_set%symbol(nset))
2204 :
2205 20 : DO iset = 1, nset
2206 16 : CALL parser_get_object(parser, nq, newline=.TRUE.)
2207 16 : CALL parser_get_object(parser, lsym)
2208 16 : CALL parser_get_object(parser, zet)
2209 16 : sto_basis_set%nq(iset) = nq
2210 16 : sto_basis_set%zet(iset) = zet
2211 16 : WRITE (nlsym, "(I2,A)") nq, TRIM(lsym)
2212 16 : sto_basis_set%symbol(iset) = TRIM(nlsym)
2213 20 : SELECT CASE (TRIM(lsym))
2214 : CASE ("S", "s")
2215 12 : sto_basis_set%lq(iset) = 0
2216 : CASE ("P", "p")
2217 4 : sto_basis_set%lq(iset) = 1
2218 : CASE ("D", "d")
2219 0 : sto_basis_set%lq(iset) = 2
2220 : CASE ("F", "f")
2221 0 : sto_basis_set%lq(iset) = 3
2222 : CASE ("G", "g")
2223 0 : sto_basis_set%lq(iset) = 4
2224 : CASE ("H", "h")
2225 0 : sto_basis_set%lq(iset) = 5
2226 : CASE ("I", "i", "J", "j")
2227 0 : sto_basis_set%lq(iset) = 6
2228 : CASE ("K", "k")
2229 0 : sto_basis_set%lq(iset) = 7
2230 : CASE ("L", "l")
2231 0 : sto_basis_set%lq(iset) = 8
2232 : CASE ("M", "m")
2233 0 : sto_basis_set%lq(iset) = 9
2234 : CASE DEFAULT
2235 : CALL cp_abort(__LOCATION__, &
2236 : "The requested basis set <"//TRIM(bsname)// &
2237 16 : "> for element <"//TRIM(symbol)//"> has an invalid component: ")
2238 : END SELECT
2239 : END DO
2240 :
2241 : basis_found = .TRUE.
2242 : EXIT search_loop
2243 : END IF
2244 : ELSE
2245 : EXIT search_loop
2246 : END IF
2247 : END DO search_loop
2248 : ELSE
2249 4 : match = .FALSE.
2250 : END IF
2251 :
2252 12 : CALL parser_release(parser)
2253 :
2254 : END DO basis_loop
2255 :
2256 4 : IF (tmp .NE. "NONE") THEN
2257 4 : IF (.NOT. basis_found) THEN
2258 0 : basis_set_file_name = ""
2259 0 : DO ibasis = 1, nbasis
2260 0 : basis_set_file_name = TRIM(basis_set_file_name)//"<"//TRIM(cbasis(ibasis))//"> "
2261 : END DO
2262 : CALL cp_abort(__LOCATION__, &
2263 : "The requested basis set <"//TRIM(bsname)// &
2264 : "> for element <"//TRIM(symbol)//"> was not "// &
2265 : "found in the basis set files "// &
2266 0 : TRIM(basis_set_file_name))
2267 : END IF
2268 : END IF
2269 4 : DEALLOCATE (cbasis)
2270 :
2271 16 : END SUBROUTINE read_sto_basis_set
2272 :
2273 : ! **************************************************************************************************
2274 : !> \brief ...
2275 : !> \param sto_basis_set ...
2276 : !> \param gto_basis_set ...
2277 : !> \param ngauss ...
2278 : !> \param ortho ...
2279 : ! **************************************************************************************************
2280 9092 : SUBROUTINE create_gto_from_sto_basis(sto_basis_set, gto_basis_set, ngauss, ortho)
2281 :
2282 : TYPE(sto_basis_set_type), INTENT(IN) :: sto_basis_set
2283 : TYPE(gto_basis_set_type), POINTER :: gto_basis_set
2284 : INTEGER, INTENT(IN), OPTIONAL :: ngauss
2285 : LOGICAL, INTENT(IN), OPTIONAL :: ortho
2286 :
2287 : INTEGER, PARAMETER :: maxng = 6
2288 :
2289 : CHARACTER(LEN=default_string_length) :: name, sng
2290 : INTEGER :: ipgf, iset, maxl, ng, nset, nshell
2291 4546 : INTEGER, DIMENSION(:), POINTER :: lq, nq
2292 : LOGICAL :: do_ortho
2293 4546 : REAL(KIND=dp), DIMENSION(:), POINTER :: zet
2294 : REAL(KIND=dp), DIMENSION(maxng) :: gcc, zetg
2295 :
2296 4546 : ng = 6
2297 4546 : IF (PRESENT(ngauss)) ng = ngauss
2298 4546 : IF (ng > maxng) CPABORT("Too many Gaussian primitives requested")
2299 4546 : do_ortho = .FALSE.
2300 4546 : IF (PRESENT(ortho)) do_ortho = ortho
2301 :
2302 4546 : CALL allocate_gto_basis_set(gto_basis_set)
2303 :
2304 : CALL get_sto_basis_set(sto_basis_set, name=name, nshell=nshell, nq=nq, &
2305 4546 : lq=lq, zet=zet)
2306 :
2307 18560 : maxl = MAXVAL(lq)
2308 4546 : CALL init_orbital_pointers(maxl)
2309 :
2310 4546 : CALL integer_to_string(ng, sng)
2311 4546 : gto_basis_set%name = TRIM(name)//"_STO-"//TRIM(sng)//"G"
2312 :
2313 4546 : nset = nshell
2314 4546 : gto_basis_set%nset = nset
2315 4546 : CALL reallocate(gto_basis_set%lmax, 1, nset)
2316 4546 : CALL reallocate(gto_basis_set%lmin, 1, nset)
2317 4546 : CALL reallocate(gto_basis_set%npgf, 1, nset)
2318 4546 : CALL reallocate(gto_basis_set%nshell, 1, nset)
2319 4546 : CALL reallocate(gto_basis_set%n, 1, 1, 1, nset)
2320 4546 : CALL reallocate(gto_basis_set%l, 1, 1, 1, nset)
2321 4546 : CALL reallocate(gto_basis_set%zet, 1, ng, 1, nset)
2322 4546 : CALL reallocate(gto_basis_set%gcc, 1, ng, 1, 1, 1, nset)
2323 :
2324 14350 : DO iset = 1, nset
2325 9804 : CALL get_sto_ng(zet(iset), ng, nq(iset), lq(iset), zetg, gcc)
2326 9804 : gto_basis_set%lmax(iset) = lq(iset)
2327 9804 : gto_basis_set%lmin(iset) = lq(iset)
2328 9804 : gto_basis_set%npgf(iset) = ng
2329 9804 : gto_basis_set%nshell(iset) = 1
2330 9804 : gto_basis_set%n(1, iset) = lq(iset) + 1
2331 9804 : gto_basis_set%l(1, iset) = lq(iset)
2332 70532 : DO ipgf = 1, ng
2333 56182 : gto_basis_set%gcc(ipgf, 1, iset) = gcc(ipgf)
2334 65986 : gto_basis_set%zet(ipgf, iset) = zetg(ipgf)
2335 : END DO
2336 : END DO
2337 :
2338 4546 : CALL process_gto_basis(gto_basis_set, do_ortho, nset, maxl)
2339 :
2340 4546 : END SUBROUTINE create_gto_from_sto_basis
2341 :
2342 : ! **************************************************************************************************
2343 : !> \brief ...
2344 : !> \param gto_basis_set ...
2345 : !> \param do_ortho ...
2346 : !> \param nset ...
2347 : !> \param maxl ...
2348 : ! **************************************************************************************************
2349 4658 : SUBROUTINE process_gto_basis(gto_basis_set, do_ortho, nset, maxl)
2350 :
2351 : TYPE(gto_basis_set_type), POINTER :: gto_basis_set
2352 : LOGICAL, INTENT(IN), OPTIONAL :: do_ortho
2353 : INTEGER, INTENT(IN) :: nset, maxl
2354 :
2355 : INTEGER :: i1, i2, ico, iset, jset, l, lshell, m, &
2356 : maxco, ncgf, ng, ngs, np, nsgf
2357 : INTEGER, DIMENSION(0:10) :: mxf
2358 4658 : REAL(KIND=dp), ALLOCATABLE, DIMENSION(:, :) :: gal, zal, zll
2359 :
2360 4658 : ng = gto_basis_set%npgf(1)
2361 14712 : DO iset = 1, nset
2362 10054 : IF ((ng /= gto_basis_set%npgf(iset)) .AND. do_ortho) &
2363 4658 : CPABORT("different number of primitves")
2364 : END DO
2365 :
2366 4658 : IF (do_ortho) THEN
2367 2088 : mxf = 0
2368 6864 : DO iset = 1, nset
2369 4776 : l = gto_basis_set%l(1, iset)
2370 6864 : mxf(l) = mxf(l) + 1
2371 : END DO
2372 25056 : m = MAXVAL(mxf)
2373 2088 : IF (m > 1) THEN
2374 4104 : ALLOCATE (gal(ng, nset), zal(ng, nset), zll(m*ng, 0:maxl))
2375 1368 : DO iset = 1, nset
2376 4560 : zal(1:ng, iset) = gto_basis_set%zet(1:ng, iset)
2377 5016 : gal(1:ng, iset) = gto_basis_set%gcc(1:ng, 1, iset)
2378 : END DO
2379 456 : CALL reallocate(gto_basis_set%zet, 1, m*ng, 1, nset)
2380 456 : CALL reallocate(gto_basis_set%gcc, 1, m*ng, 1, 1, 1, nset)
2381 1368 : DO iset = 1, nset
2382 912 : l = gto_basis_set%l(1, iset)
2383 1368 : gto_basis_set%npgf(iset) = ng*mxf(l)
2384 : END DO
2385 8664 : gto_basis_set%zet = 0.0_dp
2386 9576 : gto_basis_set%gcc = 0.0_dp
2387 4560 : zll = 0.0_dp
2388 456 : mxf = 0
2389 1368 : DO iset = 1, nset
2390 912 : l = gto_basis_set%l(1, iset)
2391 912 : mxf(l) = mxf(l) + 1
2392 912 : i1 = mxf(l)*ng - ng + 1
2393 912 : i2 = mxf(l)*ng
2394 4560 : zll(i1:i2, l) = zal(1:ng, iset)
2395 5016 : gto_basis_set%gcc(i1:i2, 1, iset) = gal(1:ng, iset)
2396 : END DO
2397 1368 : DO iset = 1, nset
2398 912 : l = gto_basis_set%l(1, iset)
2399 8664 : gto_basis_set%zet(:, iset) = zll(:, l)
2400 : END DO
2401 1368 : DO iset = 1, nset
2402 912 : l = gto_basis_set%l(1, iset)
2403 1824 : DO jset = 1, iset - 1
2404 1368 : IF (gto_basis_set%l(1, iset) == l) THEN
2405 456 : m = mxf(l)*ng
2406 : CALL orthofun(gto_basis_set%zet(1:m, iset), gto_basis_set%gcc(1:m, 1, iset), &
2407 456 : gto_basis_set%gcc(1:m, 1, jset), l)
2408 : END IF
2409 : END DO
2410 : END DO
2411 456 : DEALLOCATE (gal, zal, zll)
2412 : END IF
2413 : END IF
2414 :
2415 14712 : ngs = MAXVAL(gto_basis_set%npgf(1:nset))
2416 4658 : CALL reallocate(gto_basis_set%set_radius, 1, nset)
2417 4658 : CALL reallocate(gto_basis_set%pgf_radius, 1, ngs, 1, nset)
2418 4658 : CALL reallocate(gto_basis_set%first_cgf, 1, 1, 1, nset)
2419 4658 : CALL reallocate(gto_basis_set%first_sgf, 1, 1, 1, nset)
2420 4658 : CALL reallocate(gto_basis_set%last_cgf, 1, 1, 1, nset)
2421 4658 : CALL reallocate(gto_basis_set%last_sgf, 1, 1, 1, nset)
2422 4658 : CALL reallocate(gto_basis_set%ncgf_set, 1, nset)
2423 4658 : CALL reallocate(gto_basis_set%nsgf_set, 1, nset)
2424 :
2425 4658 : maxco = 0
2426 4658 : ncgf = 0
2427 4658 : nsgf = 0
2428 :
2429 14712 : DO iset = 1, nset
2430 10054 : gto_basis_set%ncgf_set(iset) = 0
2431 10054 : gto_basis_set%nsgf_set(iset) = 0
2432 10054 : lshell = gto_basis_set%l(1, iset)
2433 10054 : gto_basis_set%first_cgf(1, iset) = ncgf + 1
2434 10054 : ncgf = ncgf + nco(lshell)
2435 10054 : gto_basis_set%last_cgf(1, iset) = ncgf
2436 : gto_basis_set%ncgf_set(iset) = &
2437 10054 : gto_basis_set%ncgf_set(iset) + nco(lshell)
2438 10054 : gto_basis_set%first_sgf(1, iset) = nsgf + 1
2439 10054 : nsgf = nsgf + nso(lshell)
2440 10054 : gto_basis_set%last_sgf(1, iset) = nsgf
2441 : gto_basis_set%nsgf_set(iset) = &
2442 10054 : gto_basis_set%nsgf_set(iset) + nso(lshell)
2443 10054 : ngs = gto_basis_set%npgf(iset)
2444 14712 : maxco = MAX(maxco, ngs*ncoset(lshell))
2445 : END DO
2446 :
2447 4658 : gto_basis_set%ncgf = ncgf
2448 4658 : gto_basis_set%nsgf = nsgf
2449 :
2450 4658 : CALL reallocate(gto_basis_set%cphi, 1, maxco, 1, ncgf)
2451 4658 : CALL reallocate(gto_basis_set%sphi, 1, maxco, 1, nsgf)
2452 4658 : CALL reallocate(gto_basis_set%scon, 1, maxco, 1, nsgf)
2453 4658 : CALL reallocate(gto_basis_set%lx, 1, ncgf)
2454 4658 : CALL reallocate(gto_basis_set%ly, 1, ncgf)
2455 4658 : CALL reallocate(gto_basis_set%lz, 1, ncgf)
2456 4658 : CALL reallocate(gto_basis_set%m, 1, nsgf)
2457 4658 : CALL reallocate(gto_basis_set%norm_cgf, 1, ncgf)
2458 13974 : ALLOCATE (gto_basis_set%cgf_symbol(ncgf))
2459 13974 : ALLOCATE (gto_basis_set%sgf_symbol(nsgf))
2460 :
2461 4658 : ncgf = 0
2462 4658 : nsgf = 0
2463 :
2464 14712 : DO iset = 1, nset
2465 10054 : lshell = gto_basis_set%l(1, iset)
2466 10054 : np = lshell + 1
2467 33968 : DO ico = ncoset(lshell - 1) + 1, ncoset(lshell)
2468 23914 : ncgf = ncgf + 1
2469 23914 : gto_basis_set%lx(ncgf) = indco(1, ico)
2470 23914 : gto_basis_set%ly(ncgf) = indco(2, ico)
2471 23914 : gto_basis_set%lz(ncgf) = indco(3, ico)
2472 : gto_basis_set%cgf_symbol(ncgf) = &
2473 : cgf_symbol(np, (/gto_basis_set%lx(ncgf), &
2474 : gto_basis_set%ly(ncgf), &
2475 105710 : gto_basis_set%lz(ncgf)/))
2476 : END DO
2477 37274 : DO m = -lshell, lshell
2478 22562 : nsgf = nsgf + 1
2479 22562 : gto_basis_set%m(nsgf) = m
2480 32616 : gto_basis_set%sgf_symbol(nsgf) = sgf_symbol(np, lshell, m)
2481 : END DO
2482 : END DO
2483 :
2484 4658 : gto_basis_set%norm_type = -1
2485 :
2486 4658 : END SUBROUTINE process_gto_basis
2487 :
2488 : ! **************************************************************************************************
2489 : !> \brief ...
2490 : !> \param zet ...
2491 : !> \param co ...
2492 : !> \param cr ...
2493 : !> \param l ...
2494 : ! **************************************************************************************************
2495 456 : SUBROUTINE orthofun(zet, co, cr, l)
2496 : REAL(KIND=dp), DIMENSION(:), INTENT(IN) :: zet
2497 : REAL(KIND=dp), DIMENSION(:), INTENT(INOUT) :: co, cr
2498 : INTEGER, INTENT(IN) :: l
2499 :
2500 : REAL(KIND=dp) :: ss
2501 :
2502 456 : CALL aovlp(l, zet, cr, cr, ss)
2503 4104 : cr(:) = cr(:)/SQRT(ss)
2504 456 : CALL aovlp(l, zet, co, cr, ss)
2505 4104 : co(:) = co(:) - ss*cr(:)
2506 456 : CALL aovlp(l, zet, co, co, ss)
2507 4104 : co(:) = co(:)/SQRT(ss)
2508 :
2509 456 : END SUBROUTINE orthofun
2510 :
2511 : ! **************************************************************************************************
2512 : !> \brief ...
2513 : !> \param l ...
2514 : !> \param zet ...
2515 : !> \param ca ...
2516 : !> \param cb ...
2517 : !> \param ss ...
2518 : ! **************************************************************************************************
2519 1368 : SUBROUTINE aovlp(l, zet, ca, cb, ss)
2520 : INTEGER, INTENT(IN) :: l
2521 : REAL(KIND=dp), DIMENSION(:), INTENT(IN) :: zet, ca, cb
2522 : REAL(KIND=dp), INTENT(OUT) :: ss
2523 :
2524 : INTEGER :: i, j, m
2525 : REAL(KIND=dp) :: ab, ai, aj, s00, sss
2526 :
2527 : !
2528 : ! use init_norm_cgf_orb
2529 : !
2530 1368 : m = SIZE(zet)
2531 1368 : ss = 0.0_dp
2532 12312 : DO i = 1, m
2533 10944 : ai = (2.0_dp*zet(i)/pi)**0.75_dp
2534 99864 : DO j = 1, m
2535 87552 : aj = (2.0_dp*zet(j)/pi)**0.75_dp
2536 87552 : ab = 1._dp/(zet(i) + zet(j))
2537 87552 : s00 = ai*aj*(pi*ab)**1.50_dp
2538 87552 : IF (l == 0) THEN
2539 : sss = s00
2540 0 : ELSEIF (l == 1) THEN
2541 0 : sss = s00*ab*0.5_dp
2542 : ELSE
2543 0 : CPABORT("aovlp lvalue")
2544 : END IF
2545 98496 : ss = ss + sss*ca(i)*cb(j)
2546 : END DO
2547 : END DO
2548 :
2549 1368 : END SUBROUTINE aovlp
2550 :
2551 : ! **************************************************************************************************
2552 : !> \brief ...
2553 : !> \param z ...
2554 : !> \param ne ...
2555 : !> \param n ...
2556 : !> \param l ...
2557 : !> \return ...
2558 : ! **************************************************************************************************
2559 21588 : PURE FUNCTION srules(z, ne, n, l)
2560 : ! Slater rules
2561 : INTEGER, INTENT(IN) :: z
2562 : INTEGER, DIMENSION(:, :), INTENT(IN) :: ne
2563 : INTEGER, INTENT(IN) :: n, l
2564 : REAL(dp) :: srules
2565 :
2566 : REAL(dp), DIMENSION(7), PARAMETER :: &
2567 : xns = (/1.0_dp, 2.0_dp, 3.0_dp, 3.7_dp, 4.0_dp, 4.2_dp, 4.4_dp/)
2568 :
2569 : INTEGER :: i, l1, l2, m, m1, m2, nn
2570 : REAL(dp) :: s
2571 :
2572 21588 : s = 0.0_dp
2573 : ! The complete shell
2574 21588 : l1 = MIN(l + 1, 4)
2575 21588 : nn = MIN(n, 7)
2576 : IF (l1 == 1) l2 = 2
2577 21588 : IF (l1 == 2) l2 = 1
2578 15229 : IF (l1 == 3) l2 = 4
2579 20157 : IF (l1 == 4) l2 = 3
2580 : ! Rule a) no contribution from shells further out
2581 : ! Rule b) 0.35 (1s 0.3) from each other electron in the same shell
2582 21588 : IF (n == 1) THEN
2583 5958 : m = ne(1, 1)
2584 5958 : s = s + 0.3_dp*REAL(m - 1, dp)
2585 : ELSE
2586 15630 : m = ne(l1, nn) + ne(l2, nn)
2587 15630 : s = s + 0.35_dp*REAL(m - 1, dp)
2588 : END IF
2589 : ! Rule c) if (s,p) shell 0.85 from each electron with n-1, and 1.0
2590 : ! from all electrons further in
2591 21588 : IF (l1 + l2 == 3) THEN
2592 19972 : IF (nn > 1) THEN
2593 14014 : m1 = ne(1, nn - 1) + ne(2, nn - 1) + ne(3, nn - 1) + ne(4, nn - 1)
2594 14014 : m2 = 0
2595 21952 : DO i = 1, nn - 2
2596 21952 : m2 = m2 + ne(1, i) + ne(2, i) + ne(3, i) + ne(4, I)
2597 : END DO
2598 14014 : s = s + 0.85_dp*REAL(m1, dp) + 1._dp*REAL(m2, dp)
2599 : END IF
2600 : ELSE
2601 : ! Rule d) if (d,f) shell 1.0 from each electron inside
2602 : m = 0
2603 6562 : DO i = 1, nn - 1
2604 6562 : m = m + ne(1, i) + ne(2, i) + ne(3, i) + ne(4, i)
2605 : END DO
2606 1616 : s = s + 1._dp*REAL(m, dp)
2607 : END IF
2608 : ! Slater exponent is (Z-S)/NS
2609 21588 : srules = (REAL(z, dp) - s)/xns(nn)
2610 21588 : END FUNCTION srules
2611 :
2612 : ! **************************************************************************************************
2613 : !> \brief sort basis sets w.r.t. radius
2614 : !> \param basis_set ...
2615 : !> \param sort_method ...
2616 : ! **************************************************************************************************
2617 11632 : SUBROUTINE sort_gto_basis_set(basis_set, sort_method)
2618 : TYPE(gto_basis_set_type), INTENT(INOUT) :: basis_set
2619 : INTEGER, INTENT(IN) :: sort_method
2620 :
2621 11632 : CHARACTER(LEN=12), DIMENSION(:), POINTER :: cgf_symbol
2622 11632 : CHARACTER(LEN=6), DIMENSION(:), POINTER :: sgf_symbol
2623 : INTEGER :: ic, ic_max, icgf, icgf_new, icgf_old, ico, is, is_max, iset, isgf, isgf_new, &
2624 : isgf_old, ishell, lshell, maxco, maxpgf, maxshell, mm, nc, ncgf, ns, nset
2625 11632 : INTEGER, ALLOCATABLE, DIMENSION(:) :: sort_index
2626 11632 : INTEGER, ALLOCATABLE, DIMENSION(:, :) :: icgf_set, isgf_set
2627 11632 : INTEGER, DIMENSION(:), POINTER :: lx, ly, lz, m, npgf
2628 11632 : REAL(dp), ALLOCATABLE, DIMENSION(:) :: tmp
2629 11632 : REAL(dp), DIMENSION(:), POINTER :: set_radius
2630 11632 : REAL(dp), DIMENSION(:, :), POINTER :: zet
2631 11632 : REAL(KIND=dp), DIMENSION(:), POINTER :: norm_cgf
2632 11632 : REAL(KIND=dp), DIMENSION(:, :), POINTER :: cphi, scon, sphi
2633 :
2634 11632 : NULLIFY (set_radius, zet)
2635 :
2636 11632 : IF (sort_method == basis_sort_default) RETURN
2637 :
2638 : CALL get_gto_basis_set(gto_basis_set=basis_set, &
2639 : nset=nset, &
2640 : maxshell=maxshell, &
2641 : maxpgf=maxpgf, &
2642 : maxco=maxco, &
2643 : ncgf=ncgf, &
2644 : npgf=npgf, &
2645 : set_radius=set_radius, &
2646 712 : zet=zet)
2647 :
2648 2136 : ALLOCATE (sort_index(nset))
2649 2136 : ALLOCATE (tmp(nset))
2650 : SELECT CASE (sort_method)
2651 : CASE (basis_sort_zet)
2652 4502 : DO iset = 1, nset
2653 13846 : tmp(iset) = MINVAL(basis_set%zet(:npgf(iset), iset))
2654 : END DO
2655 : CASE DEFAULT
2656 712 : CPABORT("Request basis sort criterion not implemented.")
2657 : END SELECT
2658 :
2659 712 : CALL sort(tmp(1:nset), nset, sort_index)
2660 :
2661 712 : ic_max = 0
2662 712 : is_max = 0
2663 4502 : DO iset = 1, nset
2664 3790 : ic = 0
2665 3790 : is = 0
2666 10110 : DO ishell = 1, basis_set%nshell(iset)
2667 22196 : DO ico = 1, nco(basis_set%l(ishell, iset))
2668 16588 : ic = ic + 1
2669 22196 : IF (ic > ic_max) ic_max = ic
2670 : END DO
2671 5608 : lshell = basis_set%l(ishell, iset)
2672 24202 : DO mm = -lshell, lshell
2673 14804 : is = is + 1
2674 20412 : IF (is > is_max) is_max = is
2675 : END DO
2676 : END DO
2677 : END DO
2678 :
2679 712 : icgf = 0
2680 712 : isgf = 0
2681 2848 : ALLOCATE (icgf_set(nset, ic_max))
2682 42128 : icgf_set(:, :) = 0
2683 2848 : ALLOCATE (isgf_set(nset, is_max))
2684 34360 : isgf_set(:, :) = 0
2685 :
2686 4502 : DO iset = 1, nset
2687 3790 : ic = 0
2688 3790 : is = 0
2689 10110 : DO ishell = 1, basis_set%nshell(iset)
2690 22196 : DO ico = 1, nco(basis_set%l(ishell, iset))
2691 16588 : icgf = icgf + 1
2692 16588 : ic = ic + 1
2693 22196 : icgf_set(iset, ic) = icgf
2694 : END DO
2695 5608 : lshell = basis_set%l(ishell, iset)
2696 24202 : DO mm = -lshell, lshell
2697 14804 : isgf = isgf + 1
2698 14804 : is = is + 1
2699 20412 : isgf_set(iset, is) = isgf
2700 : END DO
2701 : END DO
2702 : END DO
2703 :
2704 2136 : ALLOCATE (cgf_symbol(SIZE(basis_set%cgf_symbol)))
2705 2136 : ALLOCATE (norm_cgf(SIZE(basis_set%norm_cgf)))
2706 2136 : ALLOCATE (lx(SIZE(basis_set%lx)))
2707 2136 : ALLOCATE (ly(SIZE(basis_set%ly)))
2708 2136 : ALLOCATE (lz(SIZE(basis_set%lz)))
2709 2848 : ALLOCATE (cphi(SIZE(basis_set%cphi, 1), SIZE(basis_set%cphi, 2)))
2710 606646 : cphi = 0.0_dp
2711 2848 : ALLOCATE (sphi(SIZE(basis_set%sphi, 1), SIZE(basis_set%sphi, 2)))
2712 535426 : sphi = 0.0_dp
2713 2848 : ALLOCATE (scon(SIZE(basis_set%scon, 1), SIZE(basis_set%scon, 2)))
2714 535426 : scon = 0.0_dp
2715 :
2716 2136 : ALLOCATE (sgf_symbol(SIZE(basis_set%sgf_symbol)))
2717 2136 : ALLOCATE (m(SIZE(basis_set%m)))
2718 :
2719 : icgf_new = 0
2720 : isgf_new = 0
2721 4502 : DO iset = 1, nset
2722 37596 : DO ic = 1, ic_max
2723 33806 : icgf_old = icgf_set(sort_index(iset), ic)
2724 33806 : IF (icgf_old == 0) CYCLE
2725 16588 : icgf_new = icgf_new + 1
2726 16588 : norm_cgf(icgf_new) = basis_set%norm_cgf(icgf_old)
2727 16588 : lx(icgf_new) = basis_set%lx(icgf_old)
2728 16588 : ly(icgf_new) = basis_set%ly(icgf_old)
2729 16588 : lz(icgf_new) = basis_set%lz(icgf_old)
2730 605934 : cphi(:, icgf_new) = basis_set%cphi(:, icgf_old)
2731 37596 : cgf_symbol(icgf_new) = basis_set%cgf_symbol(icgf_old)
2732 : END DO
2733 31474 : DO is = 1, is_max
2734 26972 : isgf_old = isgf_set(sort_index(iset), is)
2735 26972 : IF (isgf_old == 0) CYCLE
2736 14804 : isgf_new = isgf_new + 1
2737 14804 : m(isgf_new) = basis_set%m(isgf_old)
2738 534714 : sphi(:, isgf_new) = basis_set%sphi(:, isgf_old)
2739 534714 : scon(:, isgf_new) = basis_set%scon(:, isgf_old)
2740 30762 : sgf_symbol(isgf_new) = basis_set%sgf_symbol(isgf_old)
2741 : END DO
2742 : END DO
2743 :
2744 712 : DEALLOCATE (basis_set%cgf_symbol)
2745 712 : basis_set%cgf_symbol => cgf_symbol
2746 712 : DEALLOCATE (basis_set%norm_cgf)
2747 712 : basis_set%norm_cgf => norm_cgf
2748 712 : DEALLOCATE (basis_set%lx)
2749 712 : basis_set%lx => lx
2750 712 : DEALLOCATE (basis_set%ly)
2751 712 : basis_set%ly => ly
2752 712 : DEALLOCATE (basis_set%lz)
2753 712 : basis_set%lz => lz
2754 712 : DEALLOCATE (basis_set%cphi)
2755 712 : basis_set%cphi => cphi
2756 712 : DEALLOCATE (basis_set%sphi)
2757 712 : basis_set%sphi => sphi
2758 712 : DEALLOCATE (basis_set%scon)
2759 712 : basis_set%scon => scon
2760 :
2761 712 : DEALLOCATE (basis_set%m)
2762 712 : basis_set%m => m
2763 712 : DEALLOCATE (basis_set%sgf_symbol)
2764 712 : basis_set%sgf_symbol => sgf_symbol
2765 :
2766 8292 : basis_set%lmax = basis_set%lmax(sort_index)
2767 8292 : basis_set%lmin = basis_set%lmin(sort_index)
2768 8292 : basis_set%npgf = basis_set%npgf(sort_index)
2769 8292 : basis_set%nshell = basis_set%nshell(sort_index)
2770 8292 : basis_set%ncgf_set = basis_set%ncgf_set(sort_index)
2771 8292 : basis_set%nsgf_set = basis_set%nsgf_set(sort_index)
2772 :
2773 22184 : basis_set%n(:, :) = basis_set%n(:, sort_index)
2774 22184 : basis_set%l(:, :) = basis_set%l(:, sort_index)
2775 24944 : basis_set%zet(:, :) = basis_set%zet(:, sort_index)
2776 :
2777 92204 : basis_set%gcc(:, :, :) = basis_set%gcc(:, :, sort_index)
2778 8292 : basis_set%set_radius(:) = basis_set%set_radius(sort_index)
2779 24944 : basis_set%pgf_radius(:, :) = basis_set%pgf_radius(:, sort_index)
2780 :
2781 712 : nc = 0
2782 712 : ns = 0
2783 4502 : DO iset = 1, nset
2784 10110 : DO ishell = 1, basis_set%nshell(iset)
2785 5608 : lshell = basis_set%l(ishell, iset)
2786 5608 : basis_set%first_cgf(ishell, iset) = nc + 1
2787 5608 : nc = nc + nco(lshell)
2788 5608 : basis_set%last_cgf(ishell, iset) = nc
2789 5608 : basis_set%first_sgf(ishell, iset) = ns + 1
2790 5608 : ns = ns + nso(lshell)
2791 9398 : basis_set%last_sgf(ishell, iset) = ns
2792 : END DO
2793 : END DO
2794 :
2795 12344 : END SUBROUTINE
2796 :
2797 0 : END MODULE basis_set_types
|