Line data Source code
1 : !--------------------------------------------------------------------------------------------------!
2 : ! CP2K: A general program to perform molecular dynamics simulations !
3 : ! Copyright 2000-2024 CP2K developers group <https://cp2k.org> !
4 : ! !
5 : ! SPDX-License-Identifier: GPL-2.0-or-later !
6 : !--------------------------------------------------------------------------------------------------!
7 :
8 : ! **************************************************************************************************
9 : !> \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
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 19600 : 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 19600 : CALL deallocate_gto_basis_set(gto_basis_set)
159 :
160 19600 : ALLOCATE (gto_basis_set)
161 :
162 19600 : END SUBROUTINE allocate_gto_basis_set
163 :
164 : ! **************************************************************************************************
165 : !> \brief ...
166 : !> \param gto_basis_set ...
167 : ! **************************************************************************************************
168 39508 : 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 39508 : IF (ASSOCIATED(gto_basis_set)) THEN
177 19908 : IF (ASSOCIATED(gto_basis_set%cgf_symbol)) DEALLOCATE (gto_basis_set%cgf_symbol)
178 19908 : IF (ASSOCIATED(gto_basis_set%sgf_symbol)) DEALLOCATE (gto_basis_set%sgf_symbol)
179 19908 : IF (ASSOCIATED(gto_basis_set%norm_cgf)) DEALLOCATE (gto_basis_set%norm_cgf)
180 19908 : IF (ASSOCIATED(gto_basis_set%set_radius)) DEALLOCATE (gto_basis_set%set_radius)
181 19908 : IF (ASSOCIATED(gto_basis_set%lmax)) DEALLOCATE (gto_basis_set%lmax)
182 19908 : IF (ASSOCIATED(gto_basis_set%lmin)) DEALLOCATE (gto_basis_set%lmin)
183 19908 : IF (ASSOCIATED(gto_basis_set%lx)) DEALLOCATE (gto_basis_set%lx)
184 19908 : IF (ASSOCIATED(gto_basis_set%ly)) DEALLOCATE (gto_basis_set%ly)
185 19908 : IF (ASSOCIATED(gto_basis_set%lz)) DEALLOCATE (gto_basis_set%lz)
186 19908 : IF (ASSOCIATED(gto_basis_set%m)) DEALLOCATE (gto_basis_set%m)
187 19908 : IF (ASSOCIATED(gto_basis_set%ncgf_set)) DEALLOCATE (gto_basis_set%ncgf_set)
188 19908 : IF (ASSOCIATED(gto_basis_set%npgf)) DEALLOCATE (gto_basis_set%npgf)
189 19908 : IF (ASSOCIATED(gto_basis_set%nsgf_set)) DEALLOCATE (gto_basis_set%nsgf_set)
190 19908 : IF (ASSOCIATED(gto_basis_set%nshell)) DEALLOCATE (gto_basis_set%nshell)
191 19908 : IF (ASSOCIATED(gto_basis_set%cphi)) DEALLOCATE (gto_basis_set%cphi)
192 19908 : IF (ASSOCIATED(gto_basis_set%pgf_radius)) DEALLOCATE (gto_basis_set%pgf_radius)
193 19908 : IF (ASSOCIATED(gto_basis_set%sphi)) DEALLOCATE (gto_basis_set%sphi)
194 19908 : IF (ASSOCIATED(gto_basis_set%scon)) DEALLOCATE (gto_basis_set%scon)
195 19908 : IF (ASSOCIATED(gto_basis_set%zet)) DEALLOCATE (gto_basis_set%zet)
196 19908 : IF (ASSOCIATED(gto_basis_set%first_cgf)) DEALLOCATE (gto_basis_set%first_cgf)
197 19908 : IF (ASSOCIATED(gto_basis_set%first_sgf)) DEALLOCATE (gto_basis_set%first_sgf)
198 19908 : IF (ASSOCIATED(gto_basis_set%l)) DEALLOCATE (gto_basis_set%l)
199 19908 : IF (ASSOCIATED(gto_basis_set%last_cgf)) DEALLOCATE (gto_basis_set%last_cgf)
200 19908 : IF (ASSOCIATED(gto_basis_set%last_sgf)) DEALLOCATE (gto_basis_set%last_sgf)
201 19908 : IF (ASSOCIATED(gto_basis_set%n)) DEALLOCATE (gto_basis_set%n)
202 19908 : IF (ASSOCIATED(gto_basis_set%gcc)) DEALLOCATE (gto_basis_set%gcc)
203 19908 : DEALLOCATE (gto_basis_set)
204 : END IF
205 39508 : END SUBROUTINE deallocate_gto_basis_set
206 :
207 : ! **************************************************************************************************
208 : !> \brief ...
209 : !> \param basis_set_in ...
210 : !> \param basis_set_out ...
211 : ! **************************************************************************************************
212 2950 : 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 2950 : CALL allocate_gto_basis_set(basis_set_out)
222 :
223 2950 : basis_set_out%name = basis_set_in%name
224 2950 : basis_set_out%aliases = basis_set_in%aliases
225 2950 : basis_set_out%kind_radius = basis_set_in%kind_radius
226 2950 : basis_set_out%norm_type = basis_set_in%norm_type
227 2950 : basis_set_out%nset = basis_set_in%nset
228 2950 : basis_set_out%ncgf = basis_set_in%ncgf
229 2950 : basis_set_out%nsgf = basis_set_in%nsgf
230 2950 : nset = basis_set_in%nset
231 2950 : ncgf = basis_set_in%ncgf
232 2950 : nsgf = basis_set_in%nsgf
233 8850 : ALLOCATE (basis_set_out%cgf_symbol(ncgf))
234 8850 : ALLOCATE (basis_set_out%sgf_symbol(nsgf))
235 47600 : basis_set_out%cgf_symbol = basis_set_in%cgf_symbol
236 42452 : basis_set_out%sgf_symbol = basis_set_in%sgf_symbol
237 8850 : ALLOCATE (basis_set_out%norm_cgf(ncgf))
238 47600 : basis_set_out%norm_cgf = basis_set_in%norm_cgf
239 8850 : ALLOCATE (basis_set_out%set_radius(nset))
240 11808 : basis_set_out%set_radius = basis_set_in%set_radius
241 26550 : ALLOCATE (basis_set_out%lmax(nset), basis_set_out%lmin(nset), basis_set_out%npgf(nset), basis_set_out%nshell(nset))
242 11808 : basis_set_out%lmax = basis_set_in%lmax
243 11808 : basis_set_out%lmin = basis_set_in%lmin
244 11808 : basis_set_out%npgf = basis_set_in%npgf
245 11808 : basis_set_out%nshell = basis_set_in%nshell
246 26550 : ALLOCATE (basis_set_out%lx(ncgf), basis_set_out%ly(ncgf), basis_set_out%lz(ncgf), basis_set_out%m(nsgf))
247 47600 : basis_set_out%lx = basis_set_in%lx
248 47600 : basis_set_out%ly = basis_set_in%ly
249 47600 : basis_set_out%lz = basis_set_in%lz
250 42452 : basis_set_out%m = basis_set_in%m
251 14750 : ALLOCATE (basis_set_out%ncgf_set(nset), basis_set_out%nsgf_set(nset))
252 11808 : basis_set_out%ncgf_set = basis_set_in%ncgf_set
253 11808 : basis_set_out%nsgf_set = basis_set_in%nsgf_set
254 2950 : maxco = SIZE(basis_set_in%cphi, 1)
255 29500 : ALLOCATE (basis_set_out%cphi(maxco, ncgf), basis_set_out%sphi(maxco, nsgf), basis_set_out%scon(maxco, nsgf))
256 1240912 : basis_set_out%cphi = basis_set_in%cphi
257 1044184 : basis_set_out%sphi = basis_set_in%sphi
258 1044184 : basis_set_out%scon = basis_set_in%scon
259 11808 : maxpgf = MAXVAL(basis_set_in%npgf)
260 20650 : ALLOCATE (basis_set_out%pgf_radius(maxpgf, nset), basis_set_out%zet(maxpgf, nset))
261 44736 : basis_set_out%pgf_radius = basis_set_in%pgf_radius
262 44736 : basis_set_out%zet = basis_set_in%zet
263 11808 : maxshell = MAXVAL(basis_set_in%nshell)
264 20650 : ALLOCATE (basis_set_out%first_cgf(maxshell, nset), basis_set_out%first_sgf(maxshell, nset))
265 20650 : ALLOCATE (basis_set_out%last_cgf(maxshell, nset), basis_set_out%last_sgf(maxshell, nset))
266 32074 : basis_set_out%first_cgf = basis_set_in%first_cgf
267 32074 : basis_set_out%first_sgf = basis_set_in%first_sgf
268 32074 : basis_set_out%last_cgf = basis_set_in%last_cgf
269 32074 : basis_set_out%last_sgf = basis_set_in%last_sgf
270 20650 : ALLOCATE (basis_set_out%n(maxshell, nset), basis_set_out%l(maxshell, nset))
271 32074 : basis_set_out%n = basis_set_in%n
272 32074 : basis_set_out%l = basis_set_in%l
273 14750 : ALLOCATE (basis_set_out%gcc(maxpgf, maxshell, nset))
274 115280 : basis_set_out%gcc = basis_set_in%gcc
275 :
276 2950 : 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 : ! **************************************************************************************************
618 32154267 : SUBROUTINE get_gto_basis_set(gto_basis_set, name, aliases, norm_type, kind_radius, ncgf, &
619 : nset, nsgf, cgf_symbol, sgf_symbol, norm_cgf, set_radius, lmax, lmin, lx, ly, lz, &
620 : m, ncgf_set, npgf, nsgf_set, nshell, cphi, pgf_radius, sphi, scon, zet, first_cgf, first_sgf, l, &
621 : last_cgf, last_sgf, n, gcc, maxco, maxl, maxpgf, maxsgf_set, maxshell, maxso, nco_sum, &
622 : npgf_sum, nshell_sum, maxder, short_kind_radius)
623 :
624 : ! Get informations about a Gaussian-type orbital (GTO) basis set.
625 :
626 : ! - Creation (10.01.2002,MK)
627 :
628 : TYPE(gto_basis_set_type), INTENT(IN) :: gto_basis_set
629 : CHARACTER(LEN=default_string_length), &
630 : INTENT(OUT), OPTIONAL :: name, aliases
631 : INTEGER, INTENT(OUT), OPTIONAL :: norm_type
632 : REAL(KIND=dp), INTENT(OUT), OPTIONAL :: kind_radius
633 : INTEGER, INTENT(OUT), OPTIONAL :: ncgf, nset, nsgf
634 : CHARACTER(LEN=12), DIMENSION(:), OPTIONAL, POINTER :: cgf_symbol
635 : CHARACTER(LEN=6), DIMENSION(:), OPTIONAL, POINTER :: sgf_symbol
636 : REAL(KIND=dp), DIMENSION(:), OPTIONAL, POINTER :: norm_cgf, set_radius
637 : INTEGER, DIMENSION(:), OPTIONAL, POINTER :: lmax, lmin, lx, ly, lz, m, ncgf_set, &
638 : npgf, nsgf_set, nshell
639 : REAL(KIND=dp), DIMENSION(:, :), OPTIONAL, POINTER :: cphi, pgf_radius, sphi, scon, zet
640 : INTEGER, DIMENSION(:, :), OPTIONAL, POINTER :: first_cgf, first_sgf, l, last_cgf, &
641 : last_sgf, n
642 : REAL(KIND=dp), DIMENSION(:, :, :), OPTIONAL, &
643 : POINTER :: gcc
644 : INTEGER, INTENT(OUT), OPTIONAL :: maxco, maxl, maxpgf, maxsgf_set, &
645 : maxshell, maxso, nco_sum, npgf_sum, &
646 : nshell_sum
647 : INTEGER, INTENT(IN), OPTIONAL :: maxder
648 : REAL(KIND=dp), INTENT(OUT), OPTIONAL :: short_kind_radius
649 :
650 : INTEGER :: iset, nder
651 :
652 32154267 : IF (PRESENT(name)) name = gto_basis_set%name
653 32154267 : IF (PRESENT(aliases)) aliases = gto_basis_set%aliases
654 32154267 : IF (PRESENT(norm_type)) norm_type = gto_basis_set%norm_type
655 32154267 : IF (PRESENT(kind_radius)) kind_radius = gto_basis_set%kind_radius
656 32154267 : IF (PRESENT(short_kind_radius)) short_kind_radius = gto_basis_set%short_kind_radius
657 32154267 : IF (PRESENT(ncgf)) ncgf = gto_basis_set%ncgf
658 32154267 : IF (PRESENT(nset)) nset = gto_basis_set%nset
659 32154267 : IF (PRESENT(nsgf)) nsgf = gto_basis_set%nsgf
660 32154267 : IF (PRESENT(cgf_symbol)) cgf_symbol => gto_basis_set%cgf_symbol
661 32154267 : IF (PRESENT(sgf_symbol)) sgf_symbol => gto_basis_set%sgf_symbol
662 32154267 : IF (PRESENT(norm_cgf)) norm_cgf => gto_basis_set%norm_cgf
663 32154267 : IF (PRESENT(set_radius)) set_radius => gto_basis_set%set_radius
664 32154267 : IF (PRESENT(lmax)) lmax => gto_basis_set%lmax
665 32154267 : IF (PRESENT(lmin)) lmin => gto_basis_set%lmin
666 32154267 : IF (PRESENT(lx)) lx => gto_basis_set%lx
667 32154267 : IF (PRESENT(ly)) ly => gto_basis_set%ly
668 32154267 : IF (PRESENT(lz)) lz => gto_basis_set%lz
669 32154267 : IF (PRESENT(m)) m => gto_basis_set%m
670 32154267 : IF (PRESENT(ncgf_set)) ncgf_set => gto_basis_set%ncgf_set
671 32154267 : IF (PRESENT(npgf)) npgf => gto_basis_set%npgf
672 32154267 : IF (PRESENT(nsgf_set)) nsgf_set => gto_basis_set%nsgf_set
673 32154267 : IF (PRESENT(nshell)) nshell => gto_basis_set%nshell
674 32154267 : IF (PRESENT(cphi)) cphi => gto_basis_set%cphi
675 32154267 : IF (PRESENT(pgf_radius)) pgf_radius => gto_basis_set%pgf_radius
676 32154267 : IF (PRESENT(sphi)) sphi => gto_basis_set%sphi
677 32154267 : IF (PRESENT(scon)) scon => gto_basis_set%scon
678 32154267 : IF (PRESENT(zet)) zet => gto_basis_set%zet
679 32154267 : IF (PRESENT(first_cgf)) first_cgf => gto_basis_set%first_cgf
680 32154267 : IF (PRESENT(first_sgf)) first_sgf => gto_basis_set%first_sgf
681 32154267 : IF (PRESENT(l)) l => gto_basis_set%l
682 32154267 : IF (PRESENT(last_cgf)) last_cgf => gto_basis_set%last_cgf
683 32154267 : IF (PRESENT(last_sgf)) last_sgf => gto_basis_set%last_sgf
684 32154267 : IF (PRESENT(n)) n => gto_basis_set%n
685 32154267 : IF (PRESENT(gcc)) gcc => gto_basis_set%gcc
686 32154267 : IF (PRESENT(maxco)) THEN
687 7025465 : maxco = 0
688 7025465 : IF (PRESENT(maxder)) THEN
689 0 : nder = maxder
690 : ELSE
691 : nder = 0
692 : END IF
693 23065232 : DO iset = 1, gto_basis_set%nset
694 : maxco = MAX(maxco, gto_basis_set%npgf(iset)* &
695 23065232 : ncoset(gto_basis_set%lmax(iset) + nder))
696 : END DO
697 : END IF
698 32154267 : IF (PRESENT(maxl)) THEN
699 6907376 : maxl = -1
700 21856091 : DO iset = 1, gto_basis_set%nset
701 21856091 : maxl = MAX(maxl, gto_basis_set%lmax(iset))
702 : END DO
703 : END IF
704 32154267 : IF (PRESENT(maxpgf)) THEN
705 4246 : maxpgf = 0
706 18638 : DO iset = 1, gto_basis_set%nset
707 18638 : maxpgf = MAX(maxpgf, gto_basis_set%npgf(iset))
708 : END DO
709 : END IF
710 32154267 : IF (PRESENT(maxsgf_set)) THEN
711 439109 : maxsgf_set = 0
712 2223210 : DO iset = 1, gto_basis_set%nset
713 2223210 : maxsgf_set = MAX(maxsgf_set, gto_basis_set%nsgf_set(iset))
714 : END DO
715 : END IF
716 32154267 : IF (PRESENT(maxshell)) THEN ! MAXVAL on structure component avoided
717 2628 : maxshell = 0
718 11914 : DO iset = 1, gto_basis_set%nset
719 11914 : maxshell = MAX(maxshell, gto_basis_set%nshell(iset))
720 : END DO
721 : END IF
722 32154267 : IF (PRESENT(maxso)) THEN
723 1791781 : maxso = 0
724 6568338 : DO iset = 1, gto_basis_set%nset
725 : maxso = MAX(maxso, gto_basis_set%npgf(iset)* &
726 6568338 : nsoset(gto_basis_set%lmax(iset)))
727 : END DO
728 : END IF
729 :
730 32154267 : IF (PRESENT(nco_sum)) THEN
731 72924 : nco_sum = 0
732 257030 : DO iset = 1, gto_basis_set%nset
733 : nco_sum = nco_sum + gto_basis_set%npgf(iset)* &
734 257030 : ncoset(gto_basis_set%lmax(iset))
735 : END DO
736 : END IF
737 32168460 : IF (PRESENT(npgf_sum)) npgf_sum = SUM(gto_basis_set%npgf)
738 32168460 : IF (PRESENT(nshell_sum)) nshell_sum = SUM(gto_basis_set%nshell)
739 :
740 32154267 : END SUBROUTINE get_gto_basis_set
741 :
742 : ! **************************************************************************************************
743 : !> \brief ...
744 : !> \param gto_basis_set ...
745 : ! **************************************************************************************************
746 12 : SUBROUTINE init_aux_basis_set(gto_basis_set)
747 :
748 : ! Initialise a Gaussian-type orbital (GTO) basis set data set.
749 :
750 : ! - Creation (06.12.2000,MK)
751 :
752 : TYPE(gto_basis_set_type), POINTER :: gto_basis_set
753 :
754 : CHARACTER(len=*), PARAMETER :: routineN = 'init_aux_basis_set'
755 :
756 : INTEGER :: handle
757 :
758 : ! -------------------------------------------------------------------------
759 :
760 6 : IF (.NOT. ASSOCIATED(gto_basis_set)) RETURN
761 :
762 6 : CALL timeset(routineN, handle)
763 :
764 12 : SELECT CASE (gto_basis_set%norm_type)
765 : CASE (0)
766 : ! No normalisation requested
767 : CASE (1)
768 6 : CALL init_norm_cgf_aux_2(gto_basis_set)
769 : CASE (2)
770 : ! WARNING this was never tested
771 0 : CALL init_norm_cgf_aux(gto_basis_set)
772 : CASE DEFAULT
773 6 : CPABORT("Normalization method not specified")
774 : END SELECT
775 :
776 : ! Initialise the transformation matrices "pgf" -> "cgf"
777 6 : CALL init_cphi_and_sphi(gto_basis_set)
778 :
779 6 : CALL timestop(handle)
780 :
781 : END SUBROUTINE init_aux_basis_set
782 :
783 : ! **************************************************************************************************
784 : !> \brief ...
785 : !> \param gto_basis_set ...
786 : ! **************************************************************************************************
787 16973 : SUBROUTINE init_cphi_and_sphi(gto_basis_set)
788 :
789 : ! Initialise the matrices for the transformation of primitive Cartesian
790 : ! Gaussian-type functions to contracted Cartesian (cphi) and spherical
791 : ! (sphi) Gaussian-type functions.
792 :
793 : ! - Creation (20.09.2000,MK)
794 :
795 : TYPE(gto_basis_set_type), INTENT(INOUT) :: gto_basis_set
796 :
797 : CHARACTER(len=*), PARAMETER :: routineN = 'init_cphi_and_sphi'
798 :
799 : INTEGER :: first_cgf, first_sgf, handle, icgf, ico, &
800 : ipgf, iset, ishell, l, last_sgf, lmax, &
801 : lmin, n, n1, n2, ncgf, nn, nn1, nn2, &
802 : npgf, nsgf
803 :
804 : ! -------------------------------------------------------------------------
805 : ! Build the Cartesian transformation matrix "cphi"
806 :
807 16973 : CALL timeset(routineN, handle)
808 :
809 6460468 : gto_basis_set%cphi = 0.0_dp
810 63506 : DO iset = 1, gto_basis_set%nset
811 46533 : n = ncoset(gto_basis_set%lmax(iset))
812 139599 : DO ishell = 1, gto_basis_set%nshell(iset)
813 282323 : DO icgf = gto_basis_set%first_cgf(ishell, iset), &
814 122626 : gto_basis_set%last_cgf(ishell, iset)
815 : ico = coset(gto_basis_set%lx(icgf), &
816 : gto_basis_set%ly(icgf), &
817 206230 : gto_basis_set%lz(icgf))
818 920591 : DO ipgf = 1, gto_basis_set%npgf(iset)
819 : gto_basis_set%cphi(ico, icgf) = gto_basis_set%norm_cgf(icgf)* &
820 638268 : gto_basis_set%gcc(ipgf, ishell, iset)
821 844498 : ico = ico + n
822 : END DO
823 : END DO
824 : END DO
825 : END DO
826 :
827 : ! Build the spherical transformation matrix "sphi"
828 :
829 16973 : n = SIZE(gto_basis_set%cphi, 1)
830 :
831 5646115 : gto_basis_set%sphi = 0.0_dp
832 16973 : IF (n .GT. 0) THEN
833 16963 : lmax = -1
834 : ! Ensure proper setup of orbtramat
835 63490 : DO iset = 1, gto_basis_set%nset
836 139577 : DO ishell = 1, gto_basis_set%nshell(iset)
837 122614 : lmax = MAX(lmax, gto_basis_set%l(ishell, iset))
838 : END DO
839 : END DO
840 16963 : CALL init_spherical_harmonics(lmax, -1)
841 :
842 63490 : DO iset = 1, gto_basis_set%nset
843 139577 : DO ishell = 1, gto_basis_set%nshell(iset)
844 76087 : l = gto_basis_set%l(ishell, iset)
845 76087 : first_cgf = gto_basis_set%first_cgf(ishell, iset)
846 76087 : first_sgf = gto_basis_set%first_sgf(ishell, iset)
847 76087 : ncgf = nco(l)
848 76087 : nsgf = nso(l)
849 : CALL dgemm("N", "T", n, nsgf, ncgf, &
850 : 1.0_dp, gto_basis_set%cphi(1, first_cgf), n, &
851 : orbtramat(l)%c2s(1, 1), nsgf, &
852 122614 : 0.0_dp, gto_basis_set%sphi(1, first_sgf), n)
853 : END DO
854 : END DO
855 : END IF
856 :
857 : ! Build the reduced transformation matrix "scon"
858 : ! This matrix transforms from cartesian primitifs to spherical contracted functions
859 : ! "scon" only includes primitifs (lmin -> lmax), whereas sphi is (0 -> lmax)
860 :
861 16973 : n = SIZE(gto_basis_set%scon, 1)
862 :
863 5681317 : gto_basis_set%scon = 0.0_dp
864 16973 : IF (n .GT. 0) THEN
865 63490 : DO iset = 1, gto_basis_set%nset
866 46527 : lmin = gto_basis_set%lmin(iset)
867 46527 : lmax = gto_basis_set%lmax(iset)
868 46527 : npgf = gto_basis_set%npgf(iset)
869 46527 : nn = ncoset(lmax) - ncoset(lmin - 1)
870 139577 : DO ishell = 1, gto_basis_set%nshell(iset)
871 76087 : first_sgf = gto_basis_set%first_sgf(ishell, iset)
872 76087 : last_sgf = gto_basis_set%last_sgf(ishell, iset)
873 391257 : DO ipgf = 1, npgf
874 268643 : nn1 = (ipgf - 1)*ncoset(lmax) + ncoset(lmin - 1) + 1
875 268643 : nn2 = ipgf*ncoset(lmax)
876 268643 : n1 = (ipgf - 1)*nn + 1
877 268643 : n2 = ipgf*nn
878 4771365 : gto_basis_set%scon(n1:n2, first_sgf:last_sgf) = gto_basis_set%sphi(nn1:nn2, first_sgf:last_sgf)
879 : END DO
880 : END DO
881 : END DO
882 : END IF
883 :
884 16973 : CALL timestop(handle)
885 :
886 16973 : END SUBROUTINE init_cphi_and_sphi
887 :
888 : ! **************************************************************************************************
889 : !> \brief ...
890 : !> \param gto_basis_set ...
891 : ! **************************************************************************************************
892 0 : SUBROUTINE init_norm_cgf_aux(gto_basis_set)
893 :
894 : ! Initialise the normalization factors of the contracted Cartesian Gaussian
895 : ! functions, if the Gaussian functions represent charge distributions.
896 :
897 : ! - Creation (07.12.2000,MK)
898 :
899 : TYPE(gto_basis_set_type), INTENT(INOUT) :: gto_basis_set
900 :
901 : INTEGER :: icgf, ico, ipgf, iset, ishell, jco, &
902 : jpgf, ll, lmax, lmin, lx, ly, lz, n, &
903 : npgfa
904 : REAL(KIND=dp) :: fnorm, gcca, gccb
905 0 : REAL(KIND=dp), ALLOCATABLE, DIMENSION(:) :: ff
906 : REAL(KIND=dp), ALLOCATABLE, DIMENSION(:, :) :: gaa
907 : REAL(KIND=dp), ALLOCATABLE, DIMENSION(:, :, :) :: vv
908 0 : REAL(KIND=dp), DIMENSION(:), POINTER :: rpgfa, zeta
909 :
910 : ! -------------------------------------------------------------------------
911 :
912 0 : n = 0
913 0 : ll = 0
914 0 : DO iset = 1, gto_basis_set%nset
915 0 : n = MAX(n, gto_basis_set%npgf(iset)*ncoset(gto_basis_set%lmax(iset)))
916 0 : ll = MAX(ll, gto_basis_set%lmax(iset))
917 : END DO
918 :
919 0 : ALLOCATE (gaa(n, n))
920 0 : ALLOCATE (vv(ncoset(ll), ncoset(ll), ll + ll + 1))
921 0 : ALLOCATE (ff(0:ll + ll))
922 :
923 0 : DO iset = 1, gto_basis_set%nset
924 0 : lmax = gto_basis_set%lmax(iset)
925 0 : lmin = gto_basis_set%lmin(iset)
926 0 : n = ncoset(lmax)
927 0 : npgfa = gto_basis_set%npgf(iset)
928 0 : rpgfa => gto_basis_set%pgf_radius(1:npgfa, iset)
929 0 : zeta => gto_basis_set%zet(1:npgfa, iset)
930 : CALL coulomb2(lmax, npgfa, zeta, rpgfa, lmin, &
931 : lmax, npgfa, zeta, rpgfa, lmin, &
932 0 : (/0.0_dp, 0.0_dp, 0.0_dp/), 0.0_dp, gaa, vv, ff(0:))
933 0 : DO ishell = 1, gto_basis_set%nshell(iset)
934 0 : DO icgf = gto_basis_set%first_cgf(ishell, iset), &
935 0 : gto_basis_set%last_cgf(ishell, iset)
936 0 : lx = gto_basis_set%lx(icgf)
937 0 : ly = gto_basis_set%ly(icgf)
938 0 : lz = gto_basis_set%lz(icgf)
939 0 : ico = coset(lx, ly, lz)
940 0 : fnorm = 0.0_dp
941 0 : DO ipgf = 1, npgfa
942 0 : gcca = gto_basis_set%gcc(ipgf, ishell, iset)
943 0 : jco = coset(lx, ly, lz)
944 0 : DO jpgf = 1, npgfa
945 0 : gccb = gto_basis_set%gcc(jpgf, ishell, iset)
946 0 : fnorm = fnorm + gcca*gccb*gaa(ico, jco)
947 0 : jco = jco + n
948 : END DO
949 0 : ico = ico + n
950 : END DO
951 0 : gto_basis_set%norm_cgf(icgf) = 1.0_dp/SQRT(fnorm)
952 : END DO
953 : END DO
954 : END DO
955 :
956 0 : DEALLOCATE (vv, ff)
957 0 : DEALLOCATE (gaa)
958 :
959 0 : END SUBROUTINE init_norm_cgf_aux
960 :
961 : ! **************************************************************************************************
962 : !> \brief ...
963 : !> \param gto_basis_set ...
964 : ! **************************************************************************************************
965 6 : ELEMENTAL SUBROUTINE init_norm_cgf_aux_2(gto_basis_set)
966 :
967 : ! Initialise the normalization factors of the auxiliary Cartesian Gaussian
968 : ! functions (Kim-Gordon polarization basis) Norm = 1.
969 :
970 : ! - Creation (07.12.2000,GT)
971 :
972 : TYPE(gto_basis_set_type), INTENT(INOUT) :: gto_basis_set
973 :
974 : INTEGER :: icgf, iset, ishell
975 :
976 92 : DO iset = 1, gto_basis_set%nset
977 178 : DO ishell = 1, gto_basis_set%nshell(iset)
978 396 : DO icgf = gto_basis_set%first_cgf(ishell, iset), &
979 172 : gto_basis_set%last_cgf(ishell, iset)
980 396 : gto_basis_set%norm_cgf(icgf) = 1.0_dp
981 : END DO
982 : END DO
983 : END DO
984 :
985 6 : END SUBROUTINE init_norm_cgf_aux_2
986 :
987 : ! **************************************************************************************************
988 : !> \brief Initialise the normalization factors of the contracted Cartesian Gaussian functions.
989 : !> \param gto_basis_set ...
990 : !> \author MK
991 : ! **************************************************************************************************
992 15349 : ELEMENTAL SUBROUTINE init_norm_cgf_orb(gto_basis_set)
993 :
994 : TYPE(gto_basis_set_type), INTENT(INOUT) :: gto_basis_set
995 :
996 : INTEGER :: icgf, ipgf, iset, ishell, jpgf, l, lx, &
997 : ly, lz
998 : REAL(KIND=dp) :: expzet, fnorm, gcca, gccb, prefac, zeta, &
999 : zetb
1000 :
1001 56766 : DO iset = 1, gto_basis_set%nset
1002 124799 : DO ishell = 1, gto_basis_set%nshell(iset)
1003 :
1004 68033 : l = gto_basis_set%l(ishell, iset)
1005 :
1006 68033 : expzet = 0.5_dp*REAL(2*l + 3, dp)
1007 :
1008 68033 : fnorm = 0.0_dp
1009 :
1010 321060 : DO ipgf = 1, gto_basis_set%npgf(iset)
1011 253027 : gcca = gto_basis_set%gcc(ipgf, ishell, iset)
1012 253027 : zeta = gto_basis_set%zet(ipgf, iset)
1013 1835141 : DO jpgf = 1, gto_basis_set%npgf(iset)
1014 1514081 : gccb = gto_basis_set%gcc(jpgf, ishell, iset)
1015 1514081 : zetb = gto_basis_set%zet(jpgf, iset)
1016 1767108 : fnorm = fnorm + gcca*gccb/(zeta + zetb)**expzet
1017 : END DO
1018 : END DO
1019 :
1020 68033 : fnorm = 0.5_dp**l*pi**1.5_dp*fnorm
1021 :
1022 255105 : DO icgf = gto_basis_set%first_cgf(ishell, iset), &
1023 109450 : gto_basis_set%last_cgf(ishell, iset)
1024 187072 : lx = gto_basis_set%lx(icgf)
1025 187072 : ly = gto_basis_set%ly(icgf)
1026 187072 : lz = gto_basis_set%lz(icgf)
1027 187072 : prefac = dfac(2*lx - 1)*dfac(2*ly - 1)*dfac(2*lz - 1)
1028 255105 : gto_basis_set%norm_cgf(icgf) = 1.0_dp/SQRT(prefac*fnorm)
1029 : END DO
1030 :
1031 : END DO
1032 : END DO
1033 :
1034 15349 : END SUBROUTINE init_norm_cgf_orb
1035 :
1036 : ! **************************************************************************************************
1037 : !> \brief Initialise the normalization factors of the contracted Cartesian Gaussian
1038 : !> functions used for frozen density representation.
1039 : !> \param gto_basis_set ...
1040 : !> \author GT
1041 : ! **************************************************************************************************
1042 0 : ELEMENTAL SUBROUTINE init_norm_cgf_orb_den(gto_basis_set)
1043 :
1044 : TYPE(gto_basis_set_type), INTENT(INOUT) :: gto_basis_set
1045 :
1046 : INTEGER :: icgf, ipgf, iset, ishell, l
1047 : REAL(KIND=dp) :: expzet, gcca, prefac, zeta
1048 :
1049 0 : DO iset = 1, gto_basis_set%nset
1050 0 : DO ishell = 1, gto_basis_set%nshell(iset)
1051 0 : l = gto_basis_set%l(ishell, iset)
1052 0 : expzet = 0.5_dp*REAL(2*l + 3, dp)
1053 0 : prefac = (1.0_dp/pi)**1.5_dp
1054 0 : DO ipgf = 1, gto_basis_set%npgf(iset)
1055 0 : gcca = gto_basis_set%gcc(ipgf, ishell, iset)
1056 0 : zeta = gto_basis_set%zet(ipgf, iset)
1057 0 : gto_basis_set%gcc(ipgf, ishell, iset) = prefac*zeta**expzet*gcca
1058 : END DO
1059 0 : DO icgf = gto_basis_set%first_cgf(ishell, iset), &
1060 0 : gto_basis_set%last_cgf(ishell, iset)
1061 0 : gto_basis_set%norm_cgf(icgf) = 1.0_dp
1062 : END DO
1063 : END DO
1064 : END DO
1065 :
1066 0 : END SUBROUTINE init_norm_cgf_orb_den
1067 :
1068 : ! **************************************************************************************************
1069 : !> \brief Initialise a Gaussian-type orbital (GTO) basis set data set.
1070 : !> \param gto_basis_set ...
1071 : !> \author MK
1072 : ! **************************************************************************************************
1073 30698 : SUBROUTINE init_orb_basis_set(gto_basis_set)
1074 :
1075 : TYPE(gto_basis_set_type), POINTER :: gto_basis_set
1076 :
1077 : CHARACTER(len=*), PARAMETER :: routineN = 'init_orb_basis_set'
1078 :
1079 : INTEGER :: handle
1080 :
1081 : ! -------------------------------------------------------------------------
1082 :
1083 15349 : IF (.NOT. ASSOCIATED(gto_basis_set)) RETURN
1084 :
1085 15349 : CALL timeset(routineN, handle)
1086 :
1087 15349 : SELECT CASE (gto_basis_set%norm_type)
1088 : CASE (0)
1089 : ! No normalisation requested
1090 : CASE (1)
1091 0 : CALL init_norm_cgf_orb_den(gto_basis_set)
1092 : CASE (2)
1093 : ! Normalise the primitive Gaussian functions
1094 15349 : CALL normalise_gcc_orb(gto_basis_set)
1095 : ! Compute the normalization factors of the contracted Gaussian-type
1096 : ! functions
1097 15349 : CALL init_norm_cgf_orb(gto_basis_set)
1098 : CASE DEFAULT
1099 15349 : CPABORT("Normalization method not specified")
1100 : END SELECT
1101 :
1102 : ! Initialise the transformation matrices "pgf" -> "cgf"
1103 :
1104 15349 : CALL init_cphi_and_sphi(gto_basis_set)
1105 :
1106 15349 : CALL timestop(handle)
1107 :
1108 : END SUBROUTINE init_orb_basis_set
1109 :
1110 : ! **************************************************************************************************
1111 : !> \brief Normalise the primitive Cartesian Gaussian functions. The normalization
1112 : !> factor is included in the Gaussian contraction coefficients.
1113 : !> \param gto_basis_set ...
1114 : !> \author MK
1115 : ! **************************************************************************************************
1116 15349 : SUBROUTINE normalise_gcc_orb(gto_basis_set)
1117 :
1118 : TYPE(gto_basis_set_type), POINTER :: gto_basis_set
1119 :
1120 : INTEGER :: ipgf, iset, ishell, l
1121 : REAL(KIND=dp) :: expzet, gcca, prefac, zeta
1122 :
1123 56766 : DO iset = 1, gto_basis_set%nset
1124 124799 : DO ishell = 1, gto_basis_set%nshell(iset)
1125 68033 : l = gto_basis_set%l(ishell, iset)
1126 68033 : expzet = 0.25_dp*REAL(2*l + 3, dp)
1127 68033 : prefac = 2.0_dp**l*(2.0_dp/pi)**0.75_dp
1128 362477 : DO ipgf = 1, gto_basis_set%npgf(iset)
1129 253027 : gcca = gto_basis_set%gcc(ipgf, ishell, iset)
1130 253027 : zeta = gto_basis_set%zet(ipgf, iset)
1131 321060 : gto_basis_set%gcc(ipgf, ishell, iset) = prefac*zeta**expzet*gcca
1132 : END DO
1133 : END DO
1134 : END DO
1135 :
1136 15349 : END SUBROUTINE normalise_gcc_orb
1137 :
1138 : ! **************************************************************************************************
1139 : !> \brief Read a Gaussian-type orbital (GTO) basis set from the database file.
1140 : !> \param element_symbol ...
1141 : !> \param basis_set_name ...
1142 : !> \param gto_basis_set ...
1143 : !> \param para_env ...
1144 : !> \param dft_section ...
1145 : !> \author MK
1146 : ! **************************************************************************************************
1147 11282 : SUBROUTINE read_gto_basis_set1(element_symbol, basis_set_name, gto_basis_set, &
1148 : para_env, dft_section)
1149 :
1150 : CHARACTER(LEN=*), INTENT(IN) :: element_symbol, basis_set_name
1151 : TYPE(gto_basis_set_type), INTENT(INOUT) :: gto_basis_set
1152 : TYPE(mp_para_env_type), POINTER :: para_env
1153 : TYPE(section_vals_type), POINTER :: dft_section
1154 :
1155 : CHARACTER(LEN=240) :: line
1156 : CHARACTER(LEN=242) :: line2
1157 : CHARACTER(len=default_path_length) :: basis_set_file_name, tmp
1158 : CHARACTER(LEN=default_path_length), DIMENSION(:), &
1159 11282 : POINTER :: cbasis
1160 11282 : CHARACTER(LEN=LEN(basis_set_name)) :: bsname
1161 11282 : CHARACTER(LEN=LEN(basis_set_name)+2) :: bsname2
1162 11282 : CHARACTER(LEN=LEN(element_symbol)) :: symbol
1163 11282 : CHARACTER(LEN=LEN(element_symbol)+2) :: symbol2
1164 : INTEGER :: i, ibasis, ico, ipgf, irep, iset, ishell, lshell, m, maxco, maxl, maxpgf, &
1165 : maxshell, nbasis, ncgf, nmin, nset, nsgf, sort_method, strlen1, strlen2
1166 11282 : INTEGER, DIMENSION(:), POINTER :: lmax, lmin, npgf, nshell
1167 11282 : INTEGER, DIMENSION(:, :), POINTER :: l, n
1168 : LOGICAL :: basis_found, found, match
1169 11282 : REAL(KIND=dp), DIMENSION(:, :), POINTER :: zet
1170 11282 : REAL(KIND=dp), DIMENSION(:, :, :), POINTER :: gcc
1171 : TYPE(cp_parser_type) :: parser
1172 :
1173 11282 : line = ""
1174 11282 : line2 = ""
1175 11282 : symbol = ""
1176 11282 : symbol2 = ""
1177 11282 : bsname = ""
1178 11282 : bsname2 = ""
1179 :
1180 11282 : nbasis = 1
1181 :
1182 11282 : gto_basis_set%name = basis_set_name
1183 11282 : gto_basis_set%aliases = basis_set_name
1184 : CALL section_vals_val_get(dft_section, "BASIS_SET_FILE_NAME", &
1185 11282 : n_rep_val=nbasis)
1186 33846 : ALLOCATE (cbasis(nbasis))
1187 24480 : DO ibasis = 1, nbasis
1188 : CALL section_vals_val_get(dft_section, "BASIS_SET_FILE_NAME", &
1189 13198 : i_rep_val=ibasis, c_val=cbasis(ibasis))
1190 13198 : basis_set_file_name = cbasis(ibasis)
1191 13198 : tmp = basis_set_file_name
1192 13198 : CALL uppercase(tmp)
1193 24480 : IF (INDEX(tmp, "MOLOPT") .NE. 0) CALL cite_reference(VandeVondele2007)
1194 : END DO
1195 :
1196 : ! Search for the requested basis set in the basis set file
1197 : ! until the basis set is found or the end of file is reached
1198 :
1199 11282 : basis_found = .FALSE.
1200 23450 : basis_loop: DO ibasis = 1, nbasis
1201 13126 : IF (basis_found) EXIT basis_loop
1202 12168 : basis_set_file_name = cbasis(ibasis)
1203 12168 : CALL parser_create(parser, basis_set_file_name, para_env=para_env)
1204 :
1205 12168 : bsname = basis_set_name
1206 12168 : symbol = element_symbol
1207 12168 : irep = 0
1208 :
1209 12168 : tmp = basis_set_name
1210 12168 : CALL uppercase(tmp)
1211 12168 : IF (INDEX(tmp, "MOLOPT") .NE. 0) CALL cite_reference(VandeVondele2007)
1212 :
1213 12168 : nset = 0
1214 12168 : maxshell = 0
1215 12168 : maxpgf = 0
1216 12168 : maxco = 0
1217 12168 : ncgf = 0
1218 12168 : nsgf = 0
1219 12168 : gto_basis_set%nset = nset
1220 12168 : gto_basis_set%ncgf = ncgf
1221 12168 : gto_basis_set%nsgf = nsgf
1222 12168 : CALL reallocate(gto_basis_set%lmax, 1, nset)
1223 12168 : CALL reallocate(gto_basis_set%lmin, 1, nset)
1224 12168 : CALL reallocate(gto_basis_set%npgf, 1, nset)
1225 12168 : CALL reallocate(gto_basis_set%nshell, 1, nset)
1226 12168 : CALL reallocate(gto_basis_set%n, 1, maxshell, 1, nset)
1227 12168 : CALL reallocate(gto_basis_set%l, 1, maxshell, 1, nset)
1228 12168 : CALL reallocate(gto_basis_set%zet, 1, maxpgf, 1, nset)
1229 12168 : CALL reallocate(gto_basis_set%gcc, 1, maxpgf, 1, maxshell, 1, nset)
1230 12168 : CALL reallocate(gto_basis_set%set_radius, 1, nset)
1231 12168 : CALL reallocate(gto_basis_set%pgf_radius, 1, maxpgf, 1, nset)
1232 12168 : CALL reallocate(gto_basis_set%first_cgf, 1, maxshell, 1, nset)
1233 12168 : CALL reallocate(gto_basis_set%first_sgf, 1, maxshell, 1, nset)
1234 12168 : CALL reallocate(gto_basis_set%last_cgf, 1, maxshell, 1, nset)
1235 12168 : CALL reallocate(gto_basis_set%last_sgf, 1, maxshell, 1, nset)
1236 12168 : CALL reallocate(gto_basis_set%ncgf_set, 1, nset)
1237 12168 : CALL reallocate(gto_basis_set%nsgf_set, 1, nset)
1238 12168 : CALL reallocate(gto_basis_set%cphi, 1, maxco, 1, ncgf)
1239 12168 : CALL reallocate(gto_basis_set%sphi, 1, maxco, 1, nsgf)
1240 12168 : CALL reallocate(gto_basis_set%scon, 1, maxco, 1, nsgf)
1241 12168 : CALL reallocate(gto_basis_set%lx, 1, ncgf)
1242 12168 : CALL reallocate(gto_basis_set%ly, 1, ncgf)
1243 12168 : CALL reallocate(gto_basis_set%lz, 1, ncgf)
1244 12168 : CALL reallocate(gto_basis_set%m, 1, nsgf)
1245 12168 : CALL reallocate(gto_basis_set%norm_cgf, 1, ncgf)
1246 :
1247 12168 : IF (tmp .NE. "NONE") THEN
1248 : search_loop: DO
1249 :
1250 59979 : CALL parser_search_string(parser, TRIM(bsname), .TRUE., found, line)
1251 59979 : IF (found) THEN
1252 59093 : CALL uppercase(symbol)
1253 59093 : CALL uppercase(bsname)
1254 :
1255 59093 : match = .FALSE.
1256 59093 : CALL uppercase(line)
1257 : ! Check both the element symbol and the basis set name
1258 59093 : line2 = " "//line//" "
1259 59093 : symbol2 = " "//TRIM(symbol)//" "
1260 59093 : bsname2 = " "//TRIM(bsname)//" "
1261 59093 : strlen1 = LEN_TRIM(symbol2) + 1
1262 59093 : strlen2 = LEN_TRIM(bsname2) + 1
1263 :
1264 59093 : IF ((INDEX(line2, symbol2(:strlen1)) > 0) .AND. &
1265 11274 : (INDEX(line2, bsname2(:strlen2)) > 0)) match = .TRUE.
1266 59093 : IF (match) THEN
1267 : ! copy all names into aliases field
1268 11274 : i = INDEX(line2, symbol2(:strlen1))
1269 11274 : i = i + 1 + INDEX(line2(i + 1:), " ")
1270 11274 : gto_basis_set%aliases = line2(i:)
1271 :
1272 11274 : NULLIFY (gcc, l, lmax, lmin, n, npgf, nshell, zet)
1273 : ! Read the basis set information
1274 11274 : CALL parser_get_object(parser, nset, newline=.TRUE.)
1275 :
1276 11274 : CALL reallocate(npgf, 1, nset)
1277 11274 : CALL reallocate(nshell, 1, nset)
1278 11274 : CALL reallocate(lmax, 1, nset)
1279 11274 : CALL reallocate(lmin, 1, nset)
1280 11274 : CALL reallocate(n, 1, 1, 1, nset)
1281 :
1282 11274 : maxl = 0
1283 11274 : maxpgf = 0
1284 11274 : maxshell = 0
1285 :
1286 44258 : DO iset = 1, nset
1287 32984 : CALL parser_get_object(parser, n(1, iset), newline=.TRUE.)
1288 32984 : CALL parser_get_object(parser, lmin(iset))
1289 32984 : CALL parser_get_object(parser, lmax(iset))
1290 32984 : CALL parser_get_object(parser, npgf(iset))
1291 32984 : maxl = MAX(maxl, lmax(iset))
1292 32984 : IF (npgf(iset) > maxpgf) THEN
1293 11458 : maxpgf = npgf(iset)
1294 11458 : CALL reallocate(zet, 1, maxpgf, 1, nset)
1295 11458 : CALL reallocate(gcc, 1, maxpgf, 1, maxshell, 1, nset)
1296 : END IF
1297 32984 : nshell(iset) = 0
1298 74866 : DO lshell = lmin(iset), lmax(iset)
1299 41882 : nmin = n(1, iset) + lshell - lmin(iset)
1300 41882 : CALL parser_get_object(parser, ishell)
1301 41882 : nshell(iset) = nshell(iset) + ishell
1302 41882 : IF (nshell(iset) > maxshell) THEN
1303 17556 : maxshell = nshell(iset)
1304 17556 : CALL reallocate(n, 1, maxshell, 1, nset)
1305 17556 : CALL reallocate(l, 1, maxshell, 1, nset)
1306 17556 : CALL reallocate(gcc, 1, maxpgf, 1, maxshell, 1, nset)
1307 : END IF
1308 169065 : DO i = 1, ishell
1309 52317 : n(nshell(iset) - ishell + i, iset) = nmin + i - 1
1310 94199 : l(nshell(iset) - ishell + i, iset) = lshell
1311 : END DO
1312 : END DO
1313 112950 : DO ipgf = 1, npgf(iset)
1314 68692 : CALL parser_get_object(parser, zet(ipgf, iset), newline=.TRUE.)
1315 250001 : DO ishell = 1, nshell(iset)
1316 217017 : CALL parser_get_object(parser, gcc(ipgf, ishell, iset))
1317 : END DO
1318 : END DO
1319 : END DO
1320 :
1321 : ! Maximum angular momentum quantum number of the atomic kind
1322 :
1323 11274 : CALL init_orbital_pointers(maxl)
1324 :
1325 : ! Allocate the global variables
1326 :
1327 11274 : gto_basis_set%nset = nset
1328 11274 : CALL reallocate(gto_basis_set%lmax, 1, nset)
1329 11274 : CALL reallocate(gto_basis_set%lmin, 1, nset)
1330 11274 : CALL reallocate(gto_basis_set%npgf, 1, nset)
1331 11274 : CALL reallocate(gto_basis_set%nshell, 1, nset)
1332 11274 : CALL reallocate(gto_basis_set%n, 1, maxshell, 1, nset)
1333 11274 : CALL reallocate(gto_basis_set%l, 1, maxshell, 1, nset)
1334 11274 : CALL reallocate(gto_basis_set%zet, 1, maxpgf, 1, nset)
1335 11274 : CALL reallocate(gto_basis_set%gcc, 1, maxpgf, 1, maxshell, 1, nset)
1336 :
1337 : ! Copy the basis set information into the data structure
1338 :
1339 44258 : DO iset = 1, nset
1340 32984 : gto_basis_set%lmax(iset) = lmax(iset)
1341 32984 : gto_basis_set%lmin(iset) = lmin(iset)
1342 32984 : gto_basis_set%npgf(iset) = npgf(iset)
1343 32984 : gto_basis_set%nshell(iset) = nshell(iset)
1344 85301 : DO ishell = 1, nshell(iset)
1345 52317 : gto_basis_set%n(ishell, iset) = n(ishell, iset)
1346 52317 : gto_basis_set%l(ishell, iset) = l(ishell, iset)
1347 233626 : DO ipgf = 1, npgf(iset)
1348 200642 : gto_basis_set%gcc(ipgf, ishell, iset) = gcc(ipgf, ishell, iset)
1349 : END DO
1350 : END DO
1351 112950 : DO ipgf = 1, npgf(iset)
1352 101676 : gto_basis_set%zet(ipgf, iset) = zet(ipgf, iset)
1353 : END DO
1354 : END DO
1355 :
1356 : ! Initialise the depending atomic kind information
1357 :
1358 11274 : CALL reallocate(gto_basis_set%set_radius, 1, nset)
1359 11274 : CALL reallocate(gto_basis_set%pgf_radius, 1, maxpgf, 1, nset)
1360 11274 : CALL reallocate(gto_basis_set%first_cgf, 1, maxshell, 1, nset)
1361 11274 : CALL reallocate(gto_basis_set%first_sgf, 1, maxshell, 1, nset)
1362 11274 : CALL reallocate(gto_basis_set%last_cgf, 1, maxshell, 1, nset)
1363 11274 : CALL reallocate(gto_basis_set%last_sgf, 1, maxshell, 1, nset)
1364 11274 : CALL reallocate(gto_basis_set%ncgf_set, 1, nset)
1365 11274 : CALL reallocate(gto_basis_set%nsgf_set, 1, nset)
1366 :
1367 11274 : maxco = 0
1368 11274 : ncgf = 0
1369 11274 : nsgf = 0
1370 :
1371 44258 : DO iset = 1, nset
1372 32984 : gto_basis_set%ncgf_set(iset) = 0
1373 32984 : gto_basis_set%nsgf_set(iset) = 0
1374 85301 : DO ishell = 1, nshell(iset)
1375 52317 : lshell = gto_basis_set%l(ishell, iset)
1376 52317 : gto_basis_set%first_cgf(ishell, iset) = ncgf + 1
1377 52317 : ncgf = ncgf + nco(lshell)
1378 52317 : gto_basis_set%last_cgf(ishell, iset) = ncgf
1379 : gto_basis_set%ncgf_set(iset) = &
1380 52317 : gto_basis_set%ncgf_set(iset) + nco(lshell)
1381 52317 : gto_basis_set%first_sgf(ishell, iset) = nsgf + 1
1382 52317 : nsgf = nsgf + nso(lshell)
1383 52317 : gto_basis_set%last_sgf(ishell, iset) = nsgf
1384 : gto_basis_set%nsgf_set(iset) = &
1385 85301 : gto_basis_set%nsgf_set(iset) + nso(lshell)
1386 : END DO
1387 44258 : maxco = MAX(maxco, npgf(iset)*ncoset(lmax(iset)))
1388 : END DO
1389 :
1390 11274 : gto_basis_set%ncgf = ncgf
1391 11274 : gto_basis_set%nsgf = nsgf
1392 :
1393 11274 : CALL reallocate(gto_basis_set%cphi, 1, maxco, 1, ncgf)
1394 11274 : CALL reallocate(gto_basis_set%sphi, 1, maxco, 1, nsgf)
1395 11274 : CALL reallocate(gto_basis_set%scon, 1, maxco, 1, nsgf)
1396 11274 : CALL reallocate(gto_basis_set%lx, 1, ncgf)
1397 11274 : CALL reallocate(gto_basis_set%ly, 1, ncgf)
1398 11274 : CALL reallocate(gto_basis_set%lz, 1, ncgf)
1399 11274 : CALL reallocate(gto_basis_set%m, 1, nsgf)
1400 11274 : CALL reallocate(gto_basis_set%norm_cgf, 1, ncgf)
1401 33822 : ALLOCATE (gto_basis_set%cgf_symbol(ncgf))
1402 :
1403 33822 : ALLOCATE (gto_basis_set%sgf_symbol(nsgf))
1404 :
1405 11274 : ncgf = 0
1406 11274 : nsgf = 0
1407 :
1408 44258 : DO iset = 1, nset
1409 96575 : DO ishell = 1, nshell(iset)
1410 52317 : lshell = gto_basis_set%l(ishell, iset)
1411 195377 : DO ico = ncoset(lshell - 1) + 1, ncoset(lshell)
1412 143060 : ncgf = ncgf + 1
1413 143060 : gto_basis_set%lx(ncgf) = indco(1, ico)
1414 143060 : gto_basis_set%ly(ncgf) = indco(2, ico)
1415 143060 : gto_basis_set%lz(ncgf) = indco(3, ico)
1416 : gto_basis_set%cgf_symbol(ncgf) = &
1417 : cgf_symbol(n(ishell, iset), (/gto_basis_set%lx(ncgf), &
1418 : gto_basis_set%ly(ncgf), &
1419 624557 : gto_basis_set%lz(ncgf)/))
1420 : END DO
1421 214496 : DO m = -lshell, lshell
1422 129195 : nsgf = nsgf + 1
1423 129195 : gto_basis_set%m(nsgf) = m
1424 : gto_basis_set%sgf_symbol(nsgf) = &
1425 181512 : sgf_symbol(n(ishell, iset), lshell, m)
1426 : END DO
1427 : END DO
1428 : END DO
1429 :
1430 11274 : DEALLOCATE (gcc, l, lmax, lmin, n, npgf, nshell, zet)
1431 :
1432 11274 : basis_found = .TRUE.
1433 11274 : EXIT search_loop
1434 : END IF
1435 : ELSE
1436 : EXIT search_loop
1437 : END IF
1438 : END DO search_loop
1439 : ELSE
1440 8 : match = .FALSE.
1441 16 : ALLOCATE (gto_basis_set%cgf_symbol(ncgf))
1442 16 : ALLOCATE (gto_basis_set%sgf_symbol(nsgf))
1443 : END IF
1444 :
1445 35618 : CALL parser_release(parser)
1446 :
1447 : END DO basis_loop
1448 :
1449 11282 : IF (tmp .NE. "NONE") THEN
1450 11274 : IF (.NOT. basis_found) THEN
1451 0 : basis_set_file_name = ""
1452 0 : DO ibasis = 1, nbasis
1453 0 : basis_set_file_name = TRIM(basis_set_file_name)//"<"//TRIM(cbasis(ibasis))//"> "
1454 : END DO
1455 : CALL cp_abort(__LOCATION__, &
1456 : "The requested basis set <"//TRIM(bsname)// &
1457 : "> for element <"//TRIM(symbol)//"> was not "// &
1458 : "found in the basis set files "// &
1459 0 : TRIM(basis_set_file_name))
1460 : END IF
1461 : END IF
1462 11282 : DEALLOCATE (cbasis)
1463 :
1464 11282 : CALL section_vals_val_get(dft_section, "SORT_BASIS", i_val=sort_method)
1465 11282 : CALL sort_gto_basis_set(gto_basis_set, sort_method)
1466 :
1467 56410 : END SUBROUTINE read_gto_basis_set1
1468 :
1469 : ! **************************************************************************************************
1470 : !> \brief Read a Gaussian-type orbital (GTO) basis set from the database file.
1471 : !> \param element_symbol ...
1472 : !> \param basis_type ...
1473 : !> \param gto_basis_set ...
1474 : !> \param basis_section ...
1475 : !> \param irep ...
1476 : !> \param dft_section ...
1477 : !> \author MK
1478 : ! **************************************************************************************************
1479 174 : SUBROUTINE read_gto_basis_set2(element_symbol, basis_type, gto_basis_set, &
1480 : basis_section, irep, dft_section)
1481 :
1482 : CHARACTER(LEN=*), INTENT(IN) :: element_symbol
1483 : CHARACTER(LEN=*), INTENT(INOUT) :: basis_type
1484 : TYPE(gto_basis_set_type), INTENT(INOUT) :: gto_basis_set
1485 : TYPE(section_vals_type), OPTIONAL, POINTER :: basis_section
1486 : INTEGER :: irep
1487 : TYPE(section_vals_type), OPTIONAL, POINTER :: dft_section
1488 :
1489 : CHARACTER(len=20*default_string_length) :: line_att
1490 : CHARACTER(LEN=240) :: line
1491 : CHARACTER(LEN=242) :: line2
1492 : CHARACTER(LEN=default_path_length) :: bsname, bsname2
1493 174 : CHARACTER(LEN=LEN(element_symbol)) :: symbol
1494 174 : CHARACTER(LEN=LEN(element_symbol)+2) :: symbol2
1495 : INTEGER :: i, ico, ipgf, iset, ishell, lshell, m, &
1496 : maxco, maxl, maxpgf, maxshell, ncgf, &
1497 : nmin, nset, nsgf, sort_method
1498 174 : INTEGER, DIMENSION(:), POINTER :: lmax, lmin, npgf, nshell
1499 174 : INTEGER, DIMENSION(:, :), POINTER :: l, n
1500 : LOGICAL :: is_ok
1501 174 : REAL(KIND=dp), DIMENSION(:, :), POINTER :: zet
1502 174 : REAL(KIND=dp), DIMENSION(:, :, :), POINTER :: gcc
1503 : TYPE(cp_sll_val_type), POINTER :: list
1504 : TYPE(val_type), POINTER :: val
1505 :
1506 174 : line = ""
1507 174 : line2 = ""
1508 174 : symbol = ""
1509 174 : symbol2 = ""
1510 : bsname = ""
1511 174 : bsname2 = ""
1512 :
1513 174 : gto_basis_set%name = " "
1514 174 : gto_basis_set%aliases = " "
1515 :
1516 174 : bsname = " "
1517 174 : symbol = element_symbol
1518 :
1519 174 : nset = 0
1520 174 : maxshell = 0
1521 174 : maxpgf = 0
1522 174 : maxco = 0
1523 174 : ncgf = 0
1524 174 : nsgf = 0
1525 174 : gto_basis_set%nset = nset
1526 174 : gto_basis_set%ncgf = ncgf
1527 174 : gto_basis_set%nsgf = nsgf
1528 174 : CALL reallocate(gto_basis_set%lmax, 1, nset)
1529 174 : CALL reallocate(gto_basis_set%lmin, 1, nset)
1530 174 : CALL reallocate(gto_basis_set%npgf, 1, nset)
1531 174 : CALL reallocate(gto_basis_set%nshell, 1, nset)
1532 174 : CALL reallocate(gto_basis_set%n, 1, maxshell, 1, nset)
1533 174 : CALL reallocate(gto_basis_set%l, 1, maxshell, 1, nset)
1534 174 : CALL reallocate(gto_basis_set%zet, 1, maxpgf, 1, nset)
1535 174 : CALL reallocate(gto_basis_set%gcc, 1, maxpgf, 1, maxshell, 1, nset)
1536 174 : CALL reallocate(gto_basis_set%set_radius, 1, nset)
1537 174 : CALL reallocate(gto_basis_set%pgf_radius, 1, maxpgf, 1, nset)
1538 174 : CALL reallocate(gto_basis_set%first_cgf, 1, maxshell, 1, nset)
1539 174 : CALL reallocate(gto_basis_set%first_sgf, 1, maxshell, 1, nset)
1540 174 : CALL reallocate(gto_basis_set%last_cgf, 1, maxshell, 1, nset)
1541 174 : CALL reallocate(gto_basis_set%last_sgf, 1, maxshell, 1, nset)
1542 174 : CALL reallocate(gto_basis_set%ncgf_set, 1, nset)
1543 174 : CALL reallocate(gto_basis_set%nsgf_set, 1, nset)
1544 174 : CALL reallocate(gto_basis_set%cphi, 1, maxco, 1, ncgf)
1545 174 : CALL reallocate(gto_basis_set%sphi, 1, maxco, 1, nsgf)
1546 174 : CALL reallocate(gto_basis_set%scon, 1, maxco, 1, nsgf)
1547 174 : CALL reallocate(gto_basis_set%lx, 1, ncgf)
1548 174 : CALL reallocate(gto_basis_set%ly, 1, ncgf)
1549 174 : CALL reallocate(gto_basis_set%lz, 1, ncgf)
1550 174 : CALL reallocate(gto_basis_set%m, 1, nsgf)
1551 174 : CALL reallocate(gto_basis_set%norm_cgf, 1, ncgf)
1552 :
1553 174 : basis_type = ""
1554 174 : CALL section_vals_val_get(basis_section, "_SECTION_PARAMETERS_", i_rep_section=irep, c_val=basis_type)
1555 174 : IF (basis_type == "Orbital") basis_type = "ORB"
1556 :
1557 174 : NULLIFY (list, val)
1558 174 : CALL section_vals_list_get(basis_section, "_DEFAULT_KEYWORD_", i_rep_section=irep, list=list)
1559 174 : CALL uppercase(symbol)
1560 174 : CALL uppercase(bsname)
1561 :
1562 174 : NULLIFY (gcc, l, lmax, lmin, n, npgf, nshell, zet)
1563 : ! Read the basis set information
1564 174 : is_ok = cp_sll_val_next(list, val)
1565 174 : IF (.NOT. is_ok) CPABORT("Error reading the Basis set from input file!")
1566 174 : CALL val_get(val, c_val=line_att)
1567 174 : READ (line_att, *) nset
1568 :
1569 174 : CALL reallocate(npgf, 1, nset)
1570 174 : CALL reallocate(nshell, 1, nset)
1571 174 : CALL reallocate(lmax, 1, nset)
1572 174 : CALL reallocate(lmin, 1, nset)
1573 174 : CALL reallocate(n, 1, 1, 1, nset)
1574 :
1575 174 : maxl = 0
1576 174 : maxpgf = 0
1577 174 : maxshell = 0
1578 :
1579 500 : DO iset = 1, nset
1580 326 : is_ok = cp_sll_val_next(list, val)
1581 326 : IF (.NOT. is_ok) CPABORT("Error reading the Basis set from input file!")
1582 326 : CALL val_get(val, c_val=line_att)
1583 326 : READ (line_att, *) n(1, iset)
1584 326 : CALL remove_word(line_att)
1585 326 : READ (line_att, *) lmin(iset)
1586 326 : CALL remove_word(line_att)
1587 326 : READ (line_att, *) lmax(iset)
1588 326 : CALL remove_word(line_att)
1589 326 : READ (line_att, *) npgf(iset)
1590 326 : CALL remove_word(line_att)
1591 326 : maxl = MAX(maxl, lmax(iset))
1592 326 : IF (npgf(iset) > maxpgf) THEN
1593 174 : maxpgf = npgf(iset)
1594 174 : CALL reallocate(zet, 1, maxpgf, 1, nset)
1595 174 : CALL reallocate(gcc, 1, maxpgf, 1, maxshell, 1, nset)
1596 : END IF
1597 326 : nshell(iset) = 0
1598 1082 : DO lshell = lmin(iset), lmax(iset)
1599 756 : nmin = n(1, iset) + lshell - lmin(iset)
1600 756 : READ (line_att, *) ishell
1601 756 : CALL remove_word(line_att)
1602 756 : nshell(iset) = nshell(iset) + ishell
1603 756 : IF (nshell(iset) > maxshell) THEN
1604 334 : maxshell = nshell(iset)
1605 334 : CALL reallocate(n, 1, maxshell, 1, nset)
1606 334 : CALL reallocate(l, 1, maxshell, 1, nset)
1607 334 : CALL reallocate(gcc, 1, maxpgf, 1, maxshell, 1, nset)
1608 : END IF
1609 2126 : DO i = 1, ishell
1610 1044 : n(nshell(iset) - ishell + i, iset) = nmin + i - 1
1611 1800 : l(nshell(iset) - ishell + i, iset) = lshell
1612 : END DO
1613 : END DO
1614 326 : IF (LEN_TRIM(line_att) /= 0) &
1615 0 : CPABORT("Error reading the Basis from input file!")
1616 1308 : DO ipgf = 1, npgf(iset)
1617 808 : is_ok = cp_sll_val_next(list, val)
1618 808 : IF (.NOT. is_ok) CPABORT("Error reading the Basis set from input file!")
1619 808 : CALL val_get(val, c_val=line_att)
1620 1134 : READ (line_att, *) zet(ipgf, iset), (gcc(ipgf, ishell, iset), ishell=1, nshell(iset))
1621 : END DO
1622 : END DO
1623 :
1624 : ! Maximum angular momentum quantum number of the atomic kind
1625 :
1626 174 : CALL init_orbital_pointers(maxl)
1627 :
1628 : ! Allocate the global variables
1629 :
1630 174 : gto_basis_set%nset = nset
1631 174 : CALL reallocate(gto_basis_set%lmax, 1, nset)
1632 174 : CALL reallocate(gto_basis_set%lmin, 1, nset)
1633 174 : CALL reallocate(gto_basis_set%npgf, 1, nset)
1634 174 : CALL reallocate(gto_basis_set%nshell, 1, nset)
1635 174 : CALL reallocate(gto_basis_set%n, 1, maxshell, 1, nset)
1636 174 : CALL reallocate(gto_basis_set%l, 1, maxshell, 1, nset)
1637 174 : CALL reallocate(gto_basis_set%zet, 1, maxpgf, 1, nset)
1638 174 : CALL reallocate(gto_basis_set%gcc, 1, maxpgf, 1, maxshell, 1, nset)
1639 :
1640 : ! Copy the basis set information into the data structure
1641 :
1642 500 : DO iset = 1, nset
1643 326 : gto_basis_set%lmax(iset) = lmax(iset)
1644 326 : gto_basis_set%lmin(iset) = lmin(iset)
1645 326 : gto_basis_set%npgf(iset) = npgf(iset)
1646 326 : gto_basis_set%nshell(iset) = nshell(iset)
1647 1370 : DO ishell = 1, nshell(iset)
1648 1044 : gto_basis_set%n(ishell, iset) = n(ishell, iset)
1649 1044 : gto_basis_set%l(ishell, iset) = l(ishell, iset)
1650 5656 : DO ipgf = 1, npgf(iset)
1651 5330 : gto_basis_set%gcc(ipgf, ishell, iset) = gcc(ipgf, ishell, iset)
1652 : END DO
1653 : END DO
1654 1308 : DO ipgf = 1, npgf(iset)
1655 1134 : gto_basis_set%zet(ipgf, iset) = zet(ipgf, iset)
1656 : END DO
1657 : END DO
1658 :
1659 : ! Initialise the depending atomic kind information
1660 :
1661 174 : CALL reallocate(gto_basis_set%set_radius, 1, nset)
1662 174 : CALL reallocate(gto_basis_set%pgf_radius, 1, maxpgf, 1, nset)
1663 174 : CALL reallocate(gto_basis_set%first_cgf, 1, maxshell, 1, nset)
1664 174 : CALL reallocate(gto_basis_set%first_sgf, 1, maxshell, 1, nset)
1665 174 : CALL reallocate(gto_basis_set%last_cgf, 1, maxshell, 1, nset)
1666 174 : CALL reallocate(gto_basis_set%last_sgf, 1, maxshell, 1, nset)
1667 174 : CALL reallocate(gto_basis_set%ncgf_set, 1, nset)
1668 174 : CALL reallocate(gto_basis_set%nsgf_set, 1, nset)
1669 :
1670 174 : maxco = 0
1671 174 : ncgf = 0
1672 174 : nsgf = 0
1673 :
1674 500 : DO iset = 1, nset
1675 326 : gto_basis_set%ncgf_set(iset) = 0
1676 326 : gto_basis_set%nsgf_set(iset) = 0
1677 1370 : DO ishell = 1, nshell(iset)
1678 1044 : lshell = gto_basis_set%l(ishell, iset)
1679 1044 : gto_basis_set%first_cgf(ishell, iset) = ncgf + 1
1680 1044 : ncgf = ncgf + nco(lshell)
1681 1044 : gto_basis_set%last_cgf(ishell, iset) = ncgf
1682 : gto_basis_set%ncgf_set(iset) = &
1683 1044 : gto_basis_set%ncgf_set(iset) + nco(lshell)
1684 1044 : gto_basis_set%first_sgf(ishell, iset) = nsgf + 1
1685 1044 : nsgf = nsgf + nso(lshell)
1686 1044 : gto_basis_set%last_sgf(ishell, iset) = nsgf
1687 : gto_basis_set%nsgf_set(iset) = &
1688 1370 : gto_basis_set%nsgf_set(iset) + nso(lshell)
1689 : END DO
1690 500 : maxco = MAX(maxco, npgf(iset)*ncoset(lmax(iset)))
1691 : END DO
1692 :
1693 174 : gto_basis_set%ncgf = ncgf
1694 174 : gto_basis_set%nsgf = nsgf
1695 :
1696 174 : CALL reallocate(gto_basis_set%cphi, 1, maxco, 1, ncgf)
1697 174 : CALL reallocate(gto_basis_set%sphi, 1, maxco, 1, nsgf)
1698 174 : CALL reallocate(gto_basis_set%scon, 1, maxco, 1, nsgf)
1699 174 : CALL reallocate(gto_basis_set%lx, 1, ncgf)
1700 174 : CALL reallocate(gto_basis_set%ly, 1, ncgf)
1701 174 : CALL reallocate(gto_basis_set%lz, 1, ncgf)
1702 174 : CALL reallocate(gto_basis_set%m, 1, nsgf)
1703 174 : CALL reallocate(gto_basis_set%norm_cgf, 1, ncgf)
1704 522 : ALLOCATE (gto_basis_set%cgf_symbol(ncgf))
1705 :
1706 522 : ALLOCATE (gto_basis_set%sgf_symbol(nsgf))
1707 :
1708 174 : ncgf = 0
1709 174 : nsgf = 0
1710 :
1711 500 : DO iset = 1, nset
1712 1544 : DO ishell = 1, nshell(iset)
1713 1044 : lshell = gto_basis_set%l(ishell, iset)
1714 5040 : DO ico = ncoset(lshell - 1) + 1, ncoset(lshell)
1715 3996 : ncgf = ncgf + 1
1716 3996 : gto_basis_set%lx(ncgf) = indco(1, ico)
1717 3996 : gto_basis_set%ly(ncgf) = indco(2, ico)
1718 3996 : gto_basis_set%lz(ncgf) = indco(3, ico)
1719 : gto_basis_set%cgf_symbol(ncgf) = &
1720 : cgf_symbol(n(ishell, iset), (/gto_basis_set%lx(ncgf), &
1721 : gto_basis_set%ly(ncgf), &
1722 17028 : gto_basis_set%lz(ncgf)/))
1723 : END DO
1724 4606 : DO m = -lshell, lshell
1725 3236 : nsgf = nsgf + 1
1726 3236 : gto_basis_set%m(nsgf) = m
1727 : gto_basis_set%sgf_symbol(nsgf) = &
1728 4280 : sgf_symbol(n(ishell, iset), lshell, m)
1729 : END DO
1730 : END DO
1731 : END DO
1732 :
1733 174 : DEALLOCATE (gcc, l, lmax, lmin, n, npgf, nshell, zet)
1734 :
1735 174 : IF (PRESENT(dft_section)) THEN
1736 162 : CALL section_vals_val_get(dft_section, "SORT_BASIS", i_val=sort_method)
1737 162 : CALL sort_gto_basis_set(gto_basis_set, sort_method)
1738 : END IF
1739 :
1740 174 : END SUBROUTINE read_gto_basis_set2
1741 :
1742 : ! **************************************************************************************************
1743 : !> \brief Set the components of Gaussian-type orbital (GTO) basis set data set.
1744 : !> \param gto_basis_set ...
1745 : !> \param name ...
1746 : !> \param aliases ...
1747 : !> \param norm_type ...
1748 : !> \param kind_radius ...
1749 : !> \param ncgf ...
1750 : !> \param nset ...
1751 : !> \param nsgf ...
1752 : !> \param cgf_symbol ...
1753 : !> \param sgf_symbol ...
1754 : !> \param norm_cgf ...
1755 : !> \param set_radius ...
1756 : !> \param lmax ...
1757 : !> \param lmin ...
1758 : !> \param lx ...
1759 : !> \param ly ...
1760 : !> \param lz ...
1761 : !> \param m ...
1762 : !> \param ncgf_set ...
1763 : !> \param npgf ...
1764 : !> \param nsgf_set ...
1765 : !> \param nshell ...
1766 : !> \param cphi ...
1767 : !> \param pgf_radius ...
1768 : !> \param sphi ...
1769 : !> \param scon ...
1770 : !> \param zet ...
1771 : !> \param first_cgf ...
1772 : !> \param first_sgf ...
1773 : !> \param l ...
1774 : !> \param last_cgf ...
1775 : !> \param last_sgf ...
1776 : !> \param n ...
1777 : !> \param gcc ...
1778 : !> \param short_kind_radius ...
1779 : !> \author MK
1780 : ! **************************************************************************************************
1781 33070 : SUBROUTINE set_gto_basis_set(gto_basis_set, name, aliases, norm_type, kind_radius, ncgf, &
1782 : nset, nsgf, cgf_symbol, sgf_symbol, norm_cgf, set_radius, lmax, &
1783 : lmin, lx, ly, lz, m, ncgf_set, npgf, nsgf_set, nshell, &
1784 : cphi, pgf_radius, sphi, scon, zet, first_cgf, first_sgf, l, &
1785 : last_cgf, last_sgf, n, gcc, short_kind_radius)
1786 :
1787 : TYPE(gto_basis_set_type), INTENT(INOUT) :: gto_basis_set
1788 : CHARACTER(LEN=default_string_length), INTENT(IN), &
1789 : OPTIONAL :: name, aliases
1790 : INTEGER, INTENT(IN), OPTIONAL :: norm_type
1791 : REAL(KIND=dp), INTENT(IN), OPTIONAL :: kind_radius
1792 : INTEGER, INTENT(IN), OPTIONAL :: ncgf, nset, nsgf
1793 : CHARACTER(LEN=12), DIMENSION(:), OPTIONAL, POINTER :: cgf_symbol
1794 : CHARACTER(LEN=6), DIMENSION(:), OPTIONAL, POINTER :: sgf_symbol
1795 : REAL(KIND=dp), DIMENSION(:), OPTIONAL, POINTER :: norm_cgf, set_radius
1796 : INTEGER, DIMENSION(:), OPTIONAL, POINTER :: lmax, lmin, lx, ly, lz, m, ncgf_set, &
1797 : npgf, nsgf_set, nshell
1798 : REAL(KIND=dp), DIMENSION(:, :), OPTIONAL, POINTER :: cphi, pgf_radius, sphi, scon, zet
1799 : INTEGER, DIMENSION(:, :), OPTIONAL, POINTER :: first_cgf, first_sgf, l, last_cgf, &
1800 : last_sgf, n
1801 : REAL(KIND=dp), DIMENSION(:, :, :), OPTIONAL, &
1802 : POINTER :: gcc
1803 : REAL(KIND=dp), INTENT(IN), OPTIONAL :: short_kind_radius
1804 :
1805 33070 : IF (PRESENT(name)) gto_basis_set%name = name
1806 33070 : IF (PRESENT(aliases)) gto_basis_set%aliases = aliases
1807 33070 : IF (PRESENT(norm_type)) gto_basis_set%norm_type = norm_type
1808 33070 : IF (PRESENT(kind_radius)) gto_basis_set%kind_radius = kind_radius
1809 33070 : IF (PRESENT(short_kind_radius)) gto_basis_set%short_kind_radius = short_kind_radius
1810 33070 : IF (PRESENT(ncgf)) gto_basis_set%ncgf = ncgf
1811 33070 : IF (PRESENT(nset)) gto_basis_set%nset = nset
1812 33070 : IF (PRESENT(nsgf)) gto_basis_set%nsgf = nsgf
1813 33070 : IF (PRESENT(cgf_symbol)) gto_basis_set%cgf_symbol(:) = cgf_symbol(:)
1814 33070 : IF (PRESENT(sgf_symbol)) gto_basis_set%sgf_symbol(:) = sgf_symbol(:)
1815 33070 : IF (PRESENT(norm_cgf)) gto_basis_set%norm_cgf(:) = norm_cgf(:)
1816 149512 : IF (PRESENT(set_radius)) gto_basis_set%set_radius(:) = set_radius(:)
1817 33070 : IF (PRESENT(lmax)) gto_basis_set%lmax(:) = lmax(:)
1818 33070 : IF (PRESENT(lmin)) gto_basis_set%lmin(:) = lmin(:)
1819 33070 : IF (PRESENT(lx)) gto_basis_set%lx(:) = lx(:)
1820 33070 : IF (PRESENT(ly)) gto_basis_set%ly(:) = ly(:)
1821 33070 : IF (PRESENT(lz)) gto_basis_set%lz(:) = lz(:)
1822 33070 : IF (PRESENT(m)) gto_basis_set%m(:) = m(:)
1823 33070 : IF (PRESENT(ncgf_set)) gto_basis_set%ncgf_set(:) = ncgf_set(:)
1824 33070 : IF (PRESENT(npgf)) gto_basis_set%npgf(:) = npgf(:)
1825 33070 : IF (PRESENT(nsgf_set)) gto_basis_set%nsgf_set(:) = nsgf_set(:)
1826 33070 : IF (PRESENT(nshell)) gto_basis_set%nshell(:) = nshell(:)
1827 33070 : IF (PRESENT(cphi)) gto_basis_set%cphi(:, :) = cphi(:, :)
1828 449714 : IF (PRESENT(pgf_radius)) gto_basis_set%pgf_radius(:, :) = pgf_radius(:, :)
1829 33070 : IF (PRESENT(sphi)) gto_basis_set%sphi(:, :) = sphi(:, :)
1830 33070 : IF (PRESENT(scon)) gto_basis_set%scon(:, :) = scon(:, :)
1831 33070 : IF (PRESENT(zet)) gto_basis_set%zet(:, :) = zet(:, :)
1832 33070 : IF (PRESENT(first_cgf)) gto_basis_set%first_cgf(:, :) = first_cgf(:, :)
1833 33070 : IF (PRESENT(first_sgf)) gto_basis_set%first_sgf(:, :) = first_sgf(:, :)
1834 33070 : IF (PRESENT(l)) l(:, :) = gto_basis_set%l(:, :)
1835 33070 : IF (PRESENT(last_cgf)) gto_basis_set%last_cgf(:, :) = last_cgf(:, :)
1836 33070 : IF (PRESENT(last_sgf)) gto_basis_set%last_sgf(:, :) = last_sgf(:, :)
1837 33070 : IF (PRESENT(n)) gto_basis_set%n(:, :) = n(:, :)
1838 33070 : IF (PRESENT(gcc)) gto_basis_set%gcc(:, :, :) = gcc(:, :, :)
1839 :
1840 33070 : END SUBROUTINE set_gto_basis_set
1841 :
1842 : ! **************************************************************************************************
1843 : !> \brief Write a Gaussian-type orbital (GTO) basis set data set to the output unit.
1844 : !> \param gto_basis_set ...
1845 : !> \param output_unit ...
1846 : !> \param header ...
1847 : !> \author MK
1848 : ! **************************************************************************************************
1849 139 : SUBROUTINE write_gto_basis_set(gto_basis_set, output_unit, header)
1850 :
1851 : TYPE(gto_basis_set_type), INTENT(IN) :: gto_basis_set
1852 : INTEGER, INTENT(in) :: output_unit
1853 : CHARACTER(len=*), OPTIONAL :: header
1854 :
1855 : INTEGER :: ipgf, iset, ishell
1856 :
1857 139 : IF (output_unit > 0) THEN
1858 :
1859 139 : IF (PRESENT(header)) THEN
1860 : WRITE (UNIT=output_unit, FMT="(/,T6,A,T41,A40)") &
1861 139 : TRIM(header), TRIM(gto_basis_set%name)
1862 : END IF
1863 :
1864 : WRITE (UNIT=output_unit, FMT="(/,(T8,A,T71,I10))") &
1865 139 : "Number of orbital shell sets: ", &
1866 139 : gto_basis_set%nset, &
1867 139 : "Number of orbital shells: ", &
1868 489 : SUM(gto_basis_set%nshell(:)), &
1869 139 : "Number of primitive Cartesian functions: ", &
1870 489 : SUM(gto_basis_set%npgf(:)), &
1871 139 : "Number of Cartesian basis functions: ", &
1872 139 : gto_basis_set%ncgf, &
1873 139 : "Number of spherical basis functions: ", &
1874 139 : gto_basis_set%nsgf, &
1875 139 : "Norm type: ", &
1876 278 : gto_basis_set%norm_type
1877 :
1878 : WRITE (UNIT=output_unit, FMT="(/,T6,A,T41,A40,/,/,T25,A)") &
1879 139 : "GTO basis set information for", TRIM(gto_basis_set%name), &
1880 278 : "Set Shell n l Exponent Coefficient"
1881 :
1882 489 : DO iset = 1, gto_basis_set%nset
1883 350 : WRITE (UNIT=output_unit, FMT="(A)") ""
1884 976 : DO ishell = 1, gto_basis_set%nshell(iset)
1885 : WRITE (UNIT=output_unit, &
1886 : FMT="(T25,I3,4X,I4,4X,I2,2X,I2,(T51,2F15.6))") &
1887 487 : iset, ishell, &
1888 487 : gto_basis_set%n(ishell, iset), &
1889 487 : gto_basis_set%l(ishell, iset), &
1890 1839 : (gto_basis_set%zet(ipgf, iset), &
1891 2326 : gto_basis_set%gcc(ipgf, ishell, iset), &
1892 3650 : ipgf=1, gto_basis_set%npgf(iset))
1893 : END DO
1894 : END DO
1895 :
1896 : END IF
1897 :
1898 139 : END SUBROUTINE write_gto_basis_set
1899 :
1900 : ! **************************************************************************************************
1901 :
1902 : ! **************************************************************************************************
1903 : !> \brief Write a Gaussian-type orbital (GTO) basis set data set to the output unit.
1904 : !> \param orb_basis_set ...
1905 : !> \param output_unit ...
1906 : !> \param header ...
1907 : !> \author MK
1908 : ! **************************************************************************************************
1909 4497 : SUBROUTINE write_orb_basis_set(orb_basis_set, output_unit, header)
1910 :
1911 : TYPE(gto_basis_set_type), INTENT(IN) :: orb_basis_set
1912 : INTEGER, INTENT(in) :: output_unit
1913 : CHARACTER(len=*), OPTIONAL :: header
1914 :
1915 : INTEGER :: icgf, ico, ipgf, iset, ishell
1916 :
1917 4497 : IF (output_unit > 0) THEN
1918 4497 : IF (PRESENT(header)) THEN
1919 : WRITE (UNIT=output_unit, FMT="(/,T6,A,T41,A40)") &
1920 4497 : TRIM(header), TRIM(orb_basis_set%name)
1921 : END IF
1922 :
1923 : WRITE (UNIT=output_unit, FMT="(/,(T8,A,T71,I10))") &
1924 4497 : "Number of orbital shell sets: ", &
1925 4497 : orb_basis_set%nset, &
1926 4497 : "Number of orbital shells: ", &
1927 17524 : SUM(orb_basis_set%nshell(:)), &
1928 4497 : "Number of primitive Cartesian functions: ", &
1929 17524 : SUM(orb_basis_set%npgf(:)), &
1930 4497 : "Number of Cartesian basis functions: ", &
1931 4497 : orb_basis_set%ncgf, &
1932 4497 : "Number of spherical basis functions: ", &
1933 4497 : orb_basis_set%nsgf, &
1934 4497 : "Norm type: ", &
1935 8994 : orb_basis_set%norm_type
1936 :
1937 : WRITE (UNIT=output_unit, FMT="(/,T8,A,/,/,T25,A)") &
1938 4497 : "Normalised Cartesian orbitals:", &
1939 8994 : "Set Shell Orbital Exponent Coefficient"
1940 :
1941 4497 : icgf = 0
1942 :
1943 17524 : DO iset = 1, orb_basis_set%nset
1944 36798 : DO ishell = 1, orb_basis_set%nshell(iset)
1945 19274 : WRITE (UNIT=output_unit, FMT="(A)") ""
1946 84913 : DO ico = 1, nco(orb_basis_set%l(ishell, iset))
1947 52612 : icgf = icgf + 1
1948 : WRITE (UNIT=output_unit, &
1949 : FMT="(T25,I3,4X,I4,3X,A12,(T51,2F15.6))") &
1950 52612 : iset, ishell, orb_basis_set%cgf_symbol(icgf), &
1951 160471 : (orb_basis_set%zet(ipgf, iset), &
1952 : orb_basis_set%norm_cgf(icgf)* &
1953 213083 : orb_basis_set%gcc(ipgf, ishell, iset), &
1954 337581 : ipgf=1, orb_basis_set%npgf(iset))
1955 : END DO
1956 : END DO
1957 : END DO
1958 : END IF
1959 :
1960 4497 : END SUBROUTINE write_orb_basis_set
1961 :
1962 : ! **************************************************************************************************
1963 : !> \brief ...
1964 : !> \param sto_basis_set ...
1965 : ! **************************************************************************************************
1966 3148 : SUBROUTINE allocate_sto_basis_set(sto_basis_set)
1967 :
1968 : TYPE(sto_basis_set_type), POINTER :: sto_basis_set
1969 :
1970 3148 : CALL deallocate_sto_basis_set(sto_basis_set)
1971 :
1972 3148 : ALLOCATE (sto_basis_set)
1973 :
1974 3148 : END SUBROUTINE allocate_sto_basis_set
1975 :
1976 : ! **************************************************************************************************
1977 : !> \brief ...
1978 : !> \param sto_basis_set ...
1979 : ! **************************************************************************************************
1980 8020 : SUBROUTINE deallocate_sto_basis_set(sto_basis_set)
1981 :
1982 : TYPE(sto_basis_set_type), POINTER :: sto_basis_set
1983 :
1984 : ! -------------------------------------------------------------------------
1985 :
1986 8020 : IF (ASSOCIATED(sto_basis_set)) THEN
1987 3148 : IF (ASSOCIATED(sto_basis_set%symbol)) THEN
1988 3008 : DEALLOCATE (sto_basis_set%symbol)
1989 : END IF
1990 3148 : IF (ASSOCIATED(sto_basis_set%nq)) THEN
1991 3148 : DEALLOCATE (sto_basis_set%nq)
1992 : END IF
1993 3148 : IF (ASSOCIATED(sto_basis_set%lq)) THEN
1994 3148 : DEALLOCATE (sto_basis_set%lq)
1995 : END IF
1996 3148 : IF (ASSOCIATED(sto_basis_set%zet)) THEN
1997 3148 : DEALLOCATE (sto_basis_set%zet)
1998 : END IF
1999 :
2000 3148 : DEALLOCATE (sto_basis_set)
2001 : END IF
2002 8020 : END SUBROUTINE deallocate_sto_basis_set
2003 :
2004 : ! **************************************************************************************************
2005 : !> \brief ...
2006 : !> \param sto_basis_set ...
2007 : !> \param name ...
2008 : !> \param nshell ...
2009 : !> \param symbol ...
2010 : !> \param nq ...
2011 : !> \param lq ...
2012 : !> \param zet ...
2013 : !> \param maxlq ...
2014 : !> \param numsto ...
2015 : ! **************************************************************************************************
2016 3148 : SUBROUTINE get_sto_basis_set(sto_basis_set, name, nshell, symbol, nq, lq, zet, maxlq, numsto)
2017 :
2018 : TYPE(sto_basis_set_type), INTENT(IN) :: sto_basis_set
2019 : CHARACTER(LEN=default_string_length), &
2020 : INTENT(OUT), OPTIONAL :: name
2021 : INTEGER, INTENT(OUT), OPTIONAL :: nshell
2022 : CHARACTER(LEN=6), DIMENSION(:), OPTIONAL, POINTER :: symbol
2023 : INTEGER, DIMENSION(:), OPTIONAL, POINTER :: nq, lq
2024 : REAL(KIND=dp), DIMENSION(:), OPTIONAL, POINTER :: zet
2025 : INTEGER, INTENT(OUT), OPTIONAL :: maxlq, numsto
2026 :
2027 : INTEGER :: iset
2028 :
2029 3148 : IF (PRESENT(name)) name = sto_basis_set%name
2030 3148 : IF (PRESENT(nshell)) nshell = sto_basis_set%nshell
2031 3148 : IF (PRESENT(symbol)) symbol => sto_basis_set%symbol
2032 3148 : IF (PRESENT(nq)) nq => sto_basis_set%nq
2033 3148 : IF (PRESENT(lq)) lq => sto_basis_set%lq
2034 3148 : IF (PRESENT(zet)) zet => sto_basis_set%zet
2035 3148 : IF (PRESENT(maxlq)) THEN
2036 0 : maxlq = MAXVAL(sto_basis_set%lq(1:sto_basis_set%nshell))
2037 : END IF
2038 3148 : IF (PRESENT(numsto)) THEN
2039 0 : numsto = 0
2040 0 : DO iset = 1, sto_basis_set%nshell
2041 0 : numsto = numsto + 2*sto_basis_set%lq(iset) + 1
2042 : END DO
2043 : END IF
2044 :
2045 3148 : END SUBROUTINE get_sto_basis_set
2046 :
2047 : ! **************************************************************************************************
2048 : !> \brief ...
2049 : !> \param sto_basis_set ...
2050 : !> \param name ...
2051 : !> \param nshell ...
2052 : !> \param symbol ...
2053 : !> \param nq ...
2054 : !> \param lq ...
2055 : !> \param zet ...
2056 : ! **************************************************************************************************
2057 3144 : SUBROUTINE set_sto_basis_set(sto_basis_set, name, nshell, symbol, nq, lq, zet)
2058 :
2059 : TYPE(sto_basis_set_type), INTENT(INOUT) :: sto_basis_set
2060 : CHARACTER(LEN=default_string_length), INTENT(IN), &
2061 : OPTIONAL :: name
2062 : INTEGER, INTENT(IN), OPTIONAL :: nshell
2063 : CHARACTER(LEN=6), DIMENSION(:), OPTIONAL, POINTER :: symbol
2064 : INTEGER, DIMENSION(:), OPTIONAL, POINTER :: nq, lq
2065 : REAL(KIND=dp), DIMENSION(:), OPTIONAL, POINTER :: zet
2066 :
2067 : INTEGER :: ns
2068 :
2069 3144 : IF (PRESENT(name)) sto_basis_set%name = name
2070 3144 : IF (PRESENT(nshell)) sto_basis_set%nshell = nshell
2071 3144 : IF (PRESENT(symbol)) THEN
2072 3004 : ns = SIZE(symbol)
2073 3004 : IF (ASSOCIATED(sto_basis_set%symbol)) DEALLOCATE (sto_basis_set%symbol)
2074 9012 : ALLOCATE (sto_basis_set%symbol(1:ns))
2075 13548 : sto_basis_set%symbol(:) = symbol(:)
2076 : END IF
2077 3144 : IF (PRESENT(nq)) THEN
2078 3144 : ns = SIZE(nq)
2079 3144 : CALL reallocate(sto_basis_set%nq, 1, ns)
2080 13828 : sto_basis_set%nq = nq(:)
2081 : END IF
2082 3144 : IF (PRESENT(lq)) THEN
2083 3144 : ns = SIZE(lq)
2084 3144 : CALL reallocate(sto_basis_set%lq, 1, ns)
2085 13828 : sto_basis_set%lq = lq(:)
2086 : END IF
2087 3144 : IF (PRESENT(zet)) THEN
2088 3144 : ns = SIZE(zet)
2089 3144 : CALL reallocate(sto_basis_set%zet, 1, ns)
2090 13828 : sto_basis_set%zet = zet(:)
2091 : END IF
2092 :
2093 3144 : END SUBROUTINE set_sto_basis_set
2094 :
2095 : ! **************************************************************************************************
2096 : !> \brief ...
2097 : !> \param element_symbol ...
2098 : !> \param basis_set_name ...
2099 : !> \param sto_basis_set ...
2100 : !> \param para_env ...
2101 : !> \param dft_section ...
2102 : ! **************************************************************************************************
2103 4 : SUBROUTINE read_sto_basis_set(element_symbol, basis_set_name, sto_basis_set, para_env, dft_section)
2104 :
2105 : ! Read a Slater-type orbital (STO) basis set from the database file.
2106 :
2107 : CHARACTER(LEN=*), INTENT(IN) :: element_symbol, basis_set_name
2108 : TYPE(sto_basis_set_type), INTENT(INOUT) :: sto_basis_set
2109 : TYPE(mp_para_env_type), POINTER :: para_env
2110 : TYPE(section_vals_type), POINTER :: dft_section
2111 :
2112 : CHARACTER(LEN=10) :: nlsym
2113 : CHARACTER(LEN=2) :: lsym
2114 : CHARACTER(LEN=240) :: line
2115 : CHARACTER(LEN=242) :: line2
2116 : CHARACTER(len=default_path_length) :: basis_set_file_name, tmp
2117 : CHARACTER(LEN=default_path_length), DIMENSION(:), &
2118 4 : POINTER :: cbasis
2119 4 : CHARACTER(LEN=LEN(basis_set_name)) :: bsname
2120 4 : CHARACTER(LEN=LEN(basis_set_name)+2) :: bsname2
2121 4 : CHARACTER(LEN=LEN(element_symbol)) :: symbol
2122 4 : CHARACTER(LEN=LEN(element_symbol)+2) :: symbol2
2123 : INTEGER :: ibasis, irep, iset, nbasis, nq, nset, &
2124 : strlen1, strlen2
2125 : LOGICAL :: basis_found, found, match
2126 : REAL(KIND=dp) :: zet
2127 : TYPE(cp_parser_type) :: parser
2128 :
2129 4 : line = ""
2130 4 : line2 = ""
2131 4 : symbol = ""
2132 4 : symbol2 = ""
2133 4 : bsname = ""
2134 4 : bsname2 = ""
2135 :
2136 4 : nbasis = 1
2137 :
2138 4 : sto_basis_set%name = basis_set_name
2139 : CALL section_vals_val_get(dft_section, "BASIS_SET_FILE_NAME", &
2140 4 : n_rep_val=nbasis)
2141 12 : ALLOCATE (cbasis(nbasis))
2142 8 : DO ibasis = 1, nbasis
2143 : CALL section_vals_val_get(dft_section, "BASIS_SET_FILE_NAME", &
2144 4 : i_rep_val=ibasis, c_val=cbasis(ibasis))
2145 4 : basis_set_file_name = cbasis(ibasis)
2146 4 : tmp = basis_set_file_name
2147 8 : CALL uppercase(tmp)
2148 : END DO
2149 :
2150 : ! Search for the requested basis set in the basis set file
2151 : ! until the basis set is found or the end of file is reached
2152 :
2153 4 : basis_found = .FALSE.
2154 8 : basis_loop: DO ibasis = 1, nbasis
2155 4 : IF (basis_found) EXIT basis_loop
2156 4 : basis_set_file_name = cbasis(ibasis)
2157 4 : CALL parser_create(parser, basis_set_file_name, para_env=para_env)
2158 :
2159 4 : bsname = basis_set_name
2160 4 : symbol = element_symbol
2161 4 : irep = 0
2162 :
2163 4 : tmp = basis_set_name
2164 4 : CALL uppercase(tmp)
2165 :
2166 4 : IF (tmp .NE. "NONE") THEN
2167 : search_loop: DO
2168 :
2169 6 : CALL parser_search_string(parser, TRIM(bsname), .TRUE., found, line)
2170 6 : IF (found) THEN
2171 6 : CALL uppercase(symbol)
2172 6 : CALL uppercase(bsname)
2173 :
2174 6 : match = .FALSE.
2175 6 : CALL uppercase(line)
2176 : ! Check both the element symbol and the basis set name
2177 6 : line2 = " "//line//" "
2178 6 : symbol2 = " "//TRIM(symbol)//" "
2179 6 : bsname2 = " "//TRIM(bsname)//" "
2180 6 : strlen1 = LEN_TRIM(symbol2) + 1
2181 6 : strlen2 = LEN_TRIM(bsname2) + 1
2182 :
2183 6 : IF ((INDEX(line2, symbol2(:strlen1)) > 0) .AND. &
2184 4 : (INDEX(line2, bsname2(:strlen2)) > 0)) match = .TRUE.
2185 6 : IF (match) THEN
2186 : ! Read the basis set information
2187 4 : CALL parser_get_object(parser, nset, newline=.TRUE.)
2188 4 : sto_basis_set%nshell = nset
2189 :
2190 4 : CALL reallocate(sto_basis_set%nq, 1, nset)
2191 4 : CALL reallocate(sto_basis_set%lq, 1, nset)
2192 4 : CALL reallocate(sto_basis_set%zet, 1, nset)
2193 12 : ALLOCATE (sto_basis_set%symbol(nset))
2194 :
2195 20 : DO iset = 1, nset
2196 16 : CALL parser_get_object(parser, nq, newline=.TRUE.)
2197 16 : CALL parser_get_object(parser, lsym)
2198 16 : CALL parser_get_object(parser, zet)
2199 16 : sto_basis_set%nq(iset) = nq
2200 16 : sto_basis_set%zet(iset) = zet
2201 16 : WRITE (nlsym, "(I2,A)") nq, TRIM(lsym)
2202 16 : sto_basis_set%symbol(iset) = TRIM(nlsym)
2203 20 : SELECT CASE (TRIM(lsym))
2204 : CASE ("S", "s")
2205 12 : sto_basis_set%lq(iset) = 0
2206 : CASE ("P", "p")
2207 4 : sto_basis_set%lq(iset) = 1
2208 : CASE ("D", "d")
2209 0 : sto_basis_set%lq(iset) = 2
2210 : CASE ("F", "f")
2211 0 : sto_basis_set%lq(iset) = 3
2212 : CASE ("G", "g")
2213 0 : sto_basis_set%lq(iset) = 4
2214 : CASE ("H", "h")
2215 0 : sto_basis_set%lq(iset) = 5
2216 : CASE ("I", "i", "J", "j")
2217 0 : sto_basis_set%lq(iset) = 6
2218 : CASE ("K", "k")
2219 0 : sto_basis_set%lq(iset) = 7
2220 : CASE ("L", "l")
2221 0 : sto_basis_set%lq(iset) = 8
2222 : CASE ("M", "m")
2223 0 : sto_basis_set%lq(iset) = 9
2224 : CASE DEFAULT
2225 : CALL cp_abort(__LOCATION__, &
2226 : "The requested basis set <"//TRIM(bsname)// &
2227 16 : "> for element <"//TRIM(symbol)//"> has an invalid component: ")
2228 : END SELECT
2229 : END DO
2230 :
2231 : basis_found = .TRUE.
2232 : EXIT search_loop
2233 : END IF
2234 : ELSE
2235 : EXIT search_loop
2236 : END IF
2237 : END DO search_loop
2238 : ELSE
2239 4 : match = .FALSE.
2240 : END IF
2241 :
2242 12 : CALL parser_release(parser)
2243 :
2244 : END DO basis_loop
2245 :
2246 4 : IF (tmp .NE. "NONE") THEN
2247 4 : IF (.NOT. basis_found) THEN
2248 0 : basis_set_file_name = ""
2249 0 : DO ibasis = 1, nbasis
2250 0 : basis_set_file_name = TRIM(basis_set_file_name)//"<"//TRIM(cbasis(ibasis))//"> "
2251 : END DO
2252 : CALL cp_abort(__LOCATION__, &
2253 : "The requested basis set <"//TRIM(bsname)// &
2254 : "> for element <"//TRIM(symbol)//"> was not "// &
2255 : "found in the basis set files "// &
2256 0 : TRIM(basis_set_file_name))
2257 : END IF
2258 : END IF
2259 4 : DEALLOCATE (cbasis)
2260 :
2261 16 : END SUBROUTINE read_sto_basis_set
2262 :
2263 : ! **************************************************************************************************
2264 : !> \brief ...
2265 : !> \param sto_basis_set ...
2266 : !> \param gto_basis_set ...
2267 : !> \param ngauss ...
2268 : !> \param ortho ...
2269 : ! **************************************************************************************************
2270 3148 : SUBROUTINE create_gto_from_sto_basis(sto_basis_set, gto_basis_set, ngauss, ortho)
2271 :
2272 : TYPE(sto_basis_set_type), INTENT(IN) :: sto_basis_set
2273 : TYPE(gto_basis_set_type), POINTER :: gto_basis_set
2274 : INTEGER, INTENT(IN), OPTIONAL :: ngauss
2275 : LOGICAL, INTENT(IN), OPTIONAL :: ortho
2276 :
2277 : INTEGER, PARAMETER :: maxng = 6
2278 :
2279 : CHARACTER(LEN=default_string_length) :: name, sng
2280 : INTEGER :: i1, i2, ico, ipgf, iset, jset, l, &
2281 : lshell, m, maxco, maxl, ncgf, ng, ngs, &
2282 : np, nset, nsgf, nshell
2283 : INTEGER, DIMENSION(0:10) :: mxf
2284 3148 : INTEGER, DIMENSION(:), POINTER :: lq, nq
2285 : LOGICAL :: do_ortho
2286 3148 : REAL(KIND=dp), ALLOCATABLE, DIMENSION(:, :) :: gal, zal, zll
2287 3148 : REAL(KIND=dp), DIMENSION(:), POINTER :: zet
2288 : REAL(KIND=dp), DIMENSION(maxng) :: gcc, zetg
2289 :
2290 3148 : ng = 6
2291 3148 : IF (PRESENT(ngauss)) ng = ngauss
2292 3148 : IF (ng > maxng) CPABORT("Too many Gaussian primitives requested")
2293 3148 : do_ortho = .FALSE.
2294 3148 : IF (PRESENT(ortho)) do_ortho = ortho
2295 :
2296 3148 : CALL allocate_gto_basis_set(gto_basis_set)
2297 :
2298 : CALL get_sto_basis_set(sto_basis_set, name=name, nshell=nshell, nq=nq, &
2299 3148 : lq=lq, zet=zet)
2300 :
2301 13848 : maxl = MAXVAL(lq)
2302 3148 : CALL init_orbital_pointers(maxl)
2303 :
2304 3148 : CALL integer_to_string(ng, sng)
2305 3148 : gto_basis_set%name = TRIM(name)//"_STO-"//TRIM(sng)//"G"
2306 :
2307 3148 : nset = nshell
2308 3148 : gto_basis_set%nset = nset
2309 3148 : CALL reallocate(gto_basis_set%lmax, 1, nset)
2310 3148 : CALL reallocate(gto_basis_set%lmin, 1, nset)
2311 3148 : CALL reallocate(gto_basis_set%npgf, 1, nset)
2312 3148 : CALL reallocate(gto_basis_set%nshell, 1, nset)
2313 3148 : CALL reallocate(gto_basis_set%n, 1, 1, 1, nset)
2314 3148 : CALL reallocate(gto_basis_set%l, 1, 1, 1, nset)
2315 3148 : CALL reallocate(gto_basis_set%zet, 1, ng, 1, nset)
2316 3148 : CALL reallocate(gto_basis_set%gcc, 1, ng, 1, 1, 1, nset)
2317 :
2318 9638 : DO iset = 1, nset
2319 6490 : CALL get_sto_ng(zet(iset), ng, nq(iset), lq(iset), zetg, gcc)
2320 6490 : gto_basis_set%lmax(iset) = lq(iset)
2321 6490 : gto_basis_set%lmin(iset) = lq(iset)
2322 6490 : gto_basis_set%npgf(iset) = ng
2323 6490 : gto_basis_set%nshell(iset) = 1
2324 6490 : gto_basis_set%n(1, iset) = lq(iset) + 1
2325 6490 : gto_basis_set%l(1, iset) = lq(iset)
2326 46744 : DO ipgf = 1, ng
2327 37106 : gto_basis_set%gcc(ipgf, 1, iset) = gcc(ipgf)
2328 43596 : gto_basis_set%zet(ipgf, iset) = zetg(ipgf)
2329 : END DO
2330 : END DO
2331 :
2332 3148 : IF (do_ortho) THEN
2333 690 : mxf = 0
2334 2152 : DO iset = 1, nset
2335 1462 : l = gto_basis_set%l(1, iset)
2336 2152 : mxf(l) = mxf(l) + 1
2337 : END DO
2338 8280 : m = MAXVAL(mxf)
2339 690 : IF (m > 1) THEN
2340 2286 : ALLOCATE (gal(ng, nset), zal(ng, nset), zll(m*ng, 0:maxl))
2341 762 : DO iset = 1, nset
2342 2540 : zal(1:ng, iset) = gto_basis_set%zet(1:ng, iset)
2343 2794 : gal(1:ng, iset) = gto_basis_set%gcc(1:ng, 1, iset)
2344 : END DO
2345 254 : CALL reallocate(gto_basis_set%zet, 1, m*ng, 1, nset)
2346 254 : CALL reallocate(gto_basis_set%gcc, 1, m*ng, 1, 1, 1, nset)
2347 762 : DO iset = 1, nset
2348 508 : l = gto_basis_set%l(1, iset)
2349 762 : gto_basis_set%npgf(iset) = ng*mxf(l)
2350 : END DO
2351 4826 : gto_basis_set%zet = 0.0_dp
2352 5334 : gto_basis_set%gcc = 0.0_dp
2353 2540 : zll = 0.0_dp
2354 254 : mxf = 0
2355 762 : DO iset = 1, nset
2356 508 : l = gto_basis_set%l(1, iset)
2357 508 : mxf(l) = mxf(l) + 1
2358 508 : i1 = mxf(l)*ng - ng + 1
2359 508 : i2 = mxf(l)*ng
2360 2540 : zll(i1:i2, l) = zal(1:ng, iset)
2361 2794 : gto_basis_set%gcc(i1:i2, 1, iset) = gal(1:ng, iset)
2362 : END DO
2363 762 : DO iset = 1, nset
2364 508 : l = gto_basis_set%l(1, iset)
2365 4826 : gto_basis_set%zet(:, iset) = zll(:, l)
2366 : END DO
2367 762 : DO iset = 1, nset
2368 508 : l = gto_basis_set%l(1, iset)
2369 1016 : DO jset = 1, iset - 1
2370 762 : IF (gto_basis_set%l(1, iset) == l) THEN
2371 254 : m = mxf(l)*ng
2372 : CALL orthofun(gto_basis_set%zet(1:m, iset), gto_basis_set%gcc(1:m, 1, iset), &
2373 254 : gto_basis_set%gcc(1:m, 1, jset), l)
2374 : END IF
2375 : END DO
2376 : END DO
2377 254 : DEALLOCATE (gal, zal, zll)
2378 : END IF
2379 : END IF
2380 :
2381 9638 : ngs = MAXVAL(gto_basis_set%npgf(1:nset))
2382 3148 : CALL reallocate(gto_basis_set%set_radius, 1, nset)
2383 3148 : CALL reallocate(gto_basis_set%pgf_radius, 1, ngs, 1, nset)
2384 3148 : CALL reallocate(gto_basis_set%first_cgf, 1, 1, 1, nset)
2385 3148 : CALL reallocate(gto_basis_set%first_sgf, 1, 1, 1, nset)
2386 3148 : CALL reallocate(gto_basis_set%last_cgf, 1, 1, 1, nset)
2387 3148 : CALL reallocate(gto_basis_set%last_sgf, 1, 1, 1, nset)
2388 3148 : CALL reallocate(gto_basis_set%ncgf_set, 1, nset)
2389 3148 : CALL reallocate(gto_basis_set%nsgf_set, 1, nset)
2390 :
2391 3148 : maxco = 0
2392 3148 : ncgf = 0
2393 3148 : nsgf = 0
2394 :
2395 9638 : DO iset = 1, nset
2396 6490 : gto_basis_set%ncgf_set(iset) = 0
2397 6490 : gto_basis_set%nsgf_set(iset) = 0
2398 6490 : lshell = gto_basis_set%l(1, iset)
2399 6490 : gto_basis_set%first_cgf(1, iset) = ncgf + 1
2400 6490 : ncgf = ncgf + nco(lshell)
2401 6490 : gto_basis_set%last_cgf(1, iset) = ncgf
2402 : gto_basis_set%ncgf_set(iset) = &
2403 6490 : gto_basis_set%ncgf_set(iset) + nco(lshell)
2404 6490 : gto_basis_set%first_sgf(1, iset) = nsgf + 1
2405 6490 : nsgf = nsgf + nso(lshell)
2406 6490 : gto_basis_set%last_sgf(1, iset) = nsgf
2407 : gto_basis_set%nsgf_set(iset) = &
2408 6490 : gto_basis_set%nsgf_set(iset) + nso(lshell)
2409 6490 : ngs = gto_basis_set%npgf(iset)
2410 9638 : maxco = MAX(maxco, ngs*ncoset(lshell))
2411 : END DO
2412 :
2413 3148 : gto_basis_set%ncgf = ncgf
2414 3148 : gto_basis_set%nsgf = nsgf
2415 :
2416 3148 : CALL reallocate(gto_basis_set%cphi, 1, maxco, 1, ncgf)
2417 3148 : CALL reallocate(gto_basis_set%sphi, 1, maxco, 1, nsgf)
2418 3148 : CALL reallocate(gto_basis_set%scon, 1, maxco, 1, nsgf)
2419 3148 : CALL reallocate(gto_basis_set%lx, 1, ncgf)
2420 3148 : CALL reallocate(gto_basis_set%ly, 1, ncgf)
2421 3148 : CALL reallocate(gto_basis_set%lz, 1, ncgf)
2422 3148 : CALL reallocate(gto_basis_set%m, 1, nsgf)
2423 3148 : CALL reallocate(gto_basis_set%norm_cgf, 1, ncgf)
2424 9444 : ALLOCATE (gto_basis_set%cgf_symbol(ncgf))
2425 9444 : ALLOCATE (gto_basis_set%sgf_symbol(nsgf))
2426 :
2427 3148 : ncgf = 0
2428 3148 : nsgf = 0
2429 :
2430 9638 : DO iset = 1, nset
2431 6490 : lshell = gto_basis_set%l(1, iset)
2432 6490 : np = lshell + 1
2433 21604 : DO ico = ncoset(lshell - 1) + 1, ncoset(lshell)
2434 15114 : ncgf = ncgf + 1
2435 15114 : gto_basis_set%lx(ncgf) = indco(1, ico)
2436 15114 : gto_basis_set%ly(ncgf) = indco(2, ico)
2437 15114 : gto_basis_set%lz(ncgf) = indco(3, ico)
2438 : gto_basis_set%cgf_symbol(ncgf) = &
2439 : cgf_symbol(np, (/gto_basis_set%lx(ncgf), &
2440 : gto_basis_set%ly(ncgf), &
2441 66946 : gto_basis_set%lz(ncgf)/))
2442 : END DO
2443 23932 : DO m = -lshell, lshell
2444 14294 : nsgf = nsgf + 1
2445 14294 : gto_basis_set%m(nsgf) = m
2446 20784 : gto_basis_set%sgf_symbol(nsgf) = sgf_symbol(np, lshell, m)
2447 : END DO
2448 : END DO
2449 :
2450 3148 : gto_basis_set%norm_type = -1
2451 :
2452 6296 : END SUBROUTINE create_gto_from_sto_basis
2453 :
2454 : ! **************************************************************************************************
2455 : !> \brief ...
2456 : !> \param zet ...
2457 : !> \param co ...
2458 : !> \param cr ...
2459 : !> \param l ...
2460 : ! **************************************************************************************************
2461 254 : SUBROUTINE orthofun(zet, co, cr, l)
2462 : REAL(KIND=dp), DIMENSION(:), INTENT(IN) :: zet
2463 : REAL(KIND=dp), DIMENSION(:), INTENT(INOUT) :: co, cr
2464 : INTEGER, INTENT(IN) :: l
2465 :
2466 : REAL(KIND=dp) :: ss
2467 :
2468 254 : CALL aovlp(l, zet, cr, cr, ss)
2469 2286 : cr(:) = cr(:)/SQRT(ss)
2470 254 : CALL aovlp(l, zet, co, cr, ss)
2471 2286 : co(:) = co(:) - ss*cr(:)
2472 254 : CALL aovlp(l, zet, co, co, ss)
2473 2286 : co(:) = co(:)/SQRT(ss)
2474 :
2475 254 : END SUBROUTINE orthofun
2476 :
2477 : ! **************************************************************************************************
2478 : !> \brief ...
2479 : !> \param l ...
2480 : !> \param zet ...
2481 : !> \param ca ...
2482 : !> \param cb ...
2483 : !> \param ss ...
2484 : ! **************************************************************************************************
2485 762 : SUBROUTINE aovlp(l, zet, ca, cb, ss)
2486 : INTEGER, INTENT(IN) :: l
2487 : REAL(KIND=dp), DIMENSION(:), INTENT(IN) :: zet, ca, cb
2488 : REAL(KIND=dp), INTENT(OUT) :: ss
2489 :
2490 : INTEGER :: i, j, m
2491 : REAL(KIND=dp) :: ab, ai, aj, s00, sss
2492 :
2493 : !
2494 : ! use init_norm_cgf_orb
2495 : !
2496 762 : m = SIZE(zet)
2497 762 : ss = 0.0_dp
2498 6858 : DO i = 1, m
2499 6096 : ai = (2.0_dp*zet(i)/pi)**0.75_dp
2500 55626 : DO j = 1, m
2501 48768 : aj = (2.0_dp*zet(j)/pi)**0.75_dp
2502 48768 : ab = 1._dp/(zet(i) + zet(j))
2503 48768 : s00 = ai*aj*(pi*ab)**1.50_dp
2504 48768 : IF (l == 0) THEN
2505 : sss = s00
2506 0 : ELSEIF (l == 1) THEN
2507 0 : sss = s00*ab*0.5_dp
2508 : ELSE
2509 0 : CPABORT("aovlp lvalue")
2510 : END IF
2511 54864 : ss = ss + sss*ca(i)*cb(j)
2512 : END DO
2513 : END DO
2514 :
2515 762 : END SUBROUTINE aovlp
2516 :
2517 : ! **************************************************************************************************
2518 : !> \brief ...
2519 : !> \param z ...
2520 : !> \param ne ...
2521 : !> \param n ...
2522 : !> \param l ...
2523 : !> \return ...
2524 : ! **************************************************************************************************
2525 21570 : PURE FUNCTION srules(z, ne, n, l)
2526 : ! Slater rules
2527 : INTEGER, INTENT(IN) :: z
2528 : INTEGER, DIMENSION(:, :), INTENT(IN) :: ne
2529 : INTEGER, INTENT(IN) :: n, l
2530 : REAL(dp) :: srules
2531 :
2532 : REAL(dp), DIMENSION(7), PARAMETER :: &
2533 : xns = (/1.0_dp, 2.0_dp, 3.0_dp, 3.7_dp, 4.0_dp, 4.2_dp, 4.4_dp/)
2534 :
2535 : INTEGER :: i, l1, l2, m, m1, m2, nn
2536 : REAL(dp) :: s
2537 :
2538 21570 : s = 0.0_dp
2539 : ! The complete shell
2540 21570 : l1 = MIN(l + 1, 4)
2541 21570 : nn = MIN(n, 7)
2542 : IF (l1 == 1) l2 = 2
2543 21570 : IF (l1 == 2) l2 = 1
2544 15225 : IF (l1 == 3) l2 = 4
2545 20145 : IF (l1 == 4) l2 = 3
2546 : ! Rule a) no contribution from shells further out
2547 : ! Rule b) 0.35 (1s 0.3) from each other electron in the same shell
2548 21570 : IF (n == 1) THEN
2549 5980 : m = ne(1, 1)
2550 5980 : s = s + 0.3_dp*REAL(m - 1, dp)
2551 : ELSE
2552 15590 : m = ne(l1, nn) + ne(l2, nn)
2553 15590 : s = s + 0.35_dp*REAL(m - 1, dp)
2554 : END IF
2555 : ! Rule c) if (s,p) shell 0.85 from each electron with n-1, and 1.0
2556 : ! from all electrons further in
2557 21570 : IF (l1 + l2 == 3) THEN
2558 19960 : IF (nn > 1) THEN
2559 13980 : m1 = ne(1, nn - 1) + ne(2, nn - 1) + ne(3, nn - 1) + ne(4, nn - 1)
2560 13980 : m2 = 0
2561 21902 : DO i = 1, nn - 2
2562 21902 : m2 = m2 + ne(1, i) + ne(2, i) + ne(3, i) + ne(4, I)
2563 : END DO
2564 13980 : s = s + 0.85_dp*REAL(m1, dp) + 1._dp*REAL(m2, dp)
2565 : END IF
2566 : ELSE
2567 : ! Rule d) if (d,f) shell 1.0 from each electron inside
2568 : m = 0
2569 6532 : DO i = 1, nn - 1
2570 6532 : m = m + ne(1, i) + ne(2, i) + ne(3, i) + ne(4, i)
2571 : END DO
2572 1610 : s = s + 1._dp*REAL(m, dp)
2573 : END IF
2574 : ! Slater exponent is (Z-S)/NS
2575 21570 : srules = (REAL(z, dp) - s)/xns(nn)
2576 21570 : END FUNCTION srules
2577 :
2578 : ! **************************************************************************************************
2579 : !> \brief sort basis sets w.r.t. radius
2580 : !> \param basis_set ...
2581 : !> \param sort_method ...
2582 : ! **************************************************************************************************
2583 11630 : SUBROUTINE sort_gto_basis_set(basis_set, sort_method)
2584 : TYPE(gto_basis_set_type), INTENT(INOUT) :: basis_set
2585 : INTEGER, INTENT(IN) :: sort_method
2586 :
2587 11630 : CHARACTER(LEN=12), DIMENSION(:), POINTER :: cgf_symbol
2588 11630 : CHARACTER(LEN=6), DIMENSION(:), POINTER :: sgf_symbol
2589 : INTEGER :: ic, ic_max, icgf, icgf_new, icgf_old, ico, is, is_max, iset, isgf, isgf_new, &
2590 : isgf_old, ishell, lshell, maxco, maxpgf, maxshell, mm, nc, ncgf, ns, nset
2591 11630 : INTEGER, ALLOCATABLE, DIMENSION(:) :: sort_index
2592 11630 : INTEGER, ALLOCATABLE, DIMENSION(:, :) :: icgf_set, isgf_set
2593 11630 : INTEGER, DIMENSION(:), POINTER :: lx, ly, lz, m, npgf
2594 11630 : REAL(dp), ALLOCATABLE, DIMENSION(:) :: tmp
2595 11630 : REAL(dp), DIMENSION(:), POINTER :: set_radius
2596 11630 : REAL(dp), DIMENSION(:, :), POINTER :: zet
2597 11630 : REAL(KIND=dp), DIMENSION(:), POINTER :: norm_cgf
2598 11630 : REAL(KIND=dp), DIMENSION(:, :), POINTER :: cphi, scon, sphi
2599 :
2600 11630 : NULLIFY (set_radius, zet)
2601 :
2602 11630 : IF (sort_method == basis_sort_default) RETURN
2603 :
2604 : CALL get_gto_basis_set(gto_basis_set=basis_set, &
2605 : nset=nset, &
2606 : maxshell=maxshell, &
2607 : maxpgf=maxpgf, &
2608 : maxco=maxco, &
2609 : ncgf=ncgf, &
2610 : npgf=npgf, &
2611 : set_radius=set_radius, &
2612 728 : zet=zet)
2613 :
2614 2184 : ALLOCATE (sort_index(nset))
2615 2184 : ALLOCATE (tmp(nset))
2616 : SELECT CASE (sort_method)
2617 : CASE (basis_sort_zet)
2618 4508 : DO iset = 1, nset
2619 13856 : tmp(iset) = MINVAL(basis_set%zet(:npgf(iset), iset))
2620 : END DO
2621 : CASE DEFAULT
2622 728 : CPABORT("Request basis sort criterion not implemented.")
2623 : END SELECT
2624 :
2625 728 : CALL sort(tmp(1:nset), nset, sort_index)
2626 :
2627 728 : ic_max = 0
2628 728 : is_max = 0
2629 4508 : DO iset = 1, nset
2630 3780 : ic = 0
2631 3780 : is = 0
2632 10104 : DO ishell = 1, basis_set%nshell(iset)
2633 22072 : DO ico = 1, nco(basis_set%l(ishell, iset))
2634 16476 : ic = ic + 1
2635 22072 : IF (ic > ic_max) ic_max = ic
2636 : END DO
2637 5596 : lshell = basis_set%l(ishell, iset)
2638 24092 : DO mm = -lshell, lshell
2639 14716 : is = is + 1
2640 20312 : IF (is > is_max) is_max = is
2641 : END DO
2642 : END DO
2643 : END DO
2644 :
2645 728 : icgf = 0
2646 728 : isgf = 0
2647 2912 : ALLOCATE (icgf_set(nset, ic_max))
2648 41850 : icgf_set(:, :) = 0
2649 2912 : ALLOCATE (isgf_set(nset, is_max))
2650 34206 : isgf_set(:, :) = 0
2651 :
2652 4508 : DO iset = 1, nset
2653 3780 : ic = 0
2654 3780 : is = 0
2655 10104 : DO ishell = 1, basis_set%nshell(iset)
2656 22072 : DO ico = 1, nco(basis_set%l(ishell, iset))
2657 16476 : icgf = icgf + 1
2658 16476 : ic = ic + 1
2659 22072 : icgf_set(iset, ic) = icgf
2660 : END DO
2661 5596 : lshell = basis_set%l(ishell, iset)
2662 24092 : DO mm = -lshell, lshell
2663 14716 : isgf = isgf + 1
2664 14716 : is = is + 1
2665 20312 : isgf_set(iset, is) = isgf
2666 : END DO
2667 : END DO
2668 : END DO
2669 :
2670 2184 : ALLOCATE (cgf_symbol(SIZE(basis_set%cgf_symbol)))
2671 2184 : ALLOCATE (norm_cgf(SIZE(basis_set%norm_cgf)))
2672 2184 : ALLOCATE (lx(SIZE(basis_set%lx)))
2673 2184 : ALLOCATE (ly(SIZE(basis_set%ly)))
2674 2184 : ALLOCATE (lz(SIZE(basis_set%lz)))
2675 2912 : ALLOCATE (cphi(SIZE(basis_set%cphi, 1), SIZE(basis_set%cphi, 2)))
2676 603940 : cphi = 0.0_dp
2677 2912 : ALLOCATE (sphi(SIZE(basis_set%sphi, 1), SIZE(basis_set%sphi, 2)))
2678 533192 : sphi = 0.0_dp
2679 2912 : ALLOCATE (scon(SIZE(basis_set%scon, 1), SIZE(basis_set%scon, 2)))
2680 533192 : scon = 0.0_dp
2681 :
2682 2184 : ALLOCATE (sgf_symbol(SIZE(basis_set%sgf_symbol)))
2683 2184 : ALLOCATE (m(SIZE(basis_set%m)))
2684 :
2685 : icgf_new = 0
2686 : isgf_new = 0
2687 4508 : DO iset = 1, nset
2688 37276 : DO ic = 1, ic_max
2689 33496 : icgf_old = icgf_set(sort_index(iset), ic)
2690 33496 : IF (icgf_old == 0) CYCLE
2691 16476 : icgf_new = icgf_new + 1
2692 16476 : norm_cgf(icgf_new) = basis_set%norm_cgf(icgf_old)
2693 16476 : lx(icgf_new) = basis_set%lx(icgf_old)
2694 16476 : ly(icgf_new) = basis_set%ly(icgf_old)
2695 16476 : lz(icgf_new) = basis_set%lz(icgf_old)
2696 603212 : cphi(:, icgf_new) = basis_set%cphi(:, icgf_old)
2697 37276 : cgf_symbol(icgf_new) = basis_set%cgf_symbol(icgf_old)
2698 : END DO
2699 31284 : DO is = 1, is_max
2700 26776 : isgf_old = isgf_set(sort_index(iset), is)
2701 26776 : IF (isgf_old == 0) CYCLE
2702 14716 : isgf_new = isgf_new + 1
2703 14716 : m(isgf_new) = basis_set%m(isgf_old)
2704 532464 : sphi(:, isgf_new) = basis_set%sphi(:, isgf_old)
2705 532464 : scon(:, isgf_new) = basis_set%scon(:, isgf_old)
2706 30556 : sgf_symbol(isgf_new) = basis_set%sgf_symbol(isgf_old)
2707 : END DO
2708 : END DO
2709 :
2710 728 : DEALLOCATE (basis_set%cgf_symbol)
2711 728 : basis_set%cgf_symbol => cgf_symbol
2712 728 : DEALLOCATE (basis_set%norm_cgf)
2713 728 : basis_set%norm_cgf => norm_cgf
2714 728 : DEALLOCATE (basis_set%lx)
2715 728 : basis_set%lx => lx
2716 728 : DEALLOCATE (basis_set%ly)
2717 728 : basis_set%ly => ly
2718 728 : DEALLOCATE (basis_set%lz)
2719 728 : basis_set%lz => lz
2720 728 : DEALLOCATE (basis_set%cphi)
2721 728 : basis_set%cphi => cphi
2722 728 : DEALLOCATE (basis_set%sphi)
2723 728 : basis_set%sphi => sphi
2724 728 : DEALLOCATE (basis_set%scon)
2725 728 : basis_set%scon => scon
2726 :
2727 728 : DEALLOCATE (basis_set%m)
2728 728 : basis_set%m => m
2729 728 : DEALLOCATE (basis_set%sgf_symbol)
2730 728 : basis_set%sgf_symbol => sgf_symbol
2731 :
2732 8288 : basis_set%lmax = basis_set%lmax(sort_index)
2733 8288 : basis_set%lmin = basis_set%lmin(sort_index)
2734 8288 : basis_set%npgf = basis_set%npgf(sort_index)
2735 8288 : basis_set%nshell = basis_set%nshell(sort_index)
2736 8288 : basis_set%ncgf_set = basis_set%ncgf_set(sort_index)
2737 8288 : basis_set%nsgf_set = basis_set%nsgf_set(sort_index)
2738 :
2739 22148 : basis_set%n(:, :) = basis_set%n(:, sort_index)
2740 22148 : basis_set%l(:, :) = basis_set%l(:, sort_index)
2741 24968 : basis_set%zet(:, :) = basis_set%zet(:, sort_index)
2742 :
2743 92184 : basis_set%gcc(:, :, :) = basis_set%gcc(:, :, sort_index)
2744 8288 : basis_set%set_radius(:) = basis_set%set_radius(sort_index)
2745 24968 : basis_set%pgf_radius(:, :) = basis_set%pgf_radius(:, sort_index)
2746 :
2747 728 : nc = 0
2748 728 : ns = 0
2749 4508 : DO iset = 1, nset
2750 10104 : DO ishell = 1, basis_set%nshell(iset)
2751 5596 : lshell = basis_set%l(ishell, iset)
2752 5596 : basis_set%first_cgf(ishell, iset) = nc + 1
2753 5596 : nc = nc + nco(lshell)
2754 5596 : basis_set%last_cgf(ishell, iset) = nc
2755 5596 : basis_set%first_sgf(ishell, iset) = ns + 1
2756 5596 : ns = ns + nso(lshell)
2757 9376 : basis_set%last_sgf(ishell, iset) = ns
2758 : END DO
2759 : END DO
2760 :
2761 12358 : END SUBROUTINE
2762 :
2763 0 : END MODULE basis_set_types
|