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