58#include "./base/base_uses.f90"
64 CHARACTER(len=*),
PARAMETER,
PRIVATE :: moduleN =
'atom_types'
67 INTEGER,
PARAMETER ::
lmat = 5
74 INTEGER,
PARAMETER :: nmax = 25
80 INTEGER,
DIMENSION(0:lmat) :: nbas = 0
81 INTEGER,
DIMENSION(0:lmat) :: nprim = 0
82 REAL(kind=
dp),
DIMENSION(:, :),
POINTER :: am => null()
83 REAL(kind=
dp),
DIMENSION(:, :, :),
POINTER :: cm => null()
84 REAL(kind=
dp),
DIMENSION(:, :),
POINTER :: as => null()
85 INTEGER,
DIMENSION(:, :),
POINTER :: ns => null()
86 REAL(kind=
dp),
DIMENSION(:, :, :),
POINTER :: bf => null()
87 REAL(kind=
dp),
DIMENSION(:, :, :),
POINTER :: dbf => null()
88 REAL(kind=
dp),
DIMENSION(:, :, :),
POINTER :: ddbf => null()
89 REAL(kind=
dp) :: eps_eig = 0.0_dp
91 LOGICAL :: geometrical = .false.
92 REAL(kind=
dp) :: aval = 0.0_dp, cval = 0.0_dp
93 INTEGER,
DIMENSION(0:lmat) :: start = 0
99 CHARACTER(LEN=2) :: symbol =
""
100 CHARACTER(LEN=default_string_length) :: pname =
""
101 INTEGER,
DIMENSION(0:lmat) :: econf = 0
102 REAL(
dp) :: zion = 0.0_dp
103 REAL(
dp) :: rc = 0.0_dp
105 REAL(
dp),
DIMENSION(5) :: cl = 0.0_dp
106 INTEGER,
DIMENSION(0:lmat) :: nl = 0
107 REAL(
dp),
DIMENSION(0:lmat) :: rcnl = 0.0_dp
108 REAL(
dp),
DIMENSION(4, 4, 0:lmat) :: hnl = 0.0_dp
110 LOGICAL :: soc = .false.
111 REAL(
dp),
DIMENSION(4, 4, 0:lmat) :: knl = 0.0_dp
114 LOGICAL :: nlcc = .false.
115 INTEGER :: nexp_nlcc = 0
116 REAL(kind=
dp),
DIMENSION(10) :: alpha_nlcc = 0.0_dp
117 INTEGER,
DIMENSION(10) :: nct_nlcc = 0
118 REAL(kind=
dp),
DIMENSION(4, 10) :: cval_nlcc = 0.0_dp
120 LOGICAL :: lsdpot = .false.
121 INTEGER :: nexp_lsd = 0
122 REAL(kind=
dp),
DIMENSION(10) :: alpha_lsd = 0.0_dp
123 INTEGER,
DIMENSION(10) :: nct_lsd = 0
124 REAL(kind=
dp),
DIMENSION(4, 10) :: cval_lsd = 0.0_dp
126 LOGICAL :: lpotextended = .false.
127 INTEGER :: nexp_lpot = 0
128 REAL(kind=
dp),
DIMENSION(10) :: alpha_lpot = 0.0_dp
129 INTEGER,
DIMENSION(10) :: nct_lpot = 0
130 REAL(kind=
dp),
DIMENSION(4, 10) :: cval_lpot = 0.0_dp
134 CHARACTER(LEN=2) :: symbol =
""
135 CHARACTER(LEN=default_string_length) :: pname =
""
136 INTEGER,
DIMENSION(0:lmat) :: econf = 0
137 REAL(
dp) :: zion = 0.0_dp
140 INTEGER,
DIMENSION(1:15) :: nrloc = 0
141 REAL(
dp),
DIMENSION(1:15) :: aloc = 0.0_dp
142 REAL(
dp),
DIMENSION(1:15) :: bloc = 0.0_dp
143 INTEGER,
DIMENSION(0:10) :: npot = 0
144 INTEGER,
DIMENSION(1:15, 0:10) :: nrpot = 0
145 REAL(
dp),
DIMENSION(1:15, 0:10) :: apot = 0.0_dp
146 REAL(
dp),
DIMENSION(1:15, 0:10) :: bpot = 0.0_dp
150 CHARACTER(LEN=2) :: symbol =
""
151 CHARACTER(LEN=default_string_length) :: pname =
""
152 INTEGER,
DIMENSION(0:lmat) :: econf = 0
153 REAL(
dp) :: zion = 0.0_dp
155 LOGICAL :: has_nonlocal = .false.
156 INTEGER :: n_nonlocal = 0
157 LOGICAL,
DIMENSION(0:5) :: is_nonlocal = .false.
158 REAL(kind=
dp),
DIMENSION(nmax) :: a_nonlocal = 0.0_dp
159 REAL(kind=
dp),
DIMENSION(nmax, 0:lmat) :: h_nonlocal = 0.0_dp
160 REAL(kind=
dp),
DIMENSION(nmax, nmax, 0:lmat) :: c_nonlocal = 0.0_dp
161 INTEGER :: n_local = 0
162 REAL(kind=
dp) :: ac_local = 0.0_dp
163 REAL(kind=
dp),
DIMENSION(nmax) :: a_local = 0.0_dp
164 REAL(kind=
dp),
DIMENSION(nmax) :: c_local = 0.0_dp
165 LOGICAL :: has_nlcc = .false.
166 INTEGER :: n_nlcc = 0
167 REAL(kind=
dp),
DIMENSION(nmax) :: a_nlcc = 0.0_dp
168 REAL(kind=
dp),
DIMENSION(nmax) :: c_nlcc = 0.0_dp
172 INTEGER :: ppot_type = 0
173 LOGICAL :: confinement = .false.
174 INTEGER :: conf_type = 0
175 REAL(
dp) :: acon = 0.0_dp
176 REAL(
dp) :: rcon = 0.0_dp
177 REAL(
dp) :: scon = 0.0_dp
188 REAL(kind=
dp) :: scale_coulomb = 0.0_dp
189 REAL(kind=
dp) :: scale_longrange = 0.0_dp
190 REAL(kind=
dp) :: omega = 0.0_dp
191 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:, :, :) :: kernel
192 LOGICAL :: do_gh = .false.
199 REAL(kind=
dp),
DIMENSION(0:lmat, 10) :: occ = 0.0_dp
200 REAL(kind=
dp),
DIMENSION(0:lmat, 10) :: core = 0.0_dp
201 REAL(kind=
dp),
DIMENSION(0:lmat, 10) :: occupation = 0.0_dp
202 INTEGER :: maxl_occ = 0
203 INTEGER,
DIMENSION(0:lmat) :: maxn_occ = 0
204 INTEGER :: maxl_calc = 0
205 INTEGER,
DIMENSION(0:lmat) :: maxn_calc = 0
206 INTEGER :: multiplicity = 0
207 REAL(kind=
dp),
DIMENSION(0:lmat, 10) :: occa = 0.0_dp, occb = 0.0_dp
213 REAL(kind=
dp),
DIMENSION(:, :),
POINTER :: int => null()
217 INTEGER :: status = 0
218 INTEGER :: ppstat = 0
219 LOGICAL :: eri_coulomb = .false.
220 LOGICAL :: eri_exchange = .false.
221 LOGICAL :: all_nu = .false.
222 INTEGER,
DIMENSION(0:lmat) :: n = 0, nne = 0
223 REAL(kind=
dp),
DIMENSION(:, :, :),
POINTER :: ovlp => null(), kin => null(), core => null(), clsd => null()
224 REAL(kind=
dp),
DIMENSION(:, :, :),
POINTER :: utrans => null(), uptrans => null()
225 REAL(kind=
dp),
DIMENSION(:, :, :),
POINTER :: hnl => null()
226 REAL(kind=
dp),
DIMENSION(:, :, :),
POINTER :: conf => null()
227 TYPE(
eri),
DIMENSION(100) :: ceri =
eri()
228 TYPE(
eri),
DIMENSION(100) :: eeri =
eri()
229 INTEGER :: dkhstat = 0
230 INTEGER :: zorastat = 0
231 REAL(kind=
dp),
DIMENSION(:, :, :),
POINTER :: tzora => null()
232 REAL(kind=
dp),
DIMENSION(:, :, :),
POINTER :: hdkh => null()
238 REAL(kind=
dp),
DIMENSION(:, :, :),
POINTER :: wfn => null(), wfna => null(), wfnb => null()
239 REAL(kind=
dp),
DIMENSION(:, :, :),
POINTER :: pmat => null(), pmata => null(), pmatb => null()
240 REAL(kind=
dp),
DIMENSION(:, :),
POINTER :: ener => null(), enera => null(), enerb => null()
241 REAL(kind=
dp),
DIMENSION(:, :, :),
POINTER :: refene => null(), refchg => null(), refnod => null()
242 REAL(kind=
dp),
DIMENSION(:, :, :),
POINTER :: wrefene => null(), wrefchg => null(), wrefnod => null()
243 REAL(kind=
dp),
DIMENSION(:, :, :),
POINTER :: crefene => null(), crefchg => null(), crefnod => null()
244 REAL(kind=
dp),
DIMENSION(:, :),
POINTER :: wpsir0 => null(), tpsir0 => null()
245 REAL(kind=
dp),
DIMENSION(:, :, :),
POINTER :: rcmax => null()
246 CHARACTER(LEN=2),
DIMENSION(:, :, :),
POINTER :: reftype => null()
252 INTEGER,
DIMENSION(0:lmat) :: n = 0
253 REAL(kind=
dp),
DIMENSION(:, :, :),
POINTER :: op => null()
259 REAL(kind=
dp),
DIMENSION(:),
POINTER :: op => null()
265 TYPE atom_energy_type
266 REAL(kind=
dp) :: etot = 0.0_dp
267 REAL(kind=
dp) :: eband = 0.0_dp
268 REAL(kind=
dp) :: ekin = 0.0_dp
269 REAL(kind=
dp) :: epot = 0.0_dp
270 REAL(kind=
dp) :: ecore = 0.0_dp
271 REAL(kind=
dp) :: elsd = 0.0_dp
272 REAL(kind=
dp) :: epseudo = 0.0_dp
273 REAL(kind=
dp) :: eploc = 0.0_dp
274 REAL(kind=
dp) :: epnl = 0.0_dp
275 REAL(kind=
dp) :: exc = 0.0_dp
276 REAL(kind=
dp) :: ecoulomb = 0.0_dp
277 REAL(kind=
dp) :: eexchange = 0.0_dp
278 REAL(kind=
dp) :: econfinement = 0.0_dp
279 END TYPE atom_energy_type
284 REAL(kind=
dp) :: damping = 0.0_dp
285 REAL(kind=
dp) :: eps_scf = 0.0_dp
286 REAL(kind=
dp) :: eps_diis = 0.0_dp
287 INTEGER :: max_iter = 0
288 INTEGER :: n_diis = 0
296 LOGICAL :: pp_calc = .false.
298 LOGICAL :: do_zmp = .false., doread = .false., read_vxc = .false., dm = .false.
299 CHARACTER(LEN=default_string_length) :: ext_file =
"", ext_vxc_file =
"", &
300 zmp_restart_file =
""
307 REAL(kind=
dp) :: lambda = 0.0_dp
308 REAL(kind=
dp) :: rho_diff_integral = 0.0_dp
309 REAL(kind=
dp) :: weight = 0.0_dp, zmpgrid_tol = 0.0_dp, zmpvxcgrid_tol = 0.0_dp
316 TYPE(atom_energy_type) :: energy = atom_energy_type()
380 INTEGER,
INTENT(IN) :: zval
381 CHARACTER(LEN=2) :: btyp
383 CHARACTER(LEN=*),
PARAMETER :: routinen =
'init_atom_basis'
384 INTEGER,
PARAMETER :: nua = 40, nup = 16
385 REAL(kind=
dp),
DIMENSION(nua),
PARAMETER :: ugbs = (/0.007299_dp, 0.013705_dp, 0.025733_dp, &
386 0.048316_dp, 0.090718_dp, 0.170333_dp, 0.319819_dp, 0.600496_dp, 1.127497_dp, 2.117000_dp,&
387 3.974902_dp, 7.463317_dp, 14.013204_dp, 26.311339_dp, 49.402449_dp, 92.758561_dp, &
388 174.164456_dp, 327.013024_dp, 614.003114_dp, 1152.858743_dp, 2164.619772_dp, &
389 4064.312984_dp, 7631.197056_dp, 14328.416324_dp, 26903.186074_dp, 50513.706789_dp, &
390 94845.070265_dp, 178082.107320_dp, 334368.848683_dp, 627814.487663_dp, 1178791.123851_dp, &
391 2213310.684886_dp, 4155735.557141_dp, 7802853.046713_dp, 14650719.428954_dp, &
392 27508345.793637_dp, 51649961.080194_dp, 96978513.342764_dp, 182087882.613702_dp, &
393 341890134.751331_dp/)
395 CHARACTER(LEN=default_string_length) :: basis_fn, basis_name
396 INTEGER :: basistype, handle, i, j, k, l, ll, m, &
397 ngp, nl, nr, nu, quadtype
398 INTEGER,
DIMENSION(0:lmat) :: starti
399 INTEGER,
DIMENSION(:),
POINTER :: nqm, num_gto, num_slater, sindex
400 REAL(kind=
dp) :: al, amax, aval, cval, ear, pf, rk
401 REAL(kind=
dp),
DIMENSION(:),
POINTER :: expo
404 CALL timeset(routinen, handle)
411 NULLIFY (basis%am, basis%cm, basis%as, basis%ns, basis%bf, basis%dbf, basis%ddbf)
418 cpabort(
"# point radial grid < 0")
421 basis%geometrical = .false.
428 SELECT CASE (basistype)
435 IF (num_gto(1) < 1)
THEN
437 IF (btyp ==
"AE")
THEN
439 ELSEIF (btyp ==
"PP")
THEN
446 ALLOCATE (basis%am(nu, 0:
lmat))
448 basis%am(1:nu, i) = ugbs(1:nu)
452 DO i = 1,
SIZE(num_gto)
453 basis%nbas(i - 1) = num_gto(i)
455 basis%nprim = basis%nbas
456 m = maxval(basis%nbas)
457 ALLOCATE (basis%am(m, 0:
lmat))
460 IF (basis%nbas(l) > 0)
THEN
464 cpabort(
"Atom Basis")
474 cpassert(
SIZE(expo) >= basis%nbas(l))
475 DO i = 1, basis%nbas(l)
476 basis%am(i, l) = expo(i)
483 m = maxval(basis%nbas)
484 ALLOCATE (basis%bf(nr, m, 0:
lmat))
485 ALLOCATE (basis%dbf(nr, m, 0:
lmat))
486 ALLOCATE (basis%ddbf(nr, m, 0:
lmat))
491 DO i = 1, basis%nbas(l)
494 rk = basis%grid%rad(k)
495 ear = exp(-al*basis%grid%rad(k)**2)
496 basis%bf(k, i, l) = rk**l*ear
497 basis%dbf(k, i, l) = (real(l,
dp)*rk**(l - 1) - 2._dp*al*rk**(l + 1))*ear
498 basis%ddbf(k, i, l) = (real(l*(l - 1),
dp)*rk**(l - 2) - &
499 2._dp*al*real(2*l + 1,
dp)*rk**(l) + 4._dp*al*rk**(l + 2))*ear
507 IF (num_gto(1) < 1)
THEN
508 IF (btyp ==
"AE")
THEN
511 ELSEIF (btyp ==
"PP")
THEN
514 ELSEIF (btyp ==
"AA")
THEN
516 amax = cval**(basis%nbas(0) - 1)
517 basis%nbas(0) = nint((log(amax)/log(1.6_dp)))
520 basis%nbas(1) = basis%nbas(0) - 4
521 basis%nbas(2) = basis%nbas(0) - 8
522 basis%nbas(3) = basis%nbas(0) - 12
523 IF (
lmat > 3) basis%nbas(4:
lmat) = 0
524 ELSEIF (btyp ==
"AP")
THEN
527 basis%nbas = nint((log(amax)/log(1.6_dp)))
534 basis%nprim = basis%nbas
537 DO i = 1,
SIZE(num_gto)
538 basis%nbas(i - 1) = num_gto(i)
540 basis%nprim = basis%nbas
544 DO i = 1,
SIZE(sindex)
545 starti(i - 1) = sindex(i)
546 cpassert(sindex(i) >= 0)
551 m = maxval(basis%nbas)
552 ALLOCATE (basis%am(m, 0:
lmat))
555 DO i = 1, basis%nbas(l)
556 ll = i - 1 + starti(l)
557 basis%am(i, l) = aval*cval**(ll)
561 basis%geometrical = .true.
568 m = maxval(basis%nbas)
569 ALLOCATE (basis%bf(nr, m, 0:
lmat))
570 ALLOCATE (basis%dbf(nr, m, 0:
lmat))
571 ALLOCATE (basis%ddbf(nr, m, 0:
lmat))
576 DO i = 1, basis%nbas(l)
579 rk = basis%grid%rad(k)
580 ear = exp(-al*basis%grid%rad(k)**2)
581 basis%bf(k, i, l) = rk**l*ear
582 basis%dbf(k, i, l) = (real(l,
dp)*rk**(l - 1) - 2._dp*al*rk**(l + 1))*ear
583 basis%ddbf(k, i, l) = (real(l*(l - 1),
dp)*rk**(l - 2) - &
584 2._dp*al*real(2*l + 1,
dp)*rk**(l) + 4._dp*al*rk**(l + 2))*ear
593 CALL read_basis_set(
ptable(zval)%symbol, basis, basis_name, basis_fn, &
598 m = maxval(basis%nbas)
599 ALLOCATE (basis%bf(nr, m, 0:
lmat))
600 ALLOCATE (basis%dbf(nr, m, 0:
lmat))
601 ALLOCATE (basis%ddbf(nr, m, 0:
lmat))
606 DO i = 1, basis%nprim(l)
609 rk = basis%grid%rad(k)
610 ear = exp(-al*basis%grid%rad(k)**2)
611 DO j = 1, basis%nbas(l)
612 basis%bf(k, j, l) = basis%bf(k, j, l) + rk**l*ear*basis%cm(i, j, l)
613 basis%dbf(k, j, l) = basis%dbf(k, j, l) &
614 + (real(l,
dp)*rk**(l - 1) - 2._dp*al*rk**(l + 1))*ear*basis%cm(i, j, l)
615 basis%ddbf(k, j, l) = basis%ddbf(k, j, l) + &
616 (real(l*(l - 1),
dp)*rk**(l - 2) - 2._dp*al*real(2*l + 1,
dp)*rk**(l) + 4._dp*al*rk**(l + 2))* &
617 ear*basis%cm(i, j, l)
626 IF (num_slater(1) < 1)
THEN
630 DO i = 1,
SIZE(num_slater)
631 basis%nbas(i - 1) = num_slater(i)
633 basis%nprim = basis%nbas
634 m = maxval(basis%nbas)
635 ALLOCATE (basis%as(m, 0:
lmat), basis%ns(m, 0:
lmat))
639 IF (basis%nbas(l) > 0)
THEN
643 cpabort(
"Atom Basis")
653 cpassert(
SIZE(expo) >= basis%nbas(l))
654 DO i = 1, basis%nbas(l)
655 basis%as(i, l) = expo(i)
660 cpabort(
"Atom Basis")
670 cpassert(
SIZE(nqm) >= basis%nbas(l))
671 DO i = 1, basis%nbas(l)
672 basis%ns(i, l) = nqm(i)
679 m = maxval(basis%nbas)
680 ALLOCATE (basis%bf(nr, m, 0:
lmat))
681 ALLOCATE (basis%dbf(nr, m, 0:
lmat))
682 ALLOCATE (basis%ddbf(nr, m, 0:
lmat))
687 DO i = 1, basis%nbas(l)
690 pf = (2._dp*al)**nl*sqrt(2._dp*al/
fac(2*nl))
692 rk = basis%grid%rad(k)
693 ear = rk**(nl - 1)*exp(-al*rk)
694 basis%bf(k, i, l) = pf*ear
695 basis%dbf(k, i, l) = pf*(real(nl - 1,
dp)/rk - al)*ear
696 basis%ddbf(k, i, l) = pf*(real((nl - 2)*(nl - 1),
dp)/rk/rk &
697 - al*real(2*(nl - 1),
dp)/rk + al*al)*ear
706 CALL timestop(handle)
717 CHARACTER(LEN=*),
PARAMETER :: routinen =
'init_atom_basis_default_pp'
718 INTEGER,
PARAMETER :: nua = 40, nup = 20
719 REAL(kind=
dp),
DIMENSION(nua),
PARAMETER :: ugbs = (/0.007299_dp, 0.013705_dp, 0.025733_dp, &
720 0.048316_dp, 0.090718_dp, 0.170333_dp, 0.319819_dp, 0.600496_dp, 1.127497_dp, 2.117000_dp,&
721 3.974902_dp, 7.463317_dp, 14.013204_dp, 26.311339_dp, 49.402449_dp, 92.758561_dp, &
722 174.164456_dp, 327.013024_dp, 614.003114_dp, 1152.858743_dp, 2164.619772_dp, &
723 4064.312984_dp, 7631.197056_dp, 14328.416324_dp, 26903.186074_dp, 50513.706789_dp, &
724 94845.070265_dp, 178082.107320_dp, 334368.848683_dp, 627814.487663_dp, 1178791.123851_dp, &
725 2213310.684886_dp, 4155735.557141_dp, 7802853.046713_dp, 14650719.428954_dp, &
726 27508345.793637_dp, 51649961.080194_dp, 96978513.342764_dp, 182087882.613702_dp, &
727 341890134.751331_dp/)
729 INTEGER :: handle, i, k, l, m, ngp, nr, nu, quadtype
730 REAL(kind=
dp) :: al, ear, rk
732 CALL timeset(routinen, handle)
734 NULLIFY (basis%am, basis%cm, basis%as, basis%ns, basis%bf, basis%dbf, basis%ddbf)
743 basis%geometrical = .false.
747 basis%eps_eig = 1.e-12_dp
753 ALLOCATE (basis%am(nu, 0:
lmat))
755 basis%am(1:nu, i) = ugbs(1:nu)
759 m = maxval(basis%nbas)
760 ALLOCATE (basis%bf(nr, m, 0:
lmat))
761 ALLOCATE (basis%dbf(nr, m, 0:
lmat))
762 ALLOCATE (basis%ddbf(nr, m, 0:
lmat))
767 DO i = 1, basis%nbas(l)
770 rk = basis%grid%rad(k)
771 ear = exp(-al*basis%grid%rad(k)**2)
772 basis%bf(k, i, l) = rk**l*ear
773 basis%dbf(k, i, l) = (real(l,
dp)*rk**(l - 1) - 2._dp*al*rk**(l + 1))*ear
774 basis%ddbf(k, i, l) = (real(l*(l - 1),
dp)*rk**(l - 2) - &
775 2._dp*al*real(2*l + 1,
dp)*rk**(l) + 4._dp*al*rk**(l + 2))*ear
780 CALL timestop(handle)
794 REAL(kind=
dp),
DIMENSION(:),
INTENT(IN) :: r, rab
796 INTEGER :: i, j, k, l, m, n1, n2, n3, ngp, nl, nr, &
798 REAL(kind=
dp) :: al, ear, pf, rk
800 NULLIFY (gbasis%am, gbasis%cm, gbasis%as, gbasis%ns, gbasis%bf, gbasis%dbf, gbasis%ddbf)
803 gbasis%basis_type = basis%basis_type
804 gbasis%nbas(0:
lmat) = basis%nbas(0:
lmat)
805 gbasis%nprim(0:
lmat) = basis%nprim(0:
lmat)
806 IF (
ASSOCIATED(basis%am))
THEN
807 n1 =
SIZE(basis%am, 1)
808 n2 =
SIZE(basis%am, 2)
809 ALLOCATE (gbasis%am(n1, 0:n2 - 1))
812 IF (
ASSOCIATED(basis%cm))
THEN
813 n1 =
SIZE(basis%cm, 1)
814 n2 =
SIZE(basis%cm, 2)
815 n3 =
SIZE(basis%cm, 3)
816 ALLOCATE (gbasis%cm(n1, n2, 0:n3 - 1))
819 IF (
ASSOCIATED(basis%as))
THEN
820 n1 =
SIZE(basis%as, 1)
821 n2 =
SIZE(basis%as, 2)
822 ALLOCATE (gbasis%as(n1, 0:n2 - 1))
825 IF (
ASSOCIATED(basis%ns))
THEN
826 n1 =
SIZE(basis%ns, 1)
827 n2 =
SIZE(basis%ns, 2)
828 ALLOCATE (gbasis%ns(n1, 0:n2 - 1))
831 gbasis%eps_eig = basis%eps_eig
832 gbasis%geometrical = basis%geometrical
833 gbasis%aval = basis%aval
834 gbasis%cval = basis%cval
835 gbasis%start(0:
lmat) = basis%start(0:
lmat)
839 NULLIFY (gbasis%grid)
844 cpabort(
"# point radial grid < 0")
847 gbasis%grid%rad(:) = r(:)
848 gbasis%grid%rad2(:) = r(:)*r(:)
849 gbasis%grid%wr(:) = rab(:)*gbasis%grid%rad2(:)
853 m = maxval(gbasis%nbas)
854 ALLOCATE (gbasis%bf(nr, m, 0:
lmat))
855 ALLOCATE (gbasis%dbf(nr, m, 0:
lmat))
856 ALLOCATE (gbasis%ddbf(nr, m, 0:
lmat))
861 SELECT CASE (gbasis%basis_type)
866 DO i = 1, gbasis%nbas(l)
869 rk = gbasis%grid%rad(k)
870 ear = exp(-al*gbasis%grid%rad(k)**2)
871 gbasis%bf(k, i, l) = rk**l*ear
872 gbasis%dbf(k, i, l) = (real(l,
dp)*rk**(l - 1) - 2._dp*al*rk**(l + 1))*ear
873 gbasis%ddbf(k, i, l) = (real(l*(l - 1),
dp)*rk**(l - 2) - &
874 2._dp*al*real(2*l + 1,
dp)*rk**(l) + 4._dp*al*rk**(l + 2))*ear
880 DO i = 1, gbasis%nprim(l)
883 rk = gbasis%grid%rad(k)
884 ear = exp(-al*gbasis%grid%rad(k)**2)
885 DO j = 1, gbasis%nbas(l)
886 gbasis%bf(k, j, l) = gbasis%bf(k, j, l) + rk**l*ear*gbasis%cm(i, j, l)
887 gbasis%dbf(k, j, l) = gbasis%dbf(k, j, l) &
888 + (real(l,
dp)*rk**(l - 1) - 2._dp*al*rk**(l + 1))*ear*gbasis%cm(i, j, l)
889 gbasis%ddbf(k, j, l) = gbasis%ddbf(k, j, l) + &
890 (real(l*(l - 1),
dp)*rk**(l - 2) - 2._dp*al*real(2*l + 1,
dp)*rk**(l) + 4._dp*al*rk**(l + 2))* &
891 ear*gbasis%cm(i, j, l)
898 DO i = 1, gbasis%nbas(l)
901 pf = (2._dp*al)**nl*sqrt(2._dp*al/
fac(2*nl))
903 rk = gbasis%grid%rad(k)
904 ear = rk**(nl - 1)*exp(-al*rk)
905 gbasis%bf(k, i, l) = pf*ear
906 gbasis%dbf(k, i, l) = pf*(real(nl - 1,
dp)/rk - al)*ear
907 gbasis%ddbf(k, i, l) = pf*(real((nl - 2)*(nl - 1),
dp)/rk/rk &
908 - al*real(2*(nl - 1),
dp)/rk + al*al)*ear
926 IF (
ASSOCIATED(basis%am))
THEN
927 DEALLOCATE (basis%am)
929 IF (
ASSOCIATED(basis%cm))
THEN
930 DEALLOCATE (basis%cm)
932 IF (
ASSOCIATED(basis%as))
THEN
933 DEALLOCATE (basis%as)
935 IF (
ASSOCIATED(basis%ns))
THEN
936 DEALLOCATE (basis%ns)
938 IF (
ASSOCIATED(basis%bf))
THEN
939 DEALLOCATE (basis%bf)
941 IF (
ASSOCIATED(basis%dbf))
THEN
942 DEALLOCATE (basis%dbf)
944 IF (
ASSOCIATED(basis%ddbf))
THEN
945 DEALLOCATE (basis%ddbf)
960 cpassert(.NOT.
ASSOCIATED(
atom))
964 NULLIFY (
atom%zmp_section)
965 NULLIFY (
atom%xc_section)
967 atom%do_zmp = .false.
968 atom%doread = .false.
969 atom%read_vxc = .false.
971 atom%hfx_pot%scale_coulomb = 0.0_dp
972 atom%hfx_pot%scale_longrange = 0.0_dp
973 atom%hfx_pot%omega = 0.0_dp
984 cpassert(
ASSOCIATED(
atom))
987 NULLIFY (
atom%integrals)
988 IF (
ASSOCIATED(
atom%state))
THEN
989 DEALLOCATE (
atom%state)
991 IF (
ASSOCIATED(
atom%orbitals))
THEN
1021 SUBROUTINE set_atom(atom, basis, state, integrals, orbitals, potential, zcore, pp_calc, do_zmp, doread, &
1022 read_vxc, method_type, relativistic, coulomb_integral_type, exchange_integral_type, fmat)
1029 INTEGER,
INTENT(IN),
OPTIONAL :: zcore
1030 LOGICAL,
INTENT(IN),
OPTIONAL :: pp_calc, do_zmp, doread, read_vxc
1031 INTEGER,
INTENT(IN),
OPTIONAL :: method_type, relativistic, &
1032 coulomb_integral_type, &
1033 exchange_integral_type
1036 cpassert(
ASSOCIATED(
atom))
1038 IF (
PRESENT(basis))
atom%basis => basis
1039 IF (
PRESENT(state))
atom%state => state
1040 IF (
PRESENT(integrals))
atom%integrals => integrals
1041 IF (
PRESENT(orbitals))
atom%orbitals => orbitals
1042 IF (
PRESENT(potential))
atom%potential => potential
1043 IF (
PRESENT(zcore))
atom%zcore = zcore
1044 IF (
PRESENT(pp_calc))
atom%pp_calc = pp_calc
1046 IF (
PRESENT(do_zmp))
atom%do_zmp = do_zmp
1047 IF (
PRESENT(doread))
atom%doread = doread
1048 IF (
PRESENT(read_vxc))
atom%read_vxc = read_vxc
1050 IF (
PRESENT(method_type))
atom%method_type = method_type
1051 IF (
PRESENT(relativistic))
atom%relativistic = relativistic
1052 IF (
PRESENT(coulomb_integral_type))
atom%coulomb_integral_type = coulomb_integral_type
1053 IF (
PRESENT(exchange_integral_type))
atom%exchange_integral_type = exchange_integral_type
1055 IF (
PRESENT(fmat))
THEN
1070 INTEGER,
INTENT(IN) :: mbas, mo
1072 cpassert(.NOT.
ASSOCIATED(orbs))
1076 ALLOCATE (orbs%wfn(mbas, mo, 0:
lmat), orbs%wfna(mbas, mo, 0:
lmat), orbs%wfnb(mbas, mo, 0:
lmat))
1081 ALLOCATE (orbs%pmat(mbas, mbas, 0:
lmat), orbs%pmata(mbas, mbas, 0:
lmat), orbs%pmatb(mbas, mbas, 0:
lmat))
1086 ALLOCATE (orbs%ener(mo, 0:
lmat), orbs%enera(mo, 0:
lmat), orbs%enerb(mo, 0:
lmat))
1091 ALLOCATE (orbs%refene(mo, 0:
lmat, 2), orbs%refchg(mo, 0:
lmat, 2), orbs%refnod(mo, 0:
lmat, 2))
1095 ALLOCATE (orbs%wrefene(mo, 0:
lmat, 2), orbs%wrefchg(mo, 0:
lmat, 2), orbs%wrefnod(mo, 0:
lmat, 2))
1096 orbs%wrefene = 0._dp
1097 orbs%wrefchg = 0._dp
1098 orbs%wrefnod = 0._dp
1099 ALLOCATE (orbs%crefene(mo, 0:
lmat, 2), orbs%crefchg(mo, 0:
lmat, 2), orbs%crefnod(mo, 0:
lmat, 2))
1100 orbs%crefene = 0._dp
1101 orbs%crefchg = 0._dp
1102 orbs%crefnod = 0._dp
1103 ALLOCATE (orbs%rcmax(mo, 0:
lmat, 2))
1105 ALLOCATE (orbs%wpsir0(mo, 2), orbs%tpsir0(mo, 2))
1108 ALLOCATE (orbs%reftype(mo, 0:
lmat, 2))
1120 cpassert(
ASSOCIATED(orbs))
1122 IF (
ASSOCIATED(orbs%wfn))
THEN
1123 DEALLOCATE (orbs%wfn, orbs%wfna, orbs%wfnb)
1125 IF (
ASSOCIATED(orbs%pmat))
THEN
1126 DEALLOCATE (orbs%pmat, orbs%pmata, orbs%pmatb)
1128 IF (
ASSOCIATED(orbs%ener))
THEN
1129 DEALLOCATE (orbs%ener, orbs%enera, orbs%enerb)
1131 IF (
ASSOCIATED(orbs%refene))
THEN
1132 DEALLOCATE (orbs%refene)
1134 IF (
ASSOCIATED(orbs%refchg))
THEN
1135 DEALLOCATE (orbs%refchg)
1137 IF (
ASSOCIATED(orbs%refnod))
THEN
1138 DEALLOCATE (orbs%refnod)
1140 IF (
ASSOCIATED(orbs%wrefene))
THEN
1141 DEALLOCATE (orbs%wrefene)
1143 IF (
ASSOCIATED(orbs%wrefchg))
THEN
1144 DEALLOCATE (orbs%wrefchg)
1146 IF (
ASSOCIATED(orbs%wrefnod))
THEN
1147 DEALLOCATE (orbs%wrefnod)
1149 IF (
ASSOCIATED(orbs%crefene))
THEN
1150 DEALLOCATE (orbs%crefene)
1152 IF (
ASSOCIATED(orbs%crefchg))
THEN
1153 DEALLOCATE (orbs%crefchg)
1155 IF (
ASSOCIATED(orbs%crefnod))
THEN
1156 DEALLOCATE (orbs%crefnod)
1158 IF (
ASSOCIATED(orbs%rcmax))
THEN
1159 DEALLOCATE (orbs%rcmax)
1161 IF (
ASSOCIATED(orbs%wpsir0))
THEN
1162 DEALLOCATE (orbs%wpsir0)
1164 IF (
ASSOCIATED(orbs%tpsir0))
THEN
1165 DEALLOCATE (orbs%tpsir0)
1167 IF (
ASSOCIATED(orbs%reftype))
THEN
1168 DEALLOCATE (orbs%reftype)
1184 REAL(kind=
dp),
INTENT(OUT) :: hf_frac
1185 LOGICAL,
INTENT(OUT) :: do_hfx
1188 INTEGER,
INTENT(IN) :: extype
1190 INTEGER :: i, j, nr, nu, pot_type
1191 REAL(kind=
dp) :: scale_coulomb, scale_longrange
1192 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:) :: abscissa, weights
1196 IF (
ASSOCIATED(
atom%xc_section))
THEN
1197 xc_section =>
atom%xc_section
1202 atom%hfx_pot%scale_longrange = 0.0_dp
1203 atom%hfx_pot%scale_coulomb = 1.0_dp
1216 SELECT CASE (pot_type)
1218 cpwarn(
"Potential not implemented, use Coulomb instead!")
1220 atom%hfx_pot%scale_longrange = 0.0_dp
1221 atom%hfx_pot%scale_coulomb = scale_coulomb
1223 atom%hfx_pot%scale_coulomb = 0.0_dp
1224 atom%hfx_pot%scale_longrange = scale_longrange
1226 atom%hfx_pot%scale_coulomb = 1.0_dp
1227 atom%hfx_pot%scale_longrange = -1.0_dp
1229 atom%hfx_pot%scale_coulomb = scale_coulomb
1230 atom%hfx_pot%scale_longrange = scale_longrange
1236 cpabort(
"Only numerical and semi-analytic lrHF exchange available!")
1239 IF (
atom%hfx_pot%scale_longrange /= 0.0_dp .AND. extype ==
do_numeric .AND. .NOT.
ALLOCATED(
atom%hfx_pot%kernel))
THEN
1242 IF (
atom%hfx_pot%do_gh)
THEN
1246 ALLOCATE (weights(
atom%hfx_pot%nr_gh), abscissa(
atom%hfx_pot%nr_gh))
1247 CALL get_gauss_hermite_weights(abscissa, weights,
atom%hfx_pot%nr_gh)
1249 nr =
atom%basis%grid%nr
1250 ALLOCATE (
atom%hfx_pot%kernel(nr,
atom%hfx_pot%nr_gh, 0:
atom%state%maxl_calc +
atom%state%maxl_occ))
1251 atom%hfx_pot%kernel = 0.0_dp
1252 DO nu = 0,
atom%state%maxl_calc +
atom%state%maxl_occ
1253 DO i = 1,
atom%hfx_pot%nr_gh
1256 *abscissa(i)*
atom%basis%grid%rad(j), nu)*sqrt(weights(i))
1264 nr =
atom%basis%grid%nr
1265 ALLOCATE (
atom%hfx_pot%kernel(nr, nr, 0:
atom%state%maxl_calc +
atom%state%maxl_occ))
1266 atom%hfx_pot%kernel = 0.0_dp
1267 DO nu = 0,
atom%state%maxl_calc +
atom%state%maxl_occ
1271 *
atom%basis%grid%rad(i)*
atom%basis%grid%rad(j), nu)
1278 NULLIFY (xc_section)
1290 SUBROUTINE get_gauss_hermite_weights(abscissa, weights, nn)
1291 REAL(kind=
dp),
DIMENSION(:),
INTENT(OUT) :: abscissa, weights
1292 INTEGER,
INTENT(IN) :: nn
1294 INTEGER :: counter, ii, info, liwork, lwork
1295 INTEGER,
ALLOCATABLE,
DIMENSION(:) :: iwork
1296 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:) :: diag, subdiag, work
1297 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:, :) :: eigenvec
1301 ALLOCATE (work(1), iwork(1), diag(2*nn), subdiag(2*nn - 1), eigenvec(2*nn, 2*nn))
1306 subdiag(ii) = sqrt(real(ii, kind=
dp)/2.0_dp)
1310 CALL dstevd(
'V', 2*nn, diag, subdiag, eigenvec, 2*nn, work, lwork, iwork, liwork, info)
1313 cpabort(
'Finding size of working matrices failed!')
1317 lwork = int(work(1))
1319 DEALLOCATE (work, iwork)
1320 ALLOCATE (work(lwork), iwork(liwork))
1323 CALL dstevd(
'V', 2*nn, diag, subdiag, eigenvec, 2*nn, work, lwork, iwork, liwork, info)
1326 cpabort(
'Eigenvalue decomposition failed!')
1329 DEALLOCATE (work, iwork, subdiag)
1335 IF (diag(ii) > 0.0_dp)
THEN
1336 counter = counter + 1
1337 abscissa(counter) = diag(ii)
1338 weights(counter) =
rootpi*eigenvec(1, ii)**2
1341 IF (counter /= nn)
THEN
1342 cpabort(
'Have not found enough or too many zeros!')
1345 END SUBROUTINE get_gauss_hermite_weights
1355 INTEGER,
DIMENSION(0:lmat),
INTENT(IN) :: n
1356 INTEGER,
INTENT(IN),
OPTIONAL :: lmax
1361 IF (
PRESENT(lmax))
THEN
1367 cpassert(.NOT.
ASSOCIATED(opmat))
1372 ALLOCATE (opmat%op(m, m, 0:lm))
1384 cpassert(
ASSOCIATED(opmat))
1387 DEALLOCATE (opmat%op)
1404 cpassert(.NOT.
ASSOCIATED(opgrid))
1412 ALLOCATE (opgrid%op(nr))
1424 cpassert(
ASSOCIATED(opgrid))
1426 NULLIFY (opgrid%grid)
1427 DEALLOCATE (opgrid%op)
1442 INTEGER,
INTENT(IN) :: zval
1443 REAL(
dp),
INTENT(OUT) :: cval, aval
1444 INTEGER,
DIMENSION(0:lmat),
INTENT(OUT) :: ngto, ival
1459 cval = 2.14774520_dp
1460 aval = 0.04850670_dp
1463 cval = 2.08932430_dp
1464 aval = 0.02031060_dp
1467 cval = 2.09753060_dp
1468 aval = 0.03207070_dp
1471 cval = 2.10343410_dp
1472 aval = 0.03591970_dp
1476 cval = 2.10662820_dp
1477 aval = 0.05292410_dp
1481 cval = 2.13743840_dp
1482 aval = 0.06291970_dp
1486 cval = 2.08687310_dp
1487 aval = 0.08350860_dp
1491 cval = 2.12318180_dp
1492 aval = 0.09899170_dp
1496 cval = 2.13164810_dp
1497 aval = 0.11485350_dp
1501 cval = 2.11413310_dp
1502 aval = 0.00922630_dp
1507 cval = 2.12183620_dp
1508 aval = 0.01215850_dp
1513 cval = 2.06073230_dp
1514 aval = 0.01449350_dp
1519 cval = 2.08563660_dp
1520 aval = 0.01861460_dp
1525 cval = 2.04879270_dp
1526 aval = 0.02147790_dp
1531 cval = 2.06216660_dp
1532 aval = 0.01978920_dp
1537 cval = 2.04628670_dp
1538 aval = 0.02451470_dp
1543 cval = 2.08675200_dp
1544 aval = 0.02635040_dp
1549 cval = 2.02715220_dp
1550 aval = 0.01822040_dp
1555 cval = 2.01465650_dp
1556 aval = 0.01646570_dp
1561 cval = 2.01605240_dp
1562 aval = 0.01254190_dp
1568 cval = 2.01800000_dp
1569 aval = 0.01195490_dp
1576 cval = 1.98803560_dp
1577 aval = 0.02492140_dp
1584 cval = 1.98984000_dp
1585 aval = 0.02568400_dp
1592 cval = 2.01694380_dp
1593 aval = 0.02664480_dp
1600 cval = 2.01824090_dp
1601 aval = 0.01355000_dp
1608 cval = 1.98359400_dp
1609 aval = 0.01702210_dp
1616 cval = 1.96797340_dp
1617 aval = 0.02163180_dp
1624 cval = 1.98955180_dp
1625 aval = 0.02304480_dp
1632 cval = 1.98074320_dp
1633 aval = 0.02754320_dp
1640 cval = 2.00551070_dp
1641 aval = 0.02005530_dp
1648 cval = 2.00000030_dp
1649 aval = 0.02003000_dp
1656 cval = 2.00609100_dp
1657 aval = 0.02055620_dp
1664 cval = 2.00701000_dp
1665 aval = 0.02230400_dp
1672 cval = 2.01508710_dp
1673 aval = 0.02685790_dp
1680 cval = 2.01960430_dp
1681 aval = 0.02960430_dp
1688 cval = 2.00031000_dp
1689 aval = 0.00768400_dp
1697 cval = 1.99563960_dp
1698 aval = 0.01401940_dp
1705 cval = 1.98971210_dp
1706 aval = 0.01558470_dp
1712 cval = 1.97976190_dp
1713 aval = 0.01705520_dp
1719 cval = 1.97989290_dp
1720 aval = 0.01527040_dp
1726 cval = 1.97909240_dp
1727 aval = 0.01879720_dp
1733 cval = 1.98508430_dp
1734 aval = 0.01497550_dp
1741 cval = 1.98515010_dp
1742 aval = 0.01856670_dp
1749 cval = 1.98502970_dp
1750 aval = 0.01487000_dp
1757 cval = 1.97672850_dp
1758 aval = 0.01762500_dp
1766 cval = 1.97862730_dp
1767 aval = 0.01863310_dp
1774 cval = 1.97990020_dp
1775 aval = 0.01347150_dp
1782 cval = 1.97979410_dp
1783 aval = 0.00890265_dp
1790 cval = 1.98001000_dp
1791 aval = 0.00895215_dp
1798 cval = 1.97979980_dp
1799 aval = 0.01490290_dp
1806 cval = 1.98009310_dp
1807 aval = 0.01490390_dp
1814 cval = 1.97794750_dp
1815 aval = 0.01425880_dp
1823 cval = 1.97784450_dp
1824 aval = 0.01430130_dp
1832 cval = 1.97784450_dp
1833 aval = 0.00499318_dp
1841 cval = 1.97764820_dp
1842 aval = 0.00500392_dp
1850 cval = 1.97765150_dp
1851 aval = 0.00557083_dp
1859 cval = 1.97768750_dp
1860 aval = 0.00547531_dp
1870 cval = 1.96986600_dp
1871 aval = 0.00813143_dp
1881 cval = 1.97765720_dp
1882 aval = 0.00489201_dp
1892 cval = 1.97768120_dp
1893 aval = 0.00499000_dp
1903 cval = 1.97745700_dp
1904 aval = 0.00615587_dp
1914 cval = 1.97570240_dp
1915 aval = 0.00769959_dp
1925 cval = 1.97629350_dp
1926 aval = 0.00706610_dp
1936 cval = 1.96900000_dp
1937 aval = 0.01019150_dp
1947 cval = 1.97350000_dp
1948 aval = 0.01334320_dp
1958 cval = 1.97493000_dp
1959 aval = 0.01331360_dp
1968 cval = 1.97597670_dp
1969 aval = 0.01434040_dp
1979 cval = 1.97809240_dp
1980 aval = 0.01529430_dp
1990 cval = 1.97644360_dp
1991 aval = 0.01312770_dp
2001 cval = 1.96998000_dp
2002 aval = 0.01745150_dp
2012 cval = 1.97223830_dp
2013 aval = 0.01639750_dp
2023 cval = 1.97462110_dp
2024 aval = 0.01603680_dp
2034 cval = 1.97756000_dp
2035 aval = 0.02030570_dp
2045 cval = 1.97645760_dp
2046 aval = 0.02057180_dp
2056 cval = 1.97725820_dp
2057 aval = 0.02058210_dp
2067 cval = 1.97749380_dp
2068 aval = 0.02219380_dp
2078 cval = 1.97946280_dp
2079 aval = 0.02216280_dp
2089 cval = 1.97852130_dp
2090 aval = 0.02168500_dp
2100 cval = 1.98045190_dp
2101 aval = 0.02177860_dp
2111 cval = 1.97000000_dp
2112 aval = 0.02275000_dp
2122 cval = 1.97713580_dp
2123 aval = 0.02317030_dp
2133 cval = 1.97537880_dp
2134 aval = 0.02672860_dp
2144 cval = 1.97545360_dp
2145 aval = 0.02745360_dp
2155 cval = 1.97338370_dp
2156 aval = 0.02616310_dp
2166 cval = 1.97294240_dp
2167 aval = 0.02429220_dp
2177 cval = 1.98000000_dp
2178 aval = 0.01400000_dp
2198 SUBROUTINE read_basis_set(element_symbol, basis, basis_set_name, basis_set_file, &
2201 CHARACTER(LEN=*),
INTENT(IN) :: element_symbol
2203 CHARACTER(LEN=*),
INTENT(IN) :: basis_set_name, basis_set_file
2206 INTEGER,
PARAMETER :: maxpri = 40, maxset = 20
2208 CHARACTER(len=20*default_string_length) :: line_att
2209 CHARACTER(LEN=240) :: line
2210 CHARACTER(LEN=242) :: line2
2211 CHARACTER(LEN=LEN(basis_set_name)) :: bsname
2212 CHARACTER(LEN=LEN(basis_set_name)+2) :: bsname2
2213 CHARACTER(LEN=LEN(element_symbol)) :: symbol
2214 CHARACTER(LEN=LEN(element_symbol)+2) :: symbol2
2215 INTEGER :: i, ii, ipgf, irep, iset, ishell, j, k, &
2216 lshell, nj, nmin, ns, nset, strlen1, &
2218 INTEGER,
DIMENSION(maxpri, maxset) :: l
2219 INTEGER,
DIMENSION(maxset) :: lmax, lmin, n, npgf, nshell
2220 LOGICAL :: found, is_ok, match, read_from_input
2221 REAL(
dp) :: expzet, gcca, prefac, zeta
2222 REAL(
dp),
DIMENSION(maxpri, maxpri, maxset) :: gcc
2223 REAL(
dp),
DIMENSION(maxpri, maxset) :: zet
2227 bsname = basis_set_name
2228 symbol = element_symbol
2240 read_from_input = .false.
2242 IF (read_from_input)
THEN
2249 CALL val_get(val, c_val=line_att)
2250 READ (line_att, *) nset
2251 cpassert(nset <= maxset)
2255 CALL val_get(val, c_val=line_att)
2256 READ (line_att, *) n(iset)
2258 READ (line_att, *) lmin(iset)
2260 READ (line_att, *) lmax(iset)
2262 READ (line_att, *) npgf(iset)
2264 cpassert(npgf(iset) <= maxpri)
2266 DO lshell = lmin(iset), lmax(iset)
2267 nmin = n(iset) + lshell - lmin(iset)
2268 READ (line_att, *) ishell
2270 nshell(iset) = nshell(iset) + ishell
2272 l(nshell(iset) - ishell + i, iset) = lshell
2275 cpassert(len_trim(line_att) == 0)
2276 DO ipgf = 1, npgf(iset)
2279 CALL val_get(val, c_val=line_att)
2280 READ (line_att, *) zet(ipgf, iset), (gcc(ipgf, ishell, iset), ishell=1, nshell(iset))
2297 line2 =
" "//line//
" "
2298 symbol2 =
" "//trim(symbol)//
" "
2299 bsname2 =
" "//trim(bsname)//
" "
2300 strlen1 = len_trim(symbol2) + 1
2301 strlen2 = len_trim(bsname2) + 1
2303 IF ((index(line2, symbol2(:strlen1)) > 0) .AND. &
2304 (index(line2, bsname2(:strlen2)) > 0)) match = .true.
2309 cpassert(nset <= maxset)
2315 cpassert(npgf(iset) <= maxpri)
2317 DO lshell = lmin(iset), lmax(iset)
2318 nmin = n(iset) + lshell - lmin(iset)
2320 nshell(iset) = nshell(iset) + ishell
2322 l(nshell(iset) - ishell + i, iset) = lshell
2325 DO ipgf = 1, npgf(iset)
2327 DO ishell = 1, nshell(iset)
2351 DO j = lmin(i), min(lmax(i),
lmat)
2352 basis%nprim(j) = basis%nprim(j) + npgf(i)
2356 IF (k <=
lmat) basis%nbas(k) = basis%nbas(k) + 1
2360 nj = maxval(basis%nprim)
2361 ns = maxval(basis%nbas)
2362 ALLOCATE (basis%am(nj, 0:
lmat))
2364 ALLOCATE (basis%cm(nj, ns, 0:
lmat))
2371 IF (j >= lmin(i) .AND. j <= lmax(i))
THEN
2372 DO ipgf = 1, npgf(i)
2373 basis%am(nj + ipgf, j) = zet(ipgf, i)
2375 DO ii = 1, nshell(i)
2376 IF (l(ii, i) == j)
THEN
2378 DO ipgf = 1, npgf(i)
2379 basis%cm(nj + ipgf, ns, j) = gcc(ipgf, ii, i)
2390 expzet = 0.25_dp*real(2*j + 3,
dp)
2391 prefac = sqrt(
rootpi/2._dp**(j + 2)*
dfac(2*j + 1))
2392 DO ipgf = 1, basis%nprim(j)
2393 DO ii = 1, basis%nbas(j)
2394 gcca = basis%cm(ipgf, ii, j)
2395 zeta = 2._dp*basis%am(ipgf, j)
2396 basis%cm(ipgf, ii, j) = zeta**expzet*gcca/prefac
2401 END SUBROUTINE read_basis_set
2410 TYPE(section_vals_type),
POINTER :: opt_section
2412 INTEGER :: miter, ndiis
2413 REAL(kind=dp) :: damp, eps_diis, eps_scf
2415 CALL section_vals_val_get(opt_section,
"MAX_ITER", i_val=miter)
2416 CALL section_vals_val_get(opt_section,
"EPS_SCF", r_val=eps_scf)
2417 CALL section_vals_val_get(opt_section,
"N_DIIS", i_val=ndiis)
2418 CALL section_vals_val_get(opt_section,
"EPS_DIIS", r_val=eps_diis)
2419 CALL section_vals_val_get(opt_section,
"DAMPING", r_val=damp)
2421 optimization%max_iter = miter
2422 optimization%eps_scf = eps_scf
2423 optimization%n_diis = ndiis
2424 optimization%eps_diis = eps_diis
2425 optimization%damping = damp
2436 TYPE(section_vals_type),
POINTER :: potential_section
2437 INTEGER,
INTENT(IN) :: zval
2439 CHARACTER(LEN=default_string_length) :: pseudo_fn, pseudo_name
2441 REAL(dp),
DIMENSION(:),
POINTER :: convals
2442 TYPE(section_vals_type),
POINTER :: ecp_potential_section, &
2443 gth_potential_section
2446 CALL section_vals_val_get(potential_section,
"PSEUDO_TYPE", i_val=potential%ppot_type)
2448 SELECT CASE (potential%ppot_type)
2450 CALL section_vals_val_get(potential_section,
"POTENTIAL_FILE_NAME", c_val=pseudo_fn)
2451 CALL section_vals_val_get(potential_section,
"POTENTIAL_NAME", c_val=pseudo_name)
2452 gth_potential_section => section_vals_get_subs_vals(potential_section,
"GTH_POTENTIAL")
2453 CALL read_gth_potential(ptable(zval)%symbol, potential%gth_pot, &
2454 pseudo_name, pseudo_fn, gth_potential_section)
2456 CALL section_vals_val_get(potential_section,
"POTENTIAL_FILE_NAME", c_val=pseudo_fn)
2457 CALL section_vals_val_get(potential_section,
"POTENTIAL_NAME", c_val=pseudo_name)
2458 ecp_potential_section => section_vals_get_subs_vals(potential_section,
"ECP")
2460 pseudo_name, pseudo_fn, ecp_potential_section)
2462 CALL section_vals_val_get(potential_section,
"POTENTIAL_FILE_NAME", c_val=pseudo_fn)
2463 CALL section_vals_val_get(potential_section,
"POTENTIAL_NAME", c_val=pseudo_name)
2464 CALL atom_read_upf(potential%upf_pot, pseudo_fn)
2465 potential%upf_pot%pname = pseudo_name
2467 cpabort(
"Not implemented")
2474 potential%ppot_type = no_pseudo
2479 CALL section_vals_val_get(potential_section,
"CONFINEMENT_TYPE", i_val=ic)
2480 potential%conf_type = ic
2481 IF (potential%conf_type == no_conf)
THEN
2482 potential%acon = 0.0_dp
2483 potential%rcon = 4.0_dp
2484 potential%scon = 2.0_dp
2485 potential%confinement = .false.
2486 ELSE IF (potential%conf_type == poly_conf)
THEN
2487 CALL section_vals_val_get(potential_section,
"CONFINEMENT", r_vals=convals)
2488 IF (
SIZE(convals) >= 1)
THEN
2489 IF (convals(1) > 0.0_dp)
THEN
2490 potential%confinement = .true.
2491 potential%acon = convals(1)
2492 IF (
SIZE(convals) >= 2)
THEN
2493 potential%rcon = convals(2)
2495 potential%rcon = 4.0_dp
2497 IF (
SIZE(convals) >= 3)
THEN
2498 potential%scon = convals(3)
2500 potential%scon = 2.0_dp
2503 potential%confinement = .false.
2506 potential%confinement = .false.
2508 ELSE IF (potential%conf_type == barrier_conf)
THEN
2509 potential%acon = 200.0_dp
2510 potential%rcon = 4.0_dp
2511 potential%scon = 12.0_dp
2512 potential%confinement = .true.
2513 CALL section_vals_val_get(potential_section,
"CONFINEMENT", r_vals=convals)
2514 IF (
SIZE(convals) >= 1)
THEN
2515 IF (convals(1) > 0.0_dp)
THEN
2516 potential%acon = convals(1)
2517 IF (
SIZE(convals) >= 2)
THEN
2518 potential%rcon = convals(2)
2520 IF (
SIZE(convals) >= 3)
THEN
2521 potential%scon = convals(3)
2524 potential%confinement = .false.
2537 potential%confinement = .false.
2539 CALL atom_release_upf(potential%upf_pot)
2550 SUBROUTINE read_gth_potential(element_symbol, potential, pseudo_name, pseudo_file, &
2553 CHARACTER(LEN=*),
INTENT(IN) :: element_symbol
2555 CHARACTER(LEN=*),
INTENT(IN) :: pseudo_name, pseudo_file
2556 TYPE(section_vals_type),
POINTER :: potential_section
2558 CHARACTER(LEN=240) :: line
2559 CHARACTER(LEN=242) :: line2
2560 CHARACTER(len=5*default_string_length) :: line_att
2561 CHARACTER(LEN=LEN(element_symbol)) :: symbol
2562 CHARACTER(LEN=LEN(element_symbol)+2) :: symbol2
2563 CHARACTER(LEN=LEN(pseudo_name)) :: apname
2564 CHARACTER(LEN=LEN(pseudo_name)+2) :: apname2
2565 INTEGER :: i, ic, ipot, j, l, nlmax, strlen1, &
2567 INTEGER,
DIMENSION(0:lmat) :: elec_conf
2568 LOGICAL :: found, is_ok, match, read_from_input
2569 TYPE(cp_sll_val_type),
POINTER ::
list
2570 TYPE(val_type),
POINTER :: val
2574 apname = pseudo_name
2575 symbol = element_symbol
2577 potential%symbol = symbol
2578 potential%pname = apname
2580 potential%rc = 0._dp
2582 potential%cl = 0._dp
2584 potential%rcnl = 0._dp
2585 potential%hnl = 0._dp
2586 potential%soc = .false.
2587 potential%knl = 0._dp
2589 potential%lpotextended = .false.
2590 potential%lsdpot = .false.
2591 potential%nlcc = .false.
2592 potential%nexp_lpot = 0
2593 potential%nexp_lsd = 0
2594 potential%nexp_nlcc = 0
2596 read_from_input = .false.
2597 CALL section_vals_get(potential_section, explicit=read_from_input)
2598 IF (read_from_input)
THEN
2599 CALL section_vals_list_get(potential_section,
"_DEFAULT_KEYWORD_",
list=
list)
2600 CALL uppercase(symbol)
2601 CALL uppercase(apname)
2604 is_ok = cp_sll_val_next(
list, val)
2606 CALL val_get(val, c_val=line_att)
2607 READ (line_att, *) elec_conf(l)
2608 CALL remove_word(line_att)
2609 DO WHILE (len_trim(line_att) /= 0)
2611 READ (line_att, *) elec_conf(l)
2612 CALL remove_word(line_att)
2614 potential%econf(0:
lmat) = elec_conf(0:
lmat)
2615 potential%zion = real(sum(elec_conf), dp)
2617 is_ok = cp_sll_val_next(
list, val)
2619 CALL val_get(val, c_val=line_att)
2620 READ (line_att, *) potential%rc
2621 CALL remove_word(line_att)
2623 READ (line_att, *) potential%ncl
2624 CALL remove_word(line_att)
2625 DO i = 1, potential%ncl
2626 READ (line_att, *) potential%cl(i)
2627 CALL remove_word(line_att)
2631 is_ok = cp_sll_val_next(
list, val)
2633 CALL val_get(val, c_val=line_att)
2634 IF (index(line_att,
"LPOT") /= 0)
THEN
2635 potential%lpotextended = .true.
2636 CALL remove_word(line_att)
2637 READ (line_att, *) potential%nexp_lpot
2638 DO ipot = 1, potential%nexp_lpot
2639 is_ok = cp_sll_val_next(
list, val)
2641 CALL val_get(val, c_val=line_att)
2642 READ (line_att, *) potential%alpha_lpot(ipot)
2643 CALL remove_word(line_att)
2644 READ (line_att, *) potential%nct_lpot(ipot)
2645 CALL remove_word(line_att)
2646 DO ic = 1, potential%nct_lpot(ipot)
2647 READ (line_att, *) potential%cval_lpot(ic, ipot)
2648 CALL remove_word(line_att)
2651 ELSEIF (index(line_att,
"NLCC") /= 0)
THEN
2652 potential%nlcc = .true.
2653 CALL remove_word(line_att)
2654 READ (line_att, *) potential%nexp_nlcc
2655 DO ipot = 1, potential%nexp_nlcc
2656 is_ok = cp_sll_val_next(
list, val)
2658 CALL val_get(val, c_val=line_att)
2659 READ (line_att, *) potential%alpha_nlcc(ipot)
2660 CALL remove_word(line_att)
2661 READ (line_att, *) potential%nct_nlcc(ipot)
2662 CALL remove_word(line_att)
2663 DO ic = 1, potential%nct_nlcc(ipot)
2664 READ (line_att, *) potential%cval_nlcc(ic, ipot)
2666 potential%cval_nlcc(ic, ipot) = potential%cval_nlcc(ic, ipot)/(4.0_dp*pi)
2667 CALL remove_word(line_att)
2670 ELSEIF (index(line_att,
"LSD") /= 0)
THEN
2671 potential%lsdpot = .true.
2672 CALL remove_word(line_att)
2673 READ (line_att, *) potential%nexp_lsd
2674 DO ipot = 1, potential%nexp_lsd
2675 is_ok = cp_sll_val_next(
list, val)
2677 CALL val_get(val, c_val=line_att)
2678 READ (line_att, *) potential%alpha_lsd(ipot)
2679 CALL remove_word(line_att)
2680 READ (line_att, *) potential%nct_lsd(ipot)
2681 CALL remove_word(line_att)
2682 DO ic = 1, potential%nct_lsd(ipot)
2683 READ (line_att, *) potential%cval_lsd(ic, ipot)
2684 CALL remove_word(line_att)
2692 READ (line_att, *) nlmax
2693 CALL remove_word(line_att)
2694 IF (index(line_att,
"SOC") /= 0) potential%soc = .true.
2698 is_ok = cp_sll_val_next(
list, val)
2700 CALL val_get(val, c_val=line_att)
2701 READ (line_att, *) potential%rcnl(l)
2702 CALL remove_word(line_att)
2703 READ (line_att, *) potential%nl(l)
2704 CALL remove_word(line_att)
2705 DO i = 1, potential%nl(l)
2707 READ (line_att, *) potential%hnl(1, 1, l)
2708 CALL remove_word(line_att)
2710 cpassert(len_trim(line_att) == 0)
2711 is_ok = cp_sll_val_next(
list, val)
2713 CALL val_get(val, c_val=line_att)
2714 READ (line_att, *) potential%hnl(i, i, l)
2715 CALL remove_word(line_att)
2717 DO j = i + 1, potential%nl(l)
2718 READ (line_att, *) potential%hnl(i, j, l)
2719 potential%hnl(j, i, l) = potential%hnl(i, j, l)
2720 CALL remove_word(line_att)
2723 IF (potential%soc .AND. l /= 0)
THEN
2724 is_ok = cp_sll_val_next(
list, val)
2726 CALL val_get(val, c_val=line_att)
2727 DO i = 1, potential%nl(l)
2729 READ (line_att, *) potential%knl(1, 1, l)
2730 CALL remove_word(line_att)
2732 cpassert(len_trim(line_att) == 0)
2733 is_ok = cp_sll_val_next(
list, val)
2735 CALL val_get(val, c_val=line_att)
2736 READ (line_att, *) potential%knl(i, i, l)
2737 CALL remove_word(line_att)
2739 DO j = i + 1, potential%nl(l)
2740 READ (line_att, *) potential%knl(i, j, l)
2741 potential%knl(j, i, l) = potential%knl(i, j, l)
2742 CALL remove_word(line_att)
2746 cpassert(len_trim(line_att) == 0)
2751 TYPE(cp_parser_type) :: parser
2752 CALL parser_create(parser, pseudo_file)
2755 CALL parser_search_string(parser, trim(apname), .true., found, line)
2757 CALL uppercase(symbol)
2758 CALL uppercase(apname)
2761 CALL uppercase(line)
2762 line2 =
" "//line//
" "
2763 symbol2 =
" "//trim(symbol)//
" "
2764 apname2 =
" "//trim(apname)//
" "
2765 strlen1 = len_trim(symbol2) + 1
2766 strlen2 = len_trim(apname2) + 1
2768 IF ((index(line2, symbol2(:strlen1)) > 0) .AND. &
2769 (index(line2, apname2(:strlen2)) > 0)) match = .true.
2774 CALL parser_get_object(parser, elec_conf(l), newline=.true.)
2775 DO WHILE (parser_test_next_token(parser) ==
"INT")
2777 CALL parser_get_object(parser, elec_conf(l))
2779 potential%econf(0:
lmat) = elec_conf(0:
lmat)
2780 potential%zion = real(sum(elec_conf), dp)
2782 CALL parser_get_object(parser, potential%rc, newline=.true.)
2784 CALL parser_get_object(parser, potential%ncl)
2785 DO i = 1, potential%ncl
2786 CALL parser_get_object(parser, potential%cl(i))
2790 CALL parser_get_next_line(parser, 1)
2791 IF (parser_test_next_token(parser) ==
"INT")
THEN
2793 ELSEIF (parser_test_next_token(parser) ==
"STR")
THEN
2794 CALL parser_get_object(parser, line)
2795 IF (index(line,
"LPOT") /= 0)
THEN
2797 potential%lpotextended = .true.
2798 CALL parser_get_object(parser, potential%nexp_lpot)
2799 DO ipot = 1, potential%nexp_lpot
2800 CALL parser_get_object(parser, potential%alpha_lpot(ipot), newline=.true.)
2801 CALL parser_get_object(parser, potential%nct_lpot(ipot))
2802 DO ic = 1, potential%nct_lpot(ipot)
2803 CALL parser_get_object(parser, potential%cval_lpot(ic, ipot))
2806 ELSEIF (index(line,
"NLCC") /= 0)
THEN
2808 potential%nlcc = .true.
2809 CALL parser_get_object(parser, potential%nexp_nlcc)
2810 DO ipot = 1, potential%nexp_nlcc
2811 CALL parser_get_object(parser, potential%alpha_nlcc(ipot), newline=.true.)
2812 CALL parser_get_object(parser, potential%nct_nlcc(ipot))
2813 DO ic = 1, potential%nct_nlcc(ipot)
2814 CALL parser_get_object(parser, potential%cval_nlcc(ic, ipot))
2816 potential%cval_nlcc(ic, ipot) = potential%cval_nlcc(ic, ipot)/(4.0_dp*pi)
2819 ELSEIF (index(line,
"LSD") /= 0)
THEN
2821 potential%lsdpot = .true.
2822 CALL parser_get_object(parser, potential%nexp_lsd)
2823 DO ipot = 1, potential%nexp_lsd
2824 CALL parser_get_object(parser, potential%alpha_lsd(ipot), newline=.true.)
2825 CALL parser_get_object(parser, potential%nct_lsd(ipot))
2826 DO ic = 1, potential%nct_lsd(ipot)
2827 CALL parser_get_object(parser, potential%cval_lsd(ic, ipot))
2838 CALL parser_get_object(parser, nlmax)
2840 IF (parser_test_next_token(parser) ==
"STR")
THEN
2841 CALL parser_get_object(parser, line)
2842 IF (index(line,
"SOC") /= 0) potential%soc = .true.
2846 CALL parser_get_object(parser, potential%rcnl(l), newline=.true.)
2847 CALL parser_get_object(parser, potential%nl(l))
2848 DO i = 1, potential%nl(l)
2850 CALL parser_get_object(parser, potential%hnl(i, i, l))
2852 CALL parser_get_object(parser, potential%hnl(i, i, l), newline=.true.)
2854 DO j = i + 1, potential%nl(l)
2855 CALL parser_get_object(parser, potential%hnl(i, j, l))
2856 potential%hnl(j, i, l) = potential%hnl(i, j, l)
2859 IF (potential%soc .AND. l /= 0)
THEN
2860 DO i = 1, potential%nl(l)
2861 CALL parser_get_object(parser, potential%knl(i, i, l), newline=.true.)
2862 DO j = i + 1, potential%nl(l)
2863 CALL parser_get_object(parser, potential%knl(i, j, l))
2864 potential%knl(j, i, l) = potential%knl(i, j, l)
2879 CALL parser_release(parser)
2883 END SUBROUTINE read_gth_potential
2895 CHARACTER(LEN=*),
INTENT(IN) :: element_symbol
2897 CHARACTER(LEN=*),
INTENT(IN) :: pseudo_name, pseudo_file
2898 TYPE(section_vals_type),
POINTER :: potential_section
2900 CHARACTER(LEN=240) :: line
2901 CHARACTER(len=5*default_string_length) :: line_att
2902 CHARACTER(LEN=LEN(element_symbol)+1) :: symbol
2903 CHARACTER(LEN=LEN(pseudo_name)) :: apname
2904 INTEGER :: i, ic, l, ncore, nel
2905 LOGICAL :: found, is_ok, read_from_input
2906 TYPE(cp_sll_val_type),
POINTER ::
list
2907 TYPE(val_type),
POINTER :: val
2909 apname = pseudo_name
2910 symbol = element_symbol
2911 CALL get_ptable_info(symbol, number=ncore)
2913 potential%symbol = symbol
2914 potential%pname = apname
2920 potential%aloc = 0.0_dp
2921 potential%bloc = 0.0_dp
2924 potential%apot = 0.0_dp
2925 potential%bpot = 0.0_dp
2927 read_from_input = .false.
2928 CALL section_vals_get(potential_section, explicit=read_from_input)
2929 IF (read_from_input)
THEN
2930 CALL section_vals_list_get(potential_section,
"_DEFAULT_KEYWORD_",
list=
list)
2932 is_ok = cp_sll_val_next(
list, val)
2934 CALL val_get(val, c_val=line_att)
2935 CALL remove_word(line_att)
2936 CALL remove_word(line_att)
2938 READ (line_att, *) nel
2939 potential%zion = real(ncore - nel, kind=dp)
2941 is_ok = cp_sll_val_next(
list, val)
2943 CALL val_get(val, c_val=line_att)
2945 IF (.NOT. cp_sll_val_next(
list, val))
EXIT
2946 CALL val_get(val, c_val=line_att)
2947 IF (index(line_att, element_symbol) == 0)
THEN
2948 potential%nloc = potential%nloc + 1
2950 READ (line_att, *) potential%nrloc(ic), potential%bloc(ic), potential%aloc(ic)
2957 CALL val_get(val, c_val=line_att)
2958 IF (index(line_att, element_symbol) == 0)
THEN
2959 potential%npot(l) = potential%npot(l) + 1
2960 ic = potential%npot(l)
2961 READ (line_att, *) potential%nrpot(ic, l), potential%bpot(ic, l), potential%apot(ic, l)
2963 potential%lmax = potential%lmax + 1
2966 IF (.NOT. cp_sll_val_next(
list, val))
EXIT
2971 TYPE(cp_parser_type) :: parser
2972 CALL parser_create(parser, pseudo_file)
2975 CALL parser_search_string(parser, trim(apname), .true., found, line)
2978 CALL parser_get_object(parser, line, newline=.true.)
2979 IF (trim(line) == element_symbol)
THEN
2980 CALL parser_get_object(parser, line, lower_to_upper=.true.)
2981 cpassert(trim(line) ==
"NELEC")
2983 CALL parser_get_object(parser, nel)
2984 potential%zion = real(ncore - nel, kind=dp)
2986 CALL parser_get_object(parser, line, newline=.true.)
2989 CALL parser_read_line(parser, 1)
2990 IF (parser_test_next_token(parser) ==
"STR")
EXIT
2991 potential%nloc = potential%nloc + 1
2993 CALL parser_get_object(parser, potential%nrloc(ic))
2994 CALL parser_get_object(parser, potential%bloc(ic))
2995 CALL parser_get_object(parser, potential%aloc(ic))
2999 CALL parser_get_object(parser, symbol)
3000 IF (symbol == element_symbol)
THEN
3002 potential%lmax = potential%lmax + 1
3004 CALL parser_read_line(parser, 1)
3005 IF (parser_test_next_token(parser) ==
"STR")
EXIT
3006 potential%npot(l) = potential%npot(l) + 1
3007 ic = potential%npot(l)
3008 CALL parser_get_object(parser, potential%nrpot(ic, l))
3009 CALL parser_get_object(parser, potential%bpot(ic, l))
3010 CALL parser_get_object(parser, potential%apot(ic, l))
3017 ELSE IF (line ==
"END")
THEN
3018 cpabort(
"Element not found in ECP library")
3022 cpabort(
"ECP type not found in library")
3027 CALL parser_release(parser)
3032 potential%econf(0:3) = ptable(ncore)%e_conv(0:3)
3035 cpabort(
"Unknown Core State")
3038 potential%econf(0:3) = potential%econf(0:3) - ptable(2)%e_conv(0:3)
3040 potential%econf(0:3) = potential%econf(0:3) - ptable(10)%e_conv(0:3)
3042 potential%econf(0:3) = potential%econf(0:3) - ptable(18)%e_conv(0:3)
3044 potential%econf(0:3) = potential%econf(0:3) - ptable(18)%e_conv(0:3)
3045 potential%econf(2) = potential%econf(2) - 10
3047 potential%econf(0:3) = potential%econf(0:3) - ptable(36)%e_conv(0:3)
3049 potential%econf(0:3) = potential%econf(0:3) - ptable(36)%e_conv(0:3)
3050 potential%econf(2) = potential%econf(2) - 10
3052 potential%econf(0:3) = potential%econf(0:3) - ptable(54)%e_conv(0:3)
3054 potential%econf(0:3) = potential%econf(0:3) - ptable(36)%e_conv(0:3)
3055 potential%econf(2) = potential%econf(2) - 10
3056 potential%econf(3) = potential%econf(3) - 14
3058 potential%econf(0:3) = potential%econf(0:3) - ptable(54)%e_conv(0:3)
3059 potential%econf(3) = potential%econf(3) - 14
3061 potential%econf(0:3) = potential%econf(0:3) - ptable(54)%e_conv(0:3)
3062 potential%econf(2) = potential%econf(2) - 10
3063 potential%econf(3) = potential%econf(3) - 14
3066 cpassert(all(potential%econf >= 0))
3076 TYPE(grid_atom_type) :: grid1, grid2
3080 REAL(kind=dp) :: dr, dw
3083 IF (grid1%nr == grid2%nr)
THEN
3085 dr = abs(grid1%rad(i) - grid2%rad(i))
3086 dw = abs(grid1%wr(i) - grid2%wr(i))
3087 IF (dr + dw > 1.0e-12_dp)
THEN
Define the atom type and its sub types.
integer, parameter, public num_basis
subroutine, public read_atom_opt_section(optimization, opt_section)
...
subroutine, public create_atom_type(atom)
...
integer, parameter, public cgto_basis
integer, parameter, public gto_basis
integer, parameter, public sto_basis
subroutine, public release_atom_type(atom)
...
subroutine, public release_opmat(opmat)
...
subroutine, public release_atom_potential(potential)
...
subroutine, public init_atom_basis(basis, basis_section, zval, btyp)
Initialize the basis for the atomic code.
logical function, public atom_compare_grids(grid1, grid2)
...
subroutine, public release_opgrid(opgrid)
...
subroutine, public release_atom_orbs(orbs)
...
integer, parameter, public lmat
subroutine, public set_atom(atom, basis, state, integrals, orbitals, potential, zcore, pp_calc, do_zmp, doread, read_vxc, method_type, relativistic, coulomb_integral_type, exchange_integral_type, fmat)
...
subroutine, public init_atom_potential(potential, potential_section, zval)
...
subroutine, public release_atom_basis(basis)
...
subroutine, public create_atom_orbs(orbs, mbas, mo)
...
subroutine, public init_atom_basis_default_pp(basis)
...
subroutine, public create_opmat(opmat, n, lmax)
...
subroutine, public atom_basis_gridrep(basis, gbasis, r, rab)
...
subroutine, public read_ecp_potential(element_symbol, potential, pseudo_name, pseudo_file, potential_section)
...
subroutine, public setup_hf_section(hf_frac, do_hfx, atom, xc_section, extype)
...
subroutine, public clementi_geobas(zval, cval, aval, ngto, ival)
...
subroutine, public create_opgrid(opgrid, grid)
...
Routines that process Quantum Espresso UPF files.
subroutine, public atom_read_upf(pot, upf_filename, read_header)
...
pure subroutine, public atom_release_upf(upfpot)
...
Calculates Bessel functions.
elemental impure real(kind=dp) function, public bessel0(x, l)
...
collects all references to literature in CP2K as new algorithms / method are included from literature...
integer, save, public limpanuparb2011
Utility routines to read data from files. Kept as close as possible to the old parser because.
subroutine, public parser_read_line(parser, nline, at_end)
Read the next line from a logical unit "unit" (I/O node only). Skip (nline-1) lines and skip also all...
subroutine, public parser_get_next_line(parser, nline, at_end)
Read the next input line and broadcast the input information. Skip (nline-1) lines and skip also all ...
character(len=3) function, public parser_test_next_token(parser, string_length)
Test next input object.
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.
Defines the basic variable types.
integer, parameter, public dp
integer, parameter, public default_string_length
An array-based list which grows on demand. When the internal array is full, a new array of twice the ...
Definition of mathematical constants and functions.
real(kind=dp), parameter, public pi
real(kind=dp), dimension(-1:2 *maxfac+1), parameter, public dfac
real(kind=dp), parameter, public rootpi
real(kind=dp), dimension(0:maxfac), parameter, public fac
Periodic Table related data definitions.
type(atom), dimension(0:nelem), public ptable
subroutine, public get_ptable_info(symbol, number, amass, ielement, covalent_radius, metallic_radius, vdw_radius, found)
Pass information about the kind given the element symbol.
subroutine, public deallocate_grid_atom(grid_atom)
Deallocate a Gaussian-type orbital (GTO) basis set data set.
subroutine, public allocate_grid_atom(grid_atom)
Initialize components of the grid_atom_type structure.
subroutine, public create_grid_atom(grid_atom, nr, na, llmax, ll, quadrature)
...
Utilities for string manipulations.
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.
Provides all information about a basis set.
Provides all information about a pseudopotential.
Provides info about hartree-fock exchange (For now, we only support potentials that can be represente...
Information on optimization procedure.
Holds atomic orbitals and energies.
Provides all information on states and occupation.
Provides all information about an atomic kind.