45#include "./base/base_uses.f90"
51 CHARACTER(len=*),
PARAMETER,
PRIVATE :: moduleN =
'atom_output'
69 INTEGER,
INTENT(IN) :: zval
70 CHARACTER(len=*),
INTENT(IN) :: info
71 INTEGER,
INTENT(IN) :: iw
73 WRITE (iw,
'(/," ",A,T40,A," [",A,"]",T62,"Atomic number:",T78,I3,/)') &
74 adjustl(trim(info)), trim(
ptable(zval)%name), trim(
ptable(zval)%symbol), zval
89 INTEGER,
INTENT(IN) :: iw
91 CHARACTER(LEN=1),
DIMENSION(0:7),
PARAMETER :: &
92 label = [
"S",
"P",
"D",
"F",
"G",
"H",
"I",
"K"]
94 INTEGER :: j, l, mc, mlc, mlo, mm(0:
lmat), mo
97 WRITE (iw,
'(/,T2,A)')
"Electronic structure"
98 WRITE (iw,
'(T5,A,T71,F10.2)')
"Total number of core electrons", sum(state%core)
99 WRITE (iw,
'(T5,A,T71,F10.2)')
"Total number of valence electrons", sum(state%occ)
100 WRITE (iw,
'(T5,A,T71,F10.2)')
"Total number of electrons", sum(state%occ + state%core)
101 SELECT CASE (state%multiplicity)
103 WRITE (iw,
'(T5,A,T68,A)')
"Multiplicity",
"not specified"
105 WRITE (iw,
'(T5,A,T72,A)')
"Multiplicity",
"high spin"
107 WRITE (iw,
'(T5,A,T73,A)')
"Multiplicity",
"low spin"
109 WRITE (iw,
'(T5,A,T74,A)')
"Multiplicity",
"singlet"
111 WRITE (iw,
'(T5,A,T74,A)')
"Multiplicity",
"doublet"
113 WRITE (iw,
'(T5,A,T74,A)')
"Multiplicity",
"triplet"
115 WRITE (iw,
'(T5,A,T74,A)')
"Multiplicity",
"quartet"
117 WRITE (iw,
'(T5,A,T74,A)')
"Multiplicity",
"quintet"
119 WRITE (iw,
'(T5,A,T75,A)')
"Multiplicity",
"sextet"
121 WRITE (iw,
'(T5,A,T75,A)')
"Multiplicity",
"septet"
129 IF (state%multiplicity == -1)
THEN
130 DO l = 0, max(mlo, mlc)
131 mo = state%maxn_occ(l)
132 IF (sum(state%core(l, :)) == 0)
THEN
133 WRITE (iw,
'(A5,T10,10F6.2)') label(l), (state%occ(l, j), j=1, mo)
136 cpassert(sum(state%occ(l, 1:mc)) == 0)
137 WRITE (iw, advance=
"no", fmt=
'(A5,T9,A1,10F6.2)') label(l),
"[", (state%core(l, j), j=1, mc)
138 WRITE (iw, fmt=
'(A1,F5.2,10F6.2)')
"]", (state%occ(l, j), j=mc + 1, mc + mo)
142 WRITE (iw,
'(T5,A)')
"Alpha Electrons"
143 DO l = 0, max(mlo, mlc)
144 mo = state%maxn_occ(l)
145 IF (sum(state%core(l, :)) == 0)
THEN
146 WRITE (iw,
'(A5,T10,10F6.2)') label(l), (state%occa(l, j), j=1, mo)
149 WRITE (iw, advance=
"no", fmt=
'(A5,T9,A1,10F6.2)') label(l),
"[", (0.5_dp*state%core(l, j), j=1, mc)
150 WRITE (iw, fmt=
'(A1,F5.2,10F6.2)')
"]", (state%occa(l, j), j=1, mo)
153 WRITE (iw,
'(T5,A)')
"Beta Electrons"
154 DO l = 0, max(mlo, mlc)
155 mo = state%maxn_occ(l)
156 IF (sum(state%core(l, :)) == 0)
THEN
157 WRITE (iw,
'(A5,T10,10F6.2)') label(l), (state%occb(l, j), j=1, mo)
160 WRITE (iw, advance=
"no", fmt=
'(A5,T9,A1,10F6.2)') label(l),
"[", (0.5_dp*state%core(l, j), j=1, mc)
161 WRITE (iw, fmt=
'(A1,F5.2,10F6.2)')
"]", (state%occb(l, j), j=1, mo)
181 INTEGER,
INTENT(IN) :: iw
184 REAL(kind=
dp) :: drho
186 WRITE (iw,
'(/,A,T36,A,T61,F20.12)')
" Energy components [Hartree]", &
187 " Total Energy ::",
atom%energy%etot
188 WRITE (iw,
'(T36,A,T61,F20.12)')
" Band Energy ::",
atom%energy%eband
189 WRITE (iw,
'(T36,A,T61,F20.12)')
" Kinetic Energy ::",
atom%energy%ekin
190 WRITE (iw,
'(T36,A,T61,F20.12)')
"Potential Energy ::",
atom%energy%epot
191 IF (
atom%energy%ekin /= 0.0_dp)
THEN
192 WRITE (iw,
'(T36,A,T61,F20.12)')
" Virial (-V/T) ::", -
atom%energy%epot/
atom%energy%ekin
194 WRITE (iw,
'(T36,A,T61,F20.12)')
" Core Energy ::",
atom%energy%ecore
195 IF (
atom%energy%exc /= 0._dp) &
196 WRITE (iw,
'(T36,A,T61,F20.12)')
" XC Energy ::",
atom%energy%exc
197 WRITE (iw,
'(T36,A,T61,F20.12)')
" Coulomb Energy ::",
atom%energy%ecoulomb
198 IF (
atom%energy%eexchange /= 0._dp) &
199 WRITE (iw,
'(T34,A,T61,F20.12)')
"HF Exchange Energy ::",
atom%energy%eexchange
201 WRITE (iw,
'(T20,A,T61,F20.12)')
" Total Pseudopotential Energy ::",
atom%energy%epseudo
202 WRITE (iw,
'(T20,A,T61,F20.12)')
" Local Pseudopotential Energy ::",
atom%energy%eploc
203 IF (
atom%energy%elsd /= 0._dp) &
204 WRITE (iw,
'(T20,A,T61,F20.12)')
" Local Spin-potential Energy ::",
atom%energy%elsd
205 WRITE (iw,
'(T20,A,T61,F20.12)')
" Nonlocal Pseudopotential Energy ::",
atom%energy%epnl
207 IF (
atom%potential%confinement)
THEN
208 WRITE (iw,
'(T36,A,T61,F20.12)')
" Confinement ::",
atom%energy%econfinement
211 IF (
atom%state%multiplicity == -1)
THEN
212 WRITE (iw,
'(/,A,T20,A,T30,A,T36,A,T49,A,T71,A,/)')
" Orbital energies", &
213 "State",
"L",
"Occupation",
"Energy[a.u.]",
"Energy[eV]"
214 DO l = 0,
atom%state%maxl_calc
215 n =
atom%state%maxn_calc(l)
217 WRITE (iw,
'(T23,I2,T30,I1,T36,F10.3,T46,F15.6,T66,F15.6)') &
218 i, l,
atom%state%occupation(l, i),
atom%orbitals%ener(i, l),
atom%orbitals%ener(i, l)*
evolt
220 IF (n > 0)
WRITE (iw, *)
223 WRITE (iw,
'(/,A,T20,A,T30,A,T36,A,T42,A,T55,A,T71,A,/)')
" Orbital energies", &
224 "State",
"Spin",
"L",
"Occupation",
"Energy[a.u.]",
"Energy[eV]"
225 DO l = 0,
atom%state%maxl_calc
226 n =
atom%state%maxn_calc(l)
228 WRITE (iw,
'(T23,I2,T29,A,T36,I1,T42,F10.3,T52,F15.6,T68,F13.6)') &
229 i,
"alpha", l,
atom%state%occa(l, i),
atom%orbitals%enera(i, l),
atom%orbitals%enera(i, l)*
evolt
232 WRITE (iw,
'(T23,I2,T29,A,T36,I1,T42,F10.3,T52,F15.6,T68,F13.6)') &
233 i,
" beta", l,
atom%state%occb(l, i),
atom%orbitals%enerb(i, l),
atom%orbitals%enerb(i, l)*
evolt
235 IF (n > 0)
WRITE (iw, *)
240 WRITE (iw,
'(/,A,T66,F15.6)')
" Total Electron Density at R=0: ", drho
253 INTEGER,
INTENT(IN) :: iter
254 REAL(
dp),
INTENT(IN) :: deps
256 INTEGER,
INTENT(IN) :: iw
259 WRITE (iw,
'(/," ",79("*"),/,T33,"Integral",T48,"Integral",/,T3,A,T16,A,T33,A,T46,A,T69,A/," ",79("*"))') &
260 "Iteration",
"Convergence",
"rho diff.",
"rho*v_xc[au]",
"Energy[au]"
262 WRITE (iw,
'(T3,I9,T15,G13.6,T30,G13.6,T46,G13.6,T61,F20.12)') iter, deps,
atom%rho_diff_integral, &
277 INTEGER,
INTENT(IN) :: iter
278 REAL(
dp),
INTENT(IN) :: deps, etot
279 INTEGER,
INTENT(IN) :: iw
282 WRITE (iw,
'(/," ",79("*"),/,T19,A,T38,A,T70,A,/," ",79("*"))') &
283 "Iteration",
"Convergence",
"Energy [au]"
285 WRITE (iw,
'(T20,i8,T34,G14.6,T61,F20.12)') iter, deps, etot
299 INTEGER,
INTENT(IN) :: iw
300 CHARACTER(len=*) :: title
304 WRITE (iw,
'(/,A)') trim(title)
308 WRITE (iw,
'(/," ",21("*"),A,22("*"))')
" Geometrical Gaussian Type Orbitals "
309 WRITE (iw,
'(A,F15.8,T41,A,F15.8)')
" Initial exponent: ",
atom_basis%aval, &
312 WRITE (iw,
'(/," ",21("*"),A,21("*"))')
" Uncontracted Gaussian Type Orbitals "
318 WRITE (iw,
'(/,T2,A,(T30,I5,T51,F30.8))') &
321 WRITE (iw,
'(/,T2,A,(T30,I5,T51,F30.8))') &
324 WRITE (iw,
'(/,T2,A,(T30,I5,T51,F30.8))') &
327 WRITE (iw,
'(/,T2,A,(T30,I5,T51,F30.8))') &
330 WRITE (iw,
'(/,T2,A,(T30,I5,T51,F30.8))') &
335 WRITE (iw,
'(" ",79("*"))')
337 WRITE (iw,
'(/," ",22("*"),A,22("*"))')
" Contracted Gaussian Type Orbitals "
340 IF (l == 0)
WRITE (iw,
'(A)')
" s Functions"
341 IF (l == 1)
WRITE (iw,
'(A)')
" p Functions"
342 IF (l == 2)
WRITE (iw,
'(A)')
" d Functions"
343 IF (l == 3)
WRITE (iw,
'(A)')
" f Functions"
344 IF (l >= 3)
WRITE (iw,
'(A)')
" x Functions"
346 WRITE (iw,
'(F15.6,5(T21,6F10.6,/))') &
351 WRITE (iw,
'(" ",79("*"))')
353 WRITE (iw,
'(/," ",28("*"),A,29("*"))')
" Slater Type Orbitals "
370 WRITE (iw,
'(" ",79("*"))')
392 REAL(kind=
dp),
DIMENSION(:, :, 0:),
OPTIONAL :: wfn
394 INTEGER :: i, im, iw, l
395 REAL(kind=
dp) :: expzet, prefac, zeta
397 CALL open_file(file_name=
"OPT_BASIS", file_status=
"UNKNOWN", file_action=
"WRITE", unit_number=iw)
401 WRITE (iw,
'(/," ",21("*"),A,22("*"))')
" Geometrical Gaussian Type Orbitals "
402 WRITE (iw,
'(A,F15.8,T41,A,F15.8)')
" Initial exponent: ",
atom_basis%aval, &
405 WRITE (iw,
'(T3,A)')
"BASIS_TYPE GAUSSIAN"
411 WRITE (iw,
'(T3,A,(T15,F20.8,:," \"))') &
414 WRITE (iw,
'(T3,A,(T15,F20.8,:," \"))') &
417 WRITE (iw,
'(T3,A,(T15,F20.8,:," \"))') &
420 WRITE (iw,
'(T3,A,(T15,F20.8,:," \"))') &
423 WRITE (iw,
'(T3,A,(T15,F20.8,:," \"))') &
431 WRITE (iw,
'(T3,A)')
"BASIS_TYPE SLATER"
436 WRITE (iw,
'(T3,A,(T15,F20.8,:," \"))') &
438 WRITE (iw,
'(T3,A,60I3)') &
441 WRITE (iw,
'(T3,A,(T15,F20.8,:," \"))') &
443 WRITE (iw,
'(T3,A,60I3)') &
446 WRITE (iw,
'(T3,A,(T15,F20.8,:," \"))') &
448 WRITE (iw,
'(T3,A,60I3)') &
451 WRITE (iw,
'(T3,A,(T15,F20.8,:," \"))') &
453 WRITE (iw,
'(T3,A,60I3)') &
456 WRITE (iw,
'(T3,A,(T15,F20.8,:," \"))') &
458 WRITE (iw,
'(T3,A,60I3)') &
469 IF (
PRESENT(wfn))
THEN
474 WRITE (iw,
'(/,T3,A)')
"ORBITAL COEFFICENTS (Quickstep normalization)"
475 im = min(6,
SIZE(wfn, 2))
478 WRITE (iw,
'(T3,A,I3)')
"L Quantum Number:", l
480 expzet = 0.25_dp*real(2*l + 3,
dp)
481 prefac = sqrt(
rootpi/2._dp**(l + 2)*
dfac(2*l + 1))
484 WRITE (iw,
'(T5,F14.8,2x,6F12.8)')
atom_basis%am(i, l), wfn(i, 1:im, l)*prefac/zeta
510 INTEGER,
INTENT(IN) :: iw
512 CHARACTER(len=160) :: shortform
513 CHARACTER(LEN=20) :: tmpstr
514 CHARACTER(len=:),
ALLOCATABLE :: reference
515 INTEGER :: ifun, il, meth, myfun, reltyp
519 NULLIFY (xc_fun, xc_fun_section, xc_section)
521 meth =
atom%method_type
523 xc_section =>
atom%xc_section
544 IF (iw > 0)
WRITE (iw, fmt=
"(/,' METHOD | Restricted Kohn-Sham Calculation')")
546 IF (iw > 0)
WRITE (iw, fmt=
"(/,' METHOD | Unrestricted Kohn-Sham Calculation')")
548 IF (iw > 0)
WRITE (iw, fmt=
"(/,' METHOD | Restricted Hartree-Fock Calculation')")
550 IF (iw > 0)
WRITE (iw, fmt=
"(/,' METHOD | Unrestricted Hartree-Fock Calculation')")
552 IF (iw > 0)
WRITE (iw, fmt=
"(/,' METHOD | Restricted Open-Shell Kohn-Sham Calculation')")
556 IF (
atom%do_zmp)
THEN
557 IF (iw > 0)
WRITE (iw, fmt=
"(' ZMP | Method on atomic radial density')")
558 IF (iw > 0)
WRITE (iw, fmt=
"(' ZMP | Lambda : ',F5.1)")
atom%lambda
559 IF (iw > 0)
WRITE (iw, fmt=
"(' ZMP | Reading external density : ',A20)")
atom%ext_file
561 IF (iw > 0)
WRITE (iw, fmt=
"(' ZMP | The file is in the form of a density matrix')")
563 IF (iw > 0)
WRITE (iw, fmt=
"(' ZMP | The file is in the form of a linear density')")
565 IF (
atom%doread)
THEN
566 IF (iw > 0)
WRITE (iw, fmt=
"(' ZMP | Restarting calculation from ',A20,' file if present')")
atom%zmp_restart_file
568 ELSE IF (
atom%read_vxc)
THEN
569 IF (iw > 0)
WRITE (iw, fmt=
"(' ZMP | Calculating density from external V_xc')")
570 IF (iw > 0)
WRITE (iw, fmt=
"(' ZMP | Reading external v_xc file : ',A20)")
atom%ext_vxc_file
573 IF (
atom%pp_calc)
THEN
574 IF (iw > 0)
WRITE (iw, fmt=
"(' METHOD | Nonrelativistic Calculation')")
576 reltyp =
atom%relativistic
582 IF (iw > 0)
WRITE (iw, fmt=
"(' METHOD | Nonrelativistic Calculation')")
584 IF (iw > 0)
WRITE (iw, fmt=
"(' METHOD | Relativistic Calculation using ZORA(MP)')")
586 IF (iw > 0)
WRITE (iw, fmt=
"(' METHOD | Relativistic Calculation using scaled ZORA(MP)')")
588 IF (iw > 0)
WRITE (iw, fmt=
"(' METHOD | Relativistic Calculation using Douglas-Kroll 0th order')")
589 IF (iw > 0)
WRITE (iw, fmt=
"(' METHOD | Relativistic Calculation using kietic energy scaling')")
591 IF (iw > 0)
WRITE (iw, fmt=
"(' METHOD | Relativistic Calculation using Douglas-Kroll 1st order')")
592 IF (iw > 0)
WRITE (iw, fmt=
"(' METHOD | Relativistic Calculation using Foldy-Wouthuysen transformation')")
594 IF (iw > 0)
WRITE (iw, fmt=
"(' METHOD | Relativistic Calculation using Douglas-Kroll 2nd order')")
596 IF (iw > 0)
WRITE (iw, fmt=
"(' METHOD | Relativistic Calculation using Douglas-Kroll 3rd order')")
604 IF (iw > 0)
WRITE (iw, fmt=
"(' FUNCTIONAL| ROUTINE=',a)") trim(tmpstr)
611 IF (.NOT.
ASSOCIATED(xc_fun))
EXIT
613 ALLOCATE (
CHARACTER(LEN=libxc_get_reference_length(xc_fun, lsd)) :: reference)
615 ALLOCATE (
CHARACTER(LEN=20*default_string_length) :: reference)
618 WRITE (iw, fmt=
"(' FUNCTIONAL| ',a,':')") &
619 trim(xc_fun%section%name)
620 DO il = 1, len_trim(reference), 67
621 WRITE (iw, fmt=
"(' FUNCTIONAL| ',a67)") reference(il:)
623 DEALLOCATE (reference)
627 IF (iw > 0)
WRITE (iw, fmt=
"(' FUNCTIONAL| NO EXCHANGE-CORRELATION FUNCTIONAL USED.')")
646 INTEGER,
INTENT(IN) :: iw
648 CHARACTER(len=60) :: pline
649 INTEGER :: i, j, k, l
651 SELECT CASE (potential%ppot_type)
653 WRITE (iw,
'(/," ",28("*"),A,27("*"))')
" All Electron Potential "
655 WRITE (iw,
'(/," ",29("*"),A,29("*"))')
" GTH Pseudopotential "
656 WRITE (iw,
'(T10,A,T76,F5.1)')
" Core Charge ", potential%gth_pot%zion
657 WRITE (iw,
'(T10,A,T66,F15.6)')
" Rc ", potential%gth_pot%rc
658 WRITE (pline,
'(5F12.6)') (potential%gth_pot%cl(i), i=1, potential%gth_pot%ncl)
659 WRITE (iw,
'(T10,A,T21,A60)')
" C1 C2 ... ", adjustr(pline)
660 IF (potential%gth_pot%lpotextended)
THEN
661 DO k = 1, potential%gth_pot%nexp_lpot
662 WRITE (iw,
'(T10,A,F10.6,T38,A,4F10.6)')
" LPot: rc=", potential%gth_pot%alpha_lpot(k), &
663 "CX=", (potential%gth_pot%cval_lpot(i, k), i=1, potential%gth_pot%nct_lpot(k))
666 IF (potential%gth_pot%nlcc)
THEN
667 DO k = 1, potential%gth_pot%nexp_nlcc
668 WRITE (iw,
'(T10,A,F10.6,T38,A,4F10.6)')
" LSDPot: rc=", potential%gth_pot%alpha_nlcc(k), &
669 "CX=", (potential%gth_pot%cval_nlcc(i, k)*4.0_dp*
pi, i=1, potential%gth_pot%nct_nlcc(k))
672 IF (potential%gth_pot%lsdpot)
THEN
673 DO k = 1, potential%gth_pot%nexp_lsd
674 WRITE (iw,
'(T10,A,F10.6,T38,A,4F10.6)')
" LSDPot: rc=", potential%gth_pot%alpha_lsd(k), &
675 "CX=", (potential%gth_pot%cval_lsd(i, k), i=1, potential%gth_pot%nct_lsd(k))
679 IF (potential%gth_pot%nl(l) > 0)
THEN
680 WRITE (iw,
'(T10,A,T76,I5)')
" Angular momentum ", l
681 WRITE (iw,
'(T10,A,T66,F15.6)')
" Rcnl ", potential%gth_pot%rcnl(l)
682 WRITE (iw,
'(T10,A,T76,I5)')
" Nl ", potential%gth_pot%nl(l)
683 WRITE (pline,
'(5F12.6)') (potential%gth_pot%hnl(1, j, l), j=1, potential%gth_pot%nl(l))
684 WRITE (iw,
'(T10,A,T21,A60)')
" Hnl ", adjustr(pline)
685 DO i = 2, potential%gth_pot%nl(l)
686 WRITE (pline,
'(T21,5F12.6)') (potential%gth_pot%hnl(i, j, l), j=i, potential%gth_pot%nl(l))
687 WRITE (iw,
'(T21,A60)') adjustr(pline)
691 IF (potential%gth_pot%soc)
THEN
692 WRITE (iw,
'(T10,A)')
" Spin-orbit coupling parameters "
694 IF (potential%gth_pot%nl(l) > 0)
THEN
695 WRITE (iw,
'(T10,A,T76,I5)')
" Angular momentum ", l
696 WRITE (iw,
'(T10,A,T66,F15.6)')
" Rcnl ", potential%gth_pot%rcnl(l)
697 WRITE (iw,
'(T10,A,T76,I5)')
" Nl ", potential%gth_pot%nl(l)
698 WRITE (pline,
'(5F12.6)') (potential%gth_pot%knl(1, j, l), j=1, potential%gth_pot%nl(l))
699 WRITE (iw,
'(T10,A,T21,A60)')
" Hnl ", adjustr(pline)
700 DO i = 2, potential%gth_pot%nl(l)
701 WRITE (pline,
'(T21,5F12.6)') (potential%gth_pot%knl(i, j, l), j=i, potential%gth_pot%nl(l))
702 WRITE (iw,
'(T21,A60)') adjustr(pline)
708 WRITE (iw,
'(/," ",29("*"),A,29("*"))')
" UPF Pseudopotential "
709 DO k = 1, potential%upf_pot%maxinfo
710 WRITE (iw,
'(A80)') potential%upf_pot%info(k)
713 WRITE (iw,
'(/," ",29("*"),A,29("*"))')
" SGP Pseudopotential "
714 WRITE (iw,
'(T10,A,T76,F5.1)')
" Core Charge ", potential%sgp_pot%zion
716 WRITE (iw,
'(/," ",26("*"),A,27("*"))')
" Effective Core Potential "
717 WRITE (iw,
'(T10,A,T76,F5.1)')
" Core Charge ", potential%ecp_pot%zion
718 DO k = 1, potential%ecp_pot%nloc
720 WRITE (iw,
'(T10,A,T40,I3,T49,2F16.8)')
" Local Potential ", potential%ecp_pot%nrloc(k), &
721 potential%ecp_pot%bloc(k), potential%ecp_pot%aloc(k)
723 WRITE (iw,
'(T40,I3,T49,2F16.8)') potential%ecp_pot%nrloc(k), &
724 potential%ecp_pot%bloc(k), potential%ecp_pot%aloc(k)
727 DO l = 0, potential%ecp_pot%lmax
728 WRITE (iw,
'(T10,A,I3)')
" ECP l-value ", l
729 DO k = 1, potential%ecp_pot%npot(l)
730 WRITE (iw,
'(T40,I3,T49,2F16.8)') potential%ecp_pot%nrpot(k, l), &
731 potential%ecp_pot%bpot(k, l), potential%ecp_pot%apot(k, l)
737 IF (potential%confinement)
THEN
738 IF (potential%conf_type ==
poly_conf)
THEN
739 WRITE (iw,
'(/,T10,A,T51,F12.6," * (R /",F6.2,")**",F6.2)') &
740 " Confinement Potential ", potential%acon, potential%rcon, potential%scon
742 WRITE (iw,
'(/,T10,A)')
" Confinement Potential s*F[(r-ron)/w] "
743 WRITE (iw,
'(T57,A,F12.6,A)')
"s =", potential%acon,
" Ha"
744 WRITE (iw,
'(T57,A,F12.6,A)')
"w =", potential%rcon,
" Bohr"
745 WRITE (iw,
'(T57,A,F12.6,A)')
"ron =", potential%scon,
" Bohr"
750 WRITE (iw,
'(/,T10,A)')
" No Confinement Potential is applied "
752 WRITE (iw,
'(" ",79("*"))')
768 INTEGER,
INTENT(IN),
OPTIONAL :: iunit
769 REAL(kind=
dp),
INTENT(IN),
OPTIONAL :: fopt
771 INTEGER :: i, iw, j, k, n
773 IF (
PRESENT(iunit))
THEN
776 CALL open_file(file_name=
"GTH-PARAMETER", file_status=
"UNKNOWN", file_action=
"WRITE", unit_number=iw)
778 IF (
PRESENT(fopt))
THEN
779 WRITE (iw,
'(A,F30.8)')
"# "//trim(adjustl(gthpot%symbol)), fopt
781 WRITE (iw,
'(A)') trim(adjustl(gthpot%symbol))//
" "//trim(adjustl(gthpot%pname))
783 WRITE (iw,
'(4I5)') gthpot%econf(0:3)
784 WRITE (iw,
'(F20.14,I8,5F20.14)') gthpot%rc, gthpot%ncl, (gthpot%cl(i), i=1, gthpot%ncl)
785 IF (gthpot%lpotextended)
THEN
786 WRITE (iw,
'(A,I5)')
" LPOT", gthpot%nexp_lpot
787 DO i = 1, gthpot%nexp_lpot
788 WRITE (iw,
'(F20.14,I8,5F20.14)') gthpot%alpha_lpot(i), gthpot%nct_lpot(i), &
789 (gthpot%cval_lpot(j, i), j=1, gthpot%nct_lpot(i))
792 IF (gthpot%lsdpot)
THEN
793 WRITE (iw,
'(A,I5)')
" LSD ", gthpot%nexp_lsd
794 DO i = 1, gthpot%nexp_lsd
795 WRITE (iw,
'(F20.14,I8,5F20.14)') gthpot%alpha_lsd(i), gthpot%nct_lsd(i), &
796 (gthpot%cval_lsd(j, i), j=1, gthpot%nct_lsd(i))
799 IF (gthpot%nlcc)
THEN
800 WRITE (iw,
'(A,I5)')
" NLCC ", gthpot%nexp_nlcc
801 DO i = 1, gthpot%nexp_nlcc
802 WRITE (iw,
'(F20.14,I8,5F20.14)') gthpot%alpha_nlcc(i), gthpot%nct_nlcc(i), &
803 (gthpot%cval_nlcc(j, i)*4.0_dp*
pi, j=1, gthpot%nct_nlcc(i))
808 IF (gthpot%nl(i) > 0)
THEN
815 WRITE (iw,
'(F20.14,I8,5F20.14)') gthpot%rcnl(i), gthpot%nl(i), (gthpot%hnl(1, k, i), k=1, gthpot%nl(i))
816 SELECT CASE (gthpot%nl(i))
818 WRITE (iw,
'(T49,F20.14)') gthpot%hnl(2, 2, i)
820 WRITE (iw,
'(T49,2F20.14)') gthpot%hnl(2, 2, i), gthpot%hnl(2, 3, i)
821 WRITE (iw,
'(T69,F20.14)') gthpot%hnl(3, 3, i)
823 DO j = 2, gthpot%nl(i)
824 WRITE (iw,
'(T29,5F20.14)') (gthpot%hnl(j, k, i), k=j, gthpot%nl(i))
830 WRITE (iw,
'(T29,5F20.14)') (gthpot%hnl(1, k, i), k=1, gthpot%nl(i))
831 SELECT CASE (gthpot%nl(i))
833 WRITE (iw,
'(T49,F20.14)') gthpot%knl(2, 2, i)
835 WRITE (iw,
'(T49,2F20.14)') gthpot%knl(2, 2, i), gthpot%knl(2, 3, i)
836 WRITE (iw,
'(T69,F20.14)') gthpot%knl(3, 3, i)
838 DO j = 2, gthpot%nl(i)
839 WRITE (iw,
'(T29,5F20.14)') (gthpot%knl(j, k, i), k=j, gthpot%nl(i))
844 IF (.NOT.
PRESENT(iunit))
CALL close_file(unit_number=iw)
858 INTEGER,
INTENT(IN) :: iw
859 LOGICAL,
INTENT(IN),
OPTIONAL ::
xmgrace
861 CHARACTER(LEN=40) :: fnbody
865 SELECT CASE (
atom%method_type)
869 CALL atom_print_orbitals_helper(
atom,
atom%orbitals%wfn,
"", iw)
871 CALL atom_print_orbitals_helper(
atom,
atom%orbitals%wfna,
"Alpha", iw)
872 CALL atom_print_orbitals_helper(
atom,
atom%orbitals%wfnb,
"Beta", iw)
874 CALL atom_print_orbitals_helper(
atom,
atom%orbitals%wfn,
"", iw)
876 CALL atom_print_orbitals_helper(
atom,
atom%orbitals%wfna,
"Alpha", iw)
877 CALL atom_print_orbitals_helper(
atom,
atom%orbitals%wfnb,
"Beta", iw)
884 IF (
graph .AND. iw > 0)
THEN
886 fnbody = trim(
ptable(z)%symbol)//
"_PPorbital"
887 SELECT CASE (
atom%method_type)
891 CALL atom_orbitals_grace(
atom,
atom%orbitals%wfn, fnbody)
893 CALL atom_orbitals_grace(
atom,
atom%orbitals%wfna, trim(fnbody)//
"alpha")
894 CALL atom_orbitals_grace(
atom,
atom%orbitals%wfnb, trim(fnbody)//
"beta")
896 CALL atom_orbitals_grace(
atom,
atom%orbitals%wfn, fnbody)
898 CALL atom_orbitals_grace(
atom,
atom%orbitals%wfna, trim(fnbody)//
"alpha")
899 CALL atom_orbitals_grace(
atom,
atom%orbitals%wfnb, trim(fnbody)//
"beta")
916 SUBROUTINE atom_print_orbitals_helper(atom, wfn, description, iw)
918 REAL(kind=
dp),
DIMENSION(:, :, 0:),
INTENT(INOUT) :: wfn
919 CHARACTER(len=*),
INTENT(IN) :: description
920 INTEGER,
INTENT(IN) :: iw
922 INTEGER :: b, l, maxl, nb, nv, v
924 WRITE (iw,
'(/,A,A,A)')
" Atomic orbital expansion coefficients [", description,
"]"
926 maxl =
atom%state%maxl_calc
929 nb =
atom%basis%nbas(l)
930 nv =
atom%state%maxn_calc(l)
931 IF (nb > 0 .AND. nv > 0)
THEN
932 nv = min(nv,
SIZE(wfn, 2))
934 WRITE (iw,
'(/," ORBITAL L = ",I1," State = ",I3)') l, v
936 WRITE (iw,
'(" ",ES23.15)') wfn(b, v, l)
941 END SUBROUTINE atom_print_orbitals_helper
951 SUBROUTINE atom_orbitals_grace(atom, wfn, fnbody)
953 REAL(kind=
dp),
DIMENSION(:, :, 0:),
INTENT(INOUT) :: wfn
954 CHARACTER(len=*),
INTENT(IN) :: fnbody
956 CHARACTER(LEN=1),
DIMENSION(0:8) :: lname
957 CHARACTER(LEN=1),
DIMENSION(1:9) :: wnum
958 CHARACTER(LEN=40) :: fname, legend
959 INTEGER :: b, i, iw, l, m, maxl, nb, nv, v
960 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:, :) :: gdata, wfnr
961 REAL(kind=
dp),
DIMENSION(4) :: world_coord
963 lname = [
's',
'p',
'd',
'f',
'g',
'h',
'j',
'k',
'l']
964 wnum = [
'1',
'2',
'3',
'4',
'5',
'6',
'7',
'8',
'9']
965 m =
atom%basis%grid%nr
966 maxl =
atom%state%maxl_calc
968 fname = trim(fnbody)//
"_"//lname(l)//
".agr"
969 nb =
atom%basis%nbas(l)
970 nv =
atom%state%maxn_calc(l)
971 IF (nb > 0 .AND. nv > 0)
THEN
972 CALL open_file(file_name=fname, file_status=
"UNKNOWN", file_action=
"WRITE", unit_number=iw)
973 nv = min(nv,
SIZE(wfn, 2))
974 ALLOCATE (wfnr(m, nv))
978 wfnr(:, v) = wfnr(:, v) + wfn(b, v, l)*
atom%basis%bf(:, b, l)
981 world_coord(1) = 0.0_dp
982 world_coord(2) = minval(wfnr) - 0.5_dp
983 world_coord(3) = 15.0_dp
984 world_coord(4) = maxval(wfnr) + 0.5_dp
989 title=
"PP Radial Wavefunction", &
990 subtitle=lname(l)//
"-Quantum Number", &
991 xlabel=
"Radius [Bohr]", &
994 legend =
"WFN "//wnum(i + 1)
997 ALLOCATE (gdata(m, 2))
998 gdata(1:m, 1) =
atom%basis%grid%rad(1:m)
1000 gdata(1:m, 2) = wfnr(1:m, i + 1)
1003 DEALLOCATE (gdata, wfnr)
1007 END SUBROUTINE atom_orbitals_grace
program graph
Program to Map on grid the hills spawned during a metadynamics run.
Routines that print various information about an atomic kind.
subroutine, public atom_print_energies(atom, iw)
Print energy components.
subroutine, public atom_print_orbitals(atom, iw, xmgrace)
Print atomic orbitals.
subroutine, public atom_print_basis(atom_basis, iw, title)
Print atomic basis set.
subroutine, public atom_print_iteration(iter, deps, etot, iw)
Print convergence information.
subroutine, public atom_print_basis_file(atom_basis, wfn)
Print the optimized atomic basis set into a file.
subroutine, public atom_print_state(state, iw)
Print information about electronic state.
subroutine, public atom_print_zmp_iteration(iter, deps, atom, iw)
Printing of the atomic iterations when ZMP is active.
subroutine, public atom_print_method(atom, iw)
Print information about the electronic structure method in use.
subroutine, public atom_write_pseudo_param(gthpot, iunit, fopt)
Print GTH pseudo-potential parameters.
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.
Define the atom type and its sub types.
integer, parameter, public num_basis
integer, parameter, public cgto_basis
integer, parameter, public gto_basis
integer, parameter, public sto_basis
integer, parameter, public lmat
Some basic routines for atomic calculations.
pure integer function, dimension(0:lmat), public get_maxn_occ(occupation)
Return the maximum principal quantum number of occupied orbitals.
subroutine, public get_rho0(atom, rho0)
Calculate the total electron density at R=0.
pure integer function, public get_maxl_occ(occupation)
Return the maximum orbital quantum number of occupied orbitals.
Utility routines to open and close files. Tracking of preconnections.
subroutine, public open_file(file_name, file_status, file_form, file_action, file_position, file_pad, unit_number, debug, skip_get_unit_number, file_access)
Opens the requested file using a free unit number.
subroutine, public close_file(unit_number, file_status, keep_preconnection)
Close an open file given by its logical unit number. Optionally, keep the file and unit preconnected.
Defines the basic variable types.
integer, parameter, public dp
integer, parameter, public default_string_length
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
Periodic Table related data definitions.
type(atom), dimension(0:nelem), public ptable
Definition of physical constants:
real(kind=dp), parameter, public evolt
subroutine, public xc_functional_get_info(functional, lsd, reference, shortform, needs, max_deriv, print_warn)
get the information about the given functional
calculates a functional from libxc and its derivatives
logical function, public libxc_check_existence_in_libxc(libxc_params)
This function checks whether a functional name belongs to LibXC.
integer function, public libxc_get_reference_length(libxc_params, lsd)
This function returns the maximum length of the reference string for a given LibXC functional.
Routines to facilitate writing XMGRACE files.
subroutine, public xm_graph_data(iw, gnum, gdata)
...
subroutine, public xm_write_frameport(iw)
...
subroutine, public xm_write_frame(iw, wcoord, title, subtitle, xlabel, ylabel)
...
subroutine, public xm_write_defaults(iw)
...
subroutine, public xm_graph_info(iw, gnum, linewidth, legend)
...
Provides all information about a basis set.
Provides all information about a pseudopotential.
Provides all information on states and occupation.
Provides all information about an atomic kind.