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