57#include "./base/base_uses.f90"
63 CHARACTER(len=*),
PARAMETER,
PRIVATE :: moduleN =
'atom_energy'
83 CHARACTER(len=*),
PARAMETER :: routinen =
'atom_energy_opt'
85 CHARACTER(LEN=2) :: elem
86 CHARACTER(LEN=default_string_length) :: filename
87 CHARACTER(LEN=default_string_length), &
88 DIMENSION(:),
POINTER :: tmpstringlist
89 INTEGER :: do_eric, do_erie, handle, i, im, in, iw, k, maxl, mb, method, mo, n_meth, n_rep, &
90 nder, nr_gh, num_gau, num_gto, num_pol, reltyp, zcore, zval, zz
91 INTEGER,
DIMENSION(0:lmat) :: maxn
92 INTEGER,
DIMENSION(:),
POINTER :: cn
93 LOGICAL :: dm, do_gh, do_zmp, doread, eri_c, eri_e, &
94 had_ae, had_pp, lcomp, lcond, pp_calc, &
96 REAL(kind=
dp) :: crad, delta, lambda
97 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:) :: ext_density, ext_vxc
98 REAL(kind=
dp),
DIMENSION(0:lmat, 10) :: pocc
103 TYPE(
atom_p_type),
DIMENSION(:, :),
POINTER :: atom_info
107 TYPE(
section_vals_type),
POINTER :: admm_section, basis_section, external_vxc_section, &
108 method_section, opt_section, potential_section, powell_section, sgp_section, xc_section, &
109 zmp_restart_section, zmp_section
111 CALL timeset(routinen, handle)
118 IF (
ptable(i)%symbol == elem)
THEN
123 IF (zz /= 1) zval = zz
126 ALLOCATE (ae_basis, pp_basis)
128 NULLIFY (ae_basis%grid)
130 NULLIFY (pp_basis%grid)
137 IF (iw > 0)
CALL atom_print_info(zval,
"Atomic Energy Calculation", iw)
147 NULLIFY (potential_section)
149 ALLOCATE (ae_pot, p_pot)
167 DO in = 1, min(
SIZE(cn), 4)
168 maxn(in - 1) = cn(in)
171 maxn(in) = min(maxn(in), ae_basis%nbas(in))
172 maxn(in) = min(maxn(in), pp_basis%nbas(in))
189 ALLOCATE (ae_int, pp_int)
191 ALLOCATE (atom_info(n_rep, n_meth))
196 NULLIFY (atom_info(in, im)%atom)
199 atom_info(in, im)%atom%optimization = optimization
201 atom_info(in, im)%atom%z = zval
203 atom_info(in, im)%atom%xc_section => xc_section
212 atom_info(in, im)%atom%do_zmp = do_zmp
214 atom_info(in, im)%atom%ext_file = filename
216 r_val=atom_info(in, im)%atom%zmpgrid_tol)
218 atom_info(in, im)%atom%lambda = lambda
220 atom_info(in, im)%atom%dm = dm
224 atom_info(in, im)%atom%doread = doread
226 atom_info(in, im)%atom%zmp_restart_file = filename
231 atom_info(in, im)%atom%read_vxc = read_vxc
233 atom_info(in, im)%atom%ext_vxc_file = filename
235 r_val=atom_info(in, im)%atom%zmpvxcgrid_tol)
241 c_vals=tmpstringlist)
249 state%maxl_calc = max(maxl, state%maxl_occ)
250 state%maxl_calc = min(
lmat, state%maxl_calc)
252 DO k = 0, state%maxl_calc
253 state%maxn_calc(k) = max(maxn(k), state%maxn_occ(k))
257 pp_calc = any(index(tmpstringlist(1:),
"CORE") /= 0)
262 zcore = zval - nint(sum(state%core))
263 CALL set_atom(atom_info(in, im)%atom, zcore=zcore, pp_calc=.true.)
266 CALL set_atom(atom_info(in, im)%atom, zcore=zval, pp_calc=.false.)
271 CALL set_atom(atom_info(in, im)%atom, method_type=method, relativistic=reltyp)
291 potential=p_pot, eri_coulomb=eri_c, eri_exchange=eri_e)
295 NULLIFY (pp_int%tzora, pp_int%hdkh)
297 CALL set_atom(atom_info(in, im)%atom, basis=pp_basis, integrals=pp_int, potential=p_pot)
298 state%maxn_calc(:) = min(state%maxn_calc(:), pp_basis%nbas(:))
299 cpassert(all(state%maxn_calc(:) >= state%maxn_occ))
304 eri_coulomb=eri_c, eri_exchange=eri_e)
310 CALL set_atom(atom_info(in, im)%atom, basis=ae_basis, integrals=ae_int, potential=ae_pot)
311 state%maxn_calc(:) = min(state%maxn_calc(:), ae_basis%nbas(:))
312 cpassert(all(state%maxn_calc(:) >= state%maxn_occ))
316 CALL set_atom(atom_info(in, im)%atom, state=state)
318 CALL set_atom(atom_info(in, im)%atom, coulomb_integral_type=do_eric, &
319 exchange_integral_type=do_erie)
320 atom_info(in, im)%atom%hfx_pot%do_gh = do_gh
321 atom_info(in, im)%atom%hfx_pot%nr_gh = nr_gh
324 mo = maxval(state%maxn_calc)
325 mb = maxval(atom_info(in, im)%atom%basis%nbas)
327 CALL set_atom(atom_info(in, im)%atom, orbitals=orbitals)
336 IF (atom_info(in, im)%atom%do_zmp)
THEN
339 ALLOCATE (ext_density(atom_info(in, im)%atom%basis%grid%nr))
343 atom=atom_info(in, im)%atom, &
345 DEALLOCATE (ext_density)
348 IF (atom_info(in, im)%atom%read_vxc)
THEN
349 ALLOCATE (ext_vxc(atom_info(in, im)%atom%basis%grid%nr))
353 atom=atom_info(in, im)%atom, &
370 CALL atom_fit_density(atom_info(in, im)%atom, num_gto, 0, iw, powell_section=powell_section)
380 CALL atom_fit_kgpot(atom_info(in, im)%atom, num_gau, num_pol, iw, powell_section=powell_section)
389 CALL atom_response_basis(atom_info(in, im)%atom, delta, nder, iw)
395 file_position=
"REWIND")
397 CALL atom_write_upf(atom_info(in, im)%atom, iw)
406 iw =
cp_print_key_unit_nr(logger, atom_section,
"PRINT%GEOMETRICAL_RESPONSE_BASIS", extension=
".log")
415 CALL section_vals_val_get(atom_section,
"PRINT%ANALYZE_BASIS%OVERLAP_CONDITION_NUMBER", l_val=lcond)
433 CALL atom_admm(atom_info, admm_section, iw)
438 iw =
cp_print_key_unit_nr(logger, atom_section,
"PRINT%SEPARABLE_GAUSSIAN_PSEUDO", extension=
".log")
467 DEALLOCATE (atom_info)
469 DEALLOCATE (ae_pot, p_pot, ae_basis, pp_basis, ae_int, pp_int)
471 CALL timestop(handle)
485 SUBROUTINE atom_response_basis(atom, delta, nder, iw)
488 REAL(kind=
dp),
INTENT(IN) :: delta
489 INTEGER,
INTENT(IN) :: nder, iw
491 INTEGER :: i, ider, j, k, l, lhomo, m, n, nhomo, &
493 REAL(kind=
dp) :: dene, emax, expzet, fhomo, o, prefac, &
495 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:, :) :: amat
496 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:, :, :) :: rbasis
497 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:, :, :, :) :: wfn
498 REAL(kind=
dp),
DIMENSION(:, :, :),
POINTER :: ovlp
501 WRITE (iw,
'(/," ",79("*"),/,T30,A,A/," ",79("*"))')
"RESPONSE BASIS for ",
ptable(
atom%z)%symbol
504 ovlp =>
atom%integrals%ovlp
510 DO l = 0, state%maxl_occ
511 DO i = 1, state%maxn_occ(l)
512 IF (
atom%orbitals%ener(i, l) > emax)
THEN
515 emax =
atom%orbitals%ener(i, l)
516 fhomo = state%occupation(l, i)
521 s1 =
SIZE(
atom%orbitals%wfn, 1)
522 s2 =
SIZE(
atom%orbitals%wfn, 2)
523 ALLOCATE (wfn(s1, s2, 0:
lmat, -nder:nder))
524 s2 = maxval(state%maxn_occ) + nder
525 ALLOCATE (rbasis(s1, s2, 0:
lmat))
528 DO ider = -nder, nder
529 dene = real(ider, kind=
dp)*delta
530 cpassert(fhomo > abs(dene))
531 state%occupation(lhomo, nhomo) = fhomo + dene
533 wfn(:, :, :, ider) =
atom%orbitals%wfn
534 state%occupation(lhomo, nhomo) = fhomo
537 DO l = 0, state%maxl_occ
539 DO i = 1, max(state%maxn_occ(l), 1)
540 rbasis(:, i, l) = wfn(:, i, l, 0)
544 i = max(state%maxn_occ(l), 1)
547 rbasis(:, i + 1, l) = 0.5_dp*(wfn(:, i, l, 1) - wfn(:, i, l, -1))/delta
549 rbasis(:, i + 2, l) = 0.25_dp*(wfn(:, i, l, 2) - 2._dp*wfn(:, i, l, 0) + wfn(:, i, l, -2))/delta**2
551 rbasis(:, i + 3, l) = 0.125_dp*(wfn(:, i, l, 3) - 3._dp*wfn(:, i, l, 1) &
552 + 3._dp*wfn(:, i, l, -1) - wfn(:, i, l, -3))/delta**3
559 n = state%maxn_occ(l) + nder
560 m =
atom%basis%nbas(l)
563 o = dot_product(rbasis(1:m, j, l), reshape(matmul(ovlp(1:m, 1:m, l), rbasis(1:m, i:i, l)), (/m/)))
564 rbasis(1:m, i, l) = rbasis(1:m, i, l) - o*rbasis(1:m, j, l)
566 o = dot_product(rbasis(1:m, i, l), reshape(matmul(ovlp(1:m, 1:m, l), rbasis(1:m, i:i, l)), (/m/)))
567 rbasis(1:m, i, l) = rbasis(1:m, i, l)/sqrt(o)
571 ALLOCATE (amat(n, n))
572 amat(1:n, 1:n) = matmul(transpose(rbasis(1:m, 1:n, l)), matmul(ovlp(1:m, 1:m, l), rbasis(1:m, 1:n, l)))
574 amat(i, i) = amat(i, i) - 1._dp
576 IF (maxval(abs(amat)) > 1.e-12)
THEN
577 WRITE (iw,
'(/,A,G20.10)')
" Orthogonality error ", maxval(abs(amat))
582 WRITE (iw,
'(/,A,T30,I3)')
" Angular momentum :", l
584 WRITE (iw,
'(/,A,I0,A,I0,A)')
" Exponent Coef.(Quickstep Normalization), first ", &
585 n - nder,
" valence ", nder,
" response"
586 expzet = 0.25_dp*real(2*l + 3,
dp)
587 prefac = sqrt(
rootpi/2._dp**(l + 2)*
dfac(2*l + 1))
589 zeta = (2._dp*
atom%basis%am(i, l))**expzet
590 WRITE (iw,
'(4X,F20.10,4X,15ES20.6)')
atom%basis%am(i, l), ((prefac*rbasis(i, k, l)/zeta), k=1, n)
595 DEALLOCATE (wfn, rbasis)
597 WRITE (iw,
'(" ",79("*"))')
599 END SUBROUTINE atom_response_basis
608 SUBROUTINE atom_write_upf(atom, iw)
611 INTEGER,
INTENT(IN) :: iw
613 CHARACTER(LEN=default_string_length) :: string
614 INTEGER :: i, ibeta, j, k, l, lmax, nbeta, nr, &
617 REAL(kind=
dp) :: pf, rl, rmax
618 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:) :: beta, corden, dens, ef, locpot, rp
619 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:, :) :: dij
622 IF (.NOT.
atom%pp_calc)
RETURN
624 pot =>
atom%potential%gth_pot
625 cpassert(.NOT. pot%lsdpot)
627 WRITE (iw,
'(A)')
'<UPF version="2.0.1">'
629 WRITE (iw,
'(T4,A)')
'<PP_INFO>'
630 WRITE (iw,
'(T8,A)')
'Converted from CP2K GTH format'
631 WRITE (iw,
'(T8,A)')
'<PP_INPUTFILE>'
633 WRITE (iw,
'(T8,A)')
'</PP_INPUTFILE>'
634 WRITE (iw,
'(T4,A)')
'</PP_INFO>'
635 WRITE (iw,
'(T4,A)')
'<PP_HEADER'
636 WRITE (iw,
'(T8,A)')
'generated="Generated in analytical, separable form"'
637 WRITE (iw,
'(T8,A)')
'author="Goedecker/Hartwigsen/Hutter/Teter"'
638 WRITE (iw,
'(T8,A)')
'date="Phys.Rev.B58, 3641 (1998); B54, 1703 (1996)"'
639 WRITE (iw,
'(T8,A)')
'comment="This file generated by CP2K ATOM code"'
640 CALL compose(string,
"element", cval=pot%symbol)
641 WRITE (iw,
'(T8,A)') trim(string)
642 WRITE (iw,
'(T8,A)')
'pseudo_type="NC"'
643 WRITE (iw,
'(T8,A)')
'relativistic="no"'
644 WRITE (iw,
'(T8,A)')
'is_ultrasoft="F"'
645 WRITE (iw,
'(T8,A)')
'is_paw="F"'
646 WRITE (iw,
'(T8,A)')
'is_coulomb="F"'
647 WRITE (iw,
'(T8,A)')
'has_so="F"'
648 WRITE (iw,
'(T8,A)')
'has_wfc="F"'
649 WRITE (iw,
'(T8,A)')
'has_gipaw="F"'
650 WRITE (iw,
'(T8,A)')
'paw_as_gipaw="F"'
652 WRITE (iw,
'(T8,A)')
'core_correction="T"'
654 WRITE (iw,
'(T8,A)')
'core_correction="F"'
656 WRITE (iw,
'(T8,A)')
'functional="DFT"'
657 CALL compose(string,
"z_valence", rval=pot%zion)
658 WRITE (iw,
'(T8,A)') trim(string)
659 CALL compose(string,
"total_psenergy", rval=2._dp*
atom%energy%etot)
660 WRITE (iw,
'(T8,A)') trim(string)
661 WRITE (iw,
'(T8,A)')
'wfc_cutoff="0.0E+00"'
662 WRITE (iw,
'(T8,A)')
'rho_cutoff="0.0E+00"'
665 IF (pot%nl(l) > 0) lmax = l
667 CALL compose(string,
"l_max", ival=lmax)
668 WRITE (iw,
'(T8,A)') trim(string)
669 WRITE (iw,
'(T8,A)')
'l_max_rho="0"'
670 WRITE (iw,
'(T8,A)')
'l_local="-3"'
671 nr =
atom%basis%grid%nr
672 CALL compose(string,
"mesh_size", ival=nr)
673 WRITE (iw,
'(T8,A)') trim(string)
674 nwfc = sum(
atom%state%maxn_occ)
675 CALL compose(string,
"number_of_wfc", ival=nwfc)
676 WRITE (iw,
'(T8,A)') trim(string)
678 CALL compose(string,
"number_of_proj", ival=nbeta)
679 WRITE (iw,
'(T8,A)') trim(string)//
'/>'
682 up =
atom%basis%grid%rad(1) <
atom%basis%grid%rad(nr)
683 WRITE (iw,
'(T4,A)')
'<PP_MESH'
684 WRITE (iw,
'(T8,A)')
'dx="1.E+00"'
685 CALL compose(string,
"mesh", ival=nr)
686 WRITE (iw,
'(T8,A)') trim(string)
687 WRITE (iw,
'(T8,A)')
'xmin="1.E+00"'
688 rmax = maxval(
atom%basis%grid%rad)
689 CALL compose(string,
"rmax", rval=rmax)
690 WRITE (iw,
'(T8,A)') trim(string)
691 WRITE (iw,
'(T8,A)')
'zmesh="1.E+00">'
692 WRITE (iw,
'(T8,A)')
'<PP_R type="real"'
693 CALL compose(string,
"size", ival=nr)
694 WRITE (iw,
'(T12,A)') trim(string)
695 WRITE (iw,
'(T12,A)')
'columns="4">'
697 WRITE (iw,
'(T12,4ES25.12E3)') (
atom%basis%grid%rad(i), i=1, nr)
699 WRITE (iw,
'(T12,4ES25.12E3)') (
atom%basis%grid%rad(i), i=nr, 1, -1)
701 WRITE (iw,
'(T8,A)')
'</PP_R>'
702 WRITE (iw,
'(T8,A)')
'<PP_RAB type="real"'
703 CALL compose(string,
"size", ival=nr)
704 WRITE (iw,
'(T12,A)') trim(string)
705 WRITE (iw,
'(T12,A)')
'columns="4">'
707 WRITE (iw,
'(T12,4ES25.12E3)') (
atom%basis%grid%wr(i)/
atom%basis%grid%rad2(i), i=1, nr)
709 WRITE (iw,
'(T12,4ES25.12E3)') (
atom%basis%grid%wr(i)/
atom%basis%grid%rad2(i), i=nr, 1, -1)
711 WRITE (iw,
'(T8,A)')
'</PP_RAB>'
712 WRITE (iw,
'(T4,A)')
'</PP_MESH>'
716 WRITE (iw,
'(T4,A)')
'<PP_NLCC type="real"'
717 CALL compose(string,
"size", ival=nr)
718 WRITE (iw,
'(T8,A)') trim(string)
719 WRITE (iw,
'(T8,A)')
'columns="4">'
720 ALLOCATE (corden(nr))
724 WRITE (iw,
'(T8,4ES25.12E3)') (corden(i), i=1, nr)
726 WRITE (iw,
'(T8,4ES25.12E3)') (corden(i), i=nr, 1, -1)
729 WRITE (iw,
'(T4,A)')
'</PP_NLCC>'
733 WRITE (iw,
'(T4,A)')
'<PP_LOCAL type="real"'
734 CALL compose(string,
"size", ival=nr)
735 WRITE (iw,
'(T8,A)') trim(string)
736 WRITE (iw,
'(T8,A)')
'columns="4">'
737 ALLOCATE (locpot(nr))
741 WRITE (iw,
'(T8,4ES25.12E3)') (2.0_dp*locpot(i), i=1, nr)
743 WRITE (iw,
'(T8,4ES25.12E3)') (2.0_dp*locpot(i), i=nr, 1, -1)
746 WRITE (iw,
'(T4,A)')
'</PP_LOCAL>'
749 WRITE (iw,
'(T4,A)')
'<PP_NONLOCAL>'
750 ALLOCATE (rp(nr), ef(nr), beta(nr))
753 IF (pot%nl(l) == 0) cycle
755 rp(:) =
atom%basis%grid%rad
756 ef(:) = exp(-0.5_dp*rp*rp/(rl*rl))
758 pf = rl**(l + 0.5_dp*(4._dp*i - 1._dp))
760 pf = sqrt(2._dp)/(pf*sqrt(
gamma1(j)))
761 beta(:) = pf*rp**(l + 2*i - 2)*ef
763 CALL compose(string,
"<PP_BETA", counter=ibeta)
764 WRITE (iw,
'(T8,A)') trim(string)
765 CALL compose(string,
"angular_momentum", ival=l)
766 WRITE (iw,
'(T12,A)') trim(string)
767 WRITE (iw,
'(T12,A)')
'type="real"'
768 CALL compose(string,
"size", ival=nr)
769 WRITE (iw,
'(T12,A)') trim(string)
770 WRITE (iw,
'(T12,A)')
'columns="4">'
773 WRITE (iw,
'(T12,4ES25.12E3)') (2._dp*beta(j), j=1, nr)
775 WRITE (iw,
'(T12,4ES25.12E3)') (2._dp*beta(j), j=nr, 1, -1)
777 CALL compose(string,
"</PP_BETA", counter=ibeta, isfinal=.true.)
778 WRITE (iw,
'(T8,A)') trim(string)
781 DEALLOCATE (rp, ef, beta)
783 ALLOCATE (dij(nbeta, nbeta))
786 IF (pot%nl(l) == 0) cycle
787 ibeta = sum(pot%nl(0:l - 1)) + 1
788 i = ibeta + pot%nl(l) - 1
789 dij(ibeta:i, ibeta:i) = pot%hnl(1:pot%nl(l), 1:pot%nl(l), l)
791 WRITE (iw,
'(T8,A)')
'<PP_DIJ type="real"'
792 WRITE (iw,
'(T12,A)')
'columns="4">'
793 WRITE (iw,
'(T12,4ES25.12E3)') ((0.5_dp*dij(i, j), j=1, nbeta), i=1, nbeta)
794 WRITE (iw,
'(T8,A)')
'</PP_DIJ>'
796 WRITE (iw,
'(T4,A)')
'</PP_NONLOCAL>'
800 WRITE (iw,
'(T4,A)')
'<PP_PSWFC>'
804 IF (abs(
atom%state%occupation(l, i)) == 0._dp) cycle
806 CALL compose(string,
"<PP_CHI", counter=nwfn)
807 WRITE (iw,
'(T8,A)') trim(string)
808 CALL compose(string,
"l", ival=l)
809 WRITE (iw,
'(T12,A)') trim(string)
810 CALL compose(string,
"occupation", rval=
atom%state%occupation(l, i))
811 WRITE (iw,
'(T12,A)') trim(string)
812 CALL compose(string,
"pseudo_energy", rval=2._dp*
atom%orbitals%ener(i, l))
813 WRITE (iw,
'(T12,A)') trim(string)
814 WRITE (iw,
'(T12,A)')
'columns="4">'
816 DO k = 1,
atom%basis%nbas(l)
817 beta(:) = beta(:) +
atom%orbitals%wfn(k, i, l)*
atom%basis%bf(:, k, l)
819 beta(:) = beta*
atom%basis%grid%rad
821 WRITE (iw,
'(T12,4ES25.12E3)') (beta(j)*
atom%basis%grid%rad(j), j=1, nr)
823 WRITE (iw,
'(T12,4ES25.12E3)') (beta(j)*
atom%basis%grid%rad(j), j=nr, 1, -1)
825 CALL compose(string,
"</PP_CHI", counter=nwfn, isfinal=.true.)
826 WRITE (iw,
'(T8,A)') trim(string)
829 WRITE (iw,
'(T4,A)')
'</PP_PSWFC>'
834 WRITE (iw,
'(T4,A)')
'<PP_RHOATOM type="real"'
835 CALL compose(string,
"size", ival=nr)
836 WRITE (iw,
'(T8,A)') trim(string)
837 WRITE (iw,
'(T8,A)')
'columns="4">'
839 "RHO",
atom%basis%grid%rad)
841 WRITE (iw,
'(T8,4ES25.12E3)') (4._dp*
pi*dens(j)*
atom%basis%grid%rad2(j), j=1, nr)
843 WRITE (iw,
'(T8,4ES25.12E3)') (4._dp*
pi*dens(j)*
atom%basis%grid%rad2(j), j=nr, 1, -1)
845 WRITE (iw,
'(T4,A)')
'</PP_RHOATOM>'
848 WRITE (iw,
'(A)')
'</UPF>'
850 END SUBROUTINE atom_write_upf
864 SUBROUTINE compose(string, tag, counter, rval, ival, cval, isfinal)
865 CHARACTER(len=*),
INTENT(OUT) :: string
866 CHARACTER(len=*),
INTENT(IN) :: tag
867 INTEGER,
INTENT(IN),
OPTIONAL :: counter
868 REAL(kind=
dp),
INTENT(IN),
OPTIONAL :: rval
869 INTEGER,
INTENT(IN),
OPTIONAL :: ival
870 CHARACTER(len=*),
INTENT(IN),
OPTIONAL :: cval
871 LOGICAL,
INTENT(IN),
OPTIONAL :: isfinal
873 CHARACTER(LEN=default_string_length) :: str
876 IF (
PRESENT(counter))
THEN
877 WRITE (str,
"(I12)") counter
878 ELSEIF (
PRESENT(rval))
THEN
879 WRITE (str,
"(G18.8)") rval
880 ELSEIF (
PRESENT(ival))
THEN
881 WRITE (str,
"(I12)") ival
882 ELSEIF (
PRESENT(cval))
THEN
883 WRITE (str,
"(A)") trim(adjustl(cval))
885 WRITE (str,
"(A)")
""
888 IF (
PRESENT(isfinal)) fin = isfinal
889 IF (
PRESENT(counter))
THEN
891 WRITE (string,
"(A,A1,A,A1)") trim(adjustl(tag)),
'.', trim(adjustl(str)),
'>'
893 WRITE (string,
"(A,A1,A)") trim(adjustl(tag)),
'.', trim(adjustl(str))
897 WRITE (string,
"(A,A2,A,A2)") trim(adjustl(tag)),
'="', trim(adjustl(str)),
'>"'
899 WRITE (string,
"(A,A2,A,A1)") trim(adjustl(tag)),
'="', trim(adjustl(str)),
'"'
902 END SUBROUTINE compose
static GRID_HOST_DEVICE orbital up(const int i, const orbital a)
Increase i'th component of given orbital angular momentum.
subroutine, public atom_admm(atom_info, admm_section, iw)
Analysis of ADMM approximation to exact exchange.
subroutine, public calculate_atom(atom, iw, noguess, converged)
General routine to perform electronic structure atomic calculations.
subroutine, public atom_energy_opt(atom_section)
Compute the atomic energy.
routines that fit parameters for /from atomic calculations
subroutine, public atom_fit_density(atom, num_gto, norder, iunit, powell_section, results)
Fit the atomic electron density using a geometrical Gaussian basis set.
subroutine, public atom_fit_kgpot(atom, num_gau, num_pol, iunit, powell_section, results)
...
subroutine, public atom_grb_construction(atom_info, atom_section, iw)
Construct geometrical response basis set.
Calculate the atomic operator matrices.
subroutine, public atom_ppint_release(integrals)
Release memory allocated for atomic integrals (core electrons).
subroutine, public atom_int_setup(integrals, basis, potential, eri_coulomb, eri_exchange, all_nu)
Set up atomic integrals.
subroutine, public atom_relint_setup(integrals, basis, reltyp, zcore, alpha)
...
subroutine, public atom_relint_release(integrals)
Release memory allocated for atomic integrals (relativistic effects).
subroutine, public atom_ppint_setup(integrals, basis, potential)
...
subroutine, public atom_int_release(integrals)
Release memory allocated for atomic integrals (valence electrons).
Routines that print various information about an atomic kind.
subroutine, public atom_print_basis(atom_basis, iw, title)
Print atomic basis set.
subroutine, public atom_write_pseudo_param(gthpot, iunit)
Print GTH pseudo-potential parameters.
subroutine, public atom_print_method(atom, iw)
Print information about the electronic structure method in use.
subroutine, public atom_print_orbitals(atom, iw)
Print atomic orbitals.
subroutine, public atom_print_potential(potential, iw)
Print information about the pseudo-potential.
subroutine, public atom_print_info(zval, info, iw)
Print an information string related to the atomic kind.
subroutine, public atom_sgp_construction(atom_info, input_section, iw)
...
Define the atom type and its sub types.
subroutine, public read_atom_opt_section(optimization, opt_section)
...
subroutine, public create_atom_type(atom)
...
subroutine, public release_atom_type(atom)
...
subroutine, public release_atom_potential(potential)
...
subroutine, public init_atom_basis(basis, basis_section, zval, btyp)
Initialize the basis for the atomic code.
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)
...
Some basic routines for atomic calculations.
subroutine, public atom_condnumber(basis, crad, iw)
Print condition numbers of the given atomic basis set.
pure subroutine, public atom_local_potential(locpot, gthpot, rr)
...
subroutine, public atom_read_external_vxc(vxc, atom, iw)
ZMP subroutine to read external v_xc in the atomic code.
pure logical function, public atom_consistent_method(method, multiplicity)
Check that the atomic multiplicity is consistent with the electronic structure method.
subroutine, public atom_core_density(corden, potential, typ, rr)
...
pure integer function, dimension(0:lmat), public get_maxn_occ(occupation)
Return the maximum principal quantum number of occupied orbitals.
subroutine, public atom_density(density, pmat, basis, maxl, typ, rr)
Map the electron density on an atomic radial grid.
subroutine, public atom_read_external_density(density, atom, iw)
ZMP subroutine to read external density from linear grid of density matrix.
subroutine, public atom_completeness(basis, zv, iw)
Estimate completeness of the given atomic basis set.
subroutine, public atom_write_zmp_restart(atom)
ZMP subroutine to write external restart file.
subroutine, public atom_set_occupation(ostring, occupation, wfnocc, multiplicity)
Set occupation of atomic orbitals.
pure integer function, public get_maxl_occ(occupation)
Return the maximum orbital quantum number of occupied orbitals.
routines that build the integrals of the Vxc potential calculated for the atomic code
subroutine, public calculate_atom_ext_vxc(vxc, atom, lprint, xcmat)
ZMP subroutine for the scf external density from the external v_xc.
subroutine, public calculate_atom_zmp(ext_density, atom, lprint, xcmat)
ZMP subroutine for the calculation of the effective constraint potential in the atomic code.
various routines to log and control the output. The idea is that decisions about where to log should ...
type(cp_logger_type) function, pointer, public cp_get_default_logger()
returns the default logger
routines to handle the output, The idea is to remove the decision of wheter to output and what to out...
integer function, public cp_print_key_unit_nr(logger, basis_section, print_key_path, extension, middle_name, local, log_filename, ignore_should_output, file_form, file_position, file_action, file_status, do_backup, on_file, is_new_file, mpi_io, fout)
...
subroutine, public cp_print_key_finished_output(unit_nr, logger, basis_section, print_key_path, local, ignore_should_output, on_file, mpi_io)
should be called after you finish working with a unit obtained with cp_print_key_unit_nr,...
Defines the basic variable types.
integer, parameter, public dp
integer, parameter, public default_string_length
Definition of mathematical constants and functions.
real(kind=dp), dimension(0:maxfac), parameter, public gamma1
real(kind=dp), parameter, public pi
real(kind=dp), dimension(-1:2 *maxfac+1), parameter, public dfac
real(kind=dp), parameter, public rootpi
Periodic Table related data definitions.
type(atom), dimension(0:nelem), public ptable
integer, parameter, public nelem
Definition of physical constants:
real(kind=dp), parameter, public bohr
Provides all information about a basis set.
Provides all information about a pseudopotential.
Information on optimization procedure.
Holds atomic orbitals and energies.
Provides all information on states and occupation.
Provides all information about an atomic kind.
type of a logger, at the moment it contains just a print level starting at which level it should be l...