(git:34ef472)
atom_output.F
Go to the documentation of this file.
1 !--------------------------------------------------------------------------------------------------!
2 ! CP2K: A general program to perform molecular dynamics simulations !
3 ! Copyright 2000-2024 CP2K developers group <https://cp2k.org> !
4 ! !
5 ! SPDX-License-Identifier: GPL-2.0-or-later !
6 !--------------------------------------------------------------------------------------------------!
7 
8 ! **************************************************************************************************
9 !> \brief Routines that print various information about an atomic kind.
10 ! **************************************************************************************************
12  USE atom_types, ONLY: &
13  atom_basis_type, atom_gthpot_type, atom_potential_type, atom_state, atom_type, cgto_basis, &
14  ecp_pseudo, gth_pseudo, gto_basis, lmat, no_pseudo, num_basis, sgp_pseudo, sto_basis, &
15  upf_pseudo
16  USE atom_utils, ONLY: get_maxl_occ,&
17  get_maxn_occ,&
18  get_rho0
19  USE cp_files, ONLY: close_file,&
20  open_file
21  USE input_constants, ONLY: &
28  section_vals_type,&
30  USE kinds, ONLY: default_string_length,&
31  dp
32  USE mathconstants, ONLY: dfac,&
33  pi,&
34  rootpi
35  USE periodic_table, ONLY: ptable
36  USE physcon, ONLY: evolt
40 #include "./base/base_uses.f90"
41 
42  IMPLICIT NONE
43 
44  PRIVATE
45 
46  CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'atom_output'
47 
52 
53 CONTAINS
54 
55 ! **************************************************************************************************
56 !> \brief Print an information string related to the atomic kind.
57 !> \param zval atomic number
58 !> \param info information string
59 !> \param iw output file unit
60 !> \par History
61 !> * 09.2008 created [Juerg Hutter]
62 ! **************************************************************************************************
63  SUBROUTINE atom_print_info(zval, info, iw)
64  INTEGER, INTENT(IN) :: zval
65  CHARACTER(len=*), INTENT(IN) :: info
66  INTEGER, INTENT(IN) :: iw
67 
68  WRITE (iw, '(/," ",A,T40,A," [",A,"]",T62,"Atomic number:",T78,I3,/)') &
69  adjustl(trim(info)), trim(ptable(zval)%name), trim(ptable(zval)%symbol), zval
70 
71  END SUBROUTINE atom_print_info
72 
73 ! **************************************************************************************************
74 !> \brief Print information about electronic state.
75 !> \param state electronic state
76 !> \param iw output file unit
77 !> \par History
78 !> * 02.2010 unrestricted KS and HF methods [Juerg Hutter]
79 !> * 11.2009 print multiplicity [Juerg Hutter]
80 !> * 08.2008 created [Juerg Hutter]
81 ! **************************************************************************************************
82  SUBROUTINE atom_print_state(state, iw)
83  TYPE(atom_state) :: state
84  INTEGER, INTENT(IN) :: iw
85 
86  CHARACTER(LEN=1), DIMENSION(0:7), PARAMETER :: &
87  label = (/"S", "P", "D", "F", "G", "H", "I", "K"/)
88 
89  INTEGER :: j, l, mc, mlc, mlo, mm(0:lmat), mo
90 
91  cpassert(lmat <= 7)
92  WRITE (iw, '(/,T2,A)') "Electronic structure"
93  WRITE (iw, '(T5,A,T71,F10.2)') "Total number of core electrons", sum(state%core)
94  WRITE (iw, '(T5,A,T71,F10.2)') "Total number of valence electrons", sum(state%occ)
95  WRITE (iw, '(T5,A,T71,F10.2)') "Total number of electrons", sum(state%occ + state%core)
96  SELECT CASE (state%multiplicity)
97  CASE (-1)
98  WRITE (iw, '(T5,A,T68,A)') "Multiplicity", "not specified"
99  CASE (-2)
100  WRITE (iw, '(T5,A,T72,A)') "Multiplicity", "high spin"
101  CASE (-3)
102  WRITE (iw, '(T5,A,T73,A)') "Multiplicity", "low spin"
103  CASE (1)
104  WRITE (iw, '(T5,A,T74,A)') "Multiplicity", "singlet"
105  CASE (2)
106  WRITE (iw, '(T5,A,T74,A)') "Multiplicity", "doublet"
107  CASE (3)
108  WRITE (iw, '(T5,A,T74,A)') "Multiplicity", "triplet"
109  CASE (4)
110  WRITE (iw, '(T5,A,T74,A)') "Multiplicity", "quartet"
111  CASE (5)
112  WRITE (iw, '(T5,A,T74,A)') "Multiplicity", "quintet"
113  CASE (6)
114  WRITE (iw, '(T5,A,T75,A)') "Multiplicity", "sextet"
115  CASE (7)
116  WRITE (iw, '(T5,A,T75,A)') "Multiplicity", "septet"
117  CASE DEFAULT
118  END SELECT
119 
120  mlo = get_maxl_occ(state%occ)
121  mlc = get_maxl_occ(state%core)
122  mm = get_maxn_occ(state%core)
123 
124  IF (state%multiplicity == -1) THEN
125  DO l = 0, max(mlo, mlc)
126  mo = state%maxn_occ(l)
127  IF (sum(state%core(l, :)) == 0) THEN
128  WRITE (iw, '(A5,T10,10F6.2)') label(l), (state%occ(l, j), j=1, mo)
129  ELSE
130  mc = mm(l)
131  cpassert(sum(state%occ(l, 1:mc)) == 0)
132  WRITE (iw, advance="no", fmt='(A5,T9,A1,10F6.2)') label(l), "[", (state%core(l, j), j=1, mc)
133  WRITE (iw, fmt='(A1,F5.2,10F6.2)') "]", (state%occ(l, j), j=mc + 1, mc + mo)
134  END IF
135  END DO
136  ELSE
137  WRITE (iw, '(T5,A)') "Alpha Electrons"
138  DO l = 0, max(mlo, mlc)
139  mo = state%maxn_occ(l)
140  IF (sum(state%core(l, :)) == 0) THEN
141  WRITE (iw, '(A5,T10,10F6.2)') label(l), (state%occa(l, j), j=1, mo)
142  ELSE
143  mc = mm(l)
144  WRITE (iw, advance="no", fmt='(A5,T9,A1,10F6.2)') label(l), "[", (0.5_dp*state%core(l, j), j=1, mc)
145  WRITE (iw, fmt='(A1,F5.2,10F6.2)') "]", (state%occa(l, j), j=1, mo)
146  END IF
147  END DO
148  WRITE (iw, '(T5,A)') "Beta Electrons"
149  DO l = 0, max(mlo, mlc)
150  mo = state%maxn_occ(l)
151  IF (sum(state%core(l, :)) == 0) THEN
152  WRITE (iw, '(A5,T10,10F6.2)') label(l), (state%occb(l, j), j=1, mo)
153  ELSE
154  mc = mm(l)
155  WRITE (iw, advance="no", fmt='(A5,T9,A1,10F6.2)') label(l), "[", (0.5_dp*state%core(l, j), j=1, mc)
156  WRITE (iw, fmt='(A1,F5.2,10F6.2)') "]", (state%occb(l, j), j=1, mo)
157  END IF
158  END DO
159  END IF
160  WRITE (iw, *)
161 
162  END SUBROUTINE atom_print_state
163 
164 ! **************************************************************************************************
165 !> \brief Print energy components.
166 !> \param atom information about the atomic kind
167 !> \param iw output file unit
168 !> \par History
169 !> * 05.2010 print virial coefficient [Juerg Hutter]
170 !> * 02.2010 unrestricted KS and HF methods [Juerg Hutter]
171 !> * 09.2008 print orbital energies [Juerg Hutter]
172 !> * 08.2008 created [Juerg Hutter]
173 ! **************************************************************************************************
174  SUBROUTINE atom_print_energies(atom, iw)
175  TYPE(atom_type) :: atom
176  INTEGER, INTENT(IN) :: iw
177 
178  INTEGER :: i, l, n
179  REAL(kind=dp) :: drho
180 
181  WRITE (iw, '(/,A,T36,A,T61,F20.12)') " Energy components [Hartree]", &
182  " Total Energy ::", atom%energy%etot
183  WRITE (iw, '(T36,A,T61,F20.12)') " Band Energy ::", atom%energy%eband
184  WRITE (iw, '(T36,A,T61,F20.12)') " Kinetic Energy ::", atom%energy%ekin
185  WRITE (iw, '(T36,A,T61,F20.12)') "Potential Energy ::", atom%energy%epot
186  IF (atom%energy%ekin /= 0.0_dp) THEN
187  WRITE (iw, '(T36,A,T61,F20.12)') " Virial (-V/T) ::", -atom%energy%epot/atom%energy%ekin
188  END IF
189  WRITE (iw, '(T36,A,T61,F20.12)') " Core Energy ::", atom%energy%ecore
190  IF (atom%energy%exc /= 0._dp) &
191  WRITE (iw, '(T36,A,T61,F20.12)') " XC Energy ::", atom%energy%exc
192  WRITE (iw, '(T36,A,T61,F20.12)') " Coulomb Energy ::", atom%energy%ecoulomb
193  IF (atom%energy%eexchange /= 0._dp) &
194  WRITE (iw, '(T34,A,T61,F20.12)') "HF Exchange Energy ::", atom%energy%eexchange
195  IF (atom%potential%ppot_type /= no_pseudo) THEN
196  WRITE (iw, '(T20,A,T61,F20.12)') " Total Pseudopotential Energy ::", atom%energy%epseudo
197  WRITE (iw, '(T20,A,T61,F20.12)') " Local Pseudopotential Energy ::", atom%energy%eploc
198  IF (atom%energy%elsd /= 0._dp) &
199  WRITE (iw, '(T20,A,T61,F20.12)') " Local Spin-potential Energy ::", atom%energy%elsd
200  WRITE (iw, '(T20,A,T61,F20.12)') " Nonlocal Pseudopotential Energy ::", atom%energy%epnl
201  END IF
202  IF (atom%potential%confinement) THEN
203  WRITE (iw, '(T36,A,T61,F20.12)') " Confinement ::", atom%energy%econfinement
204  END IF
205 
206  IF (atom%state%multiplicity == -1) THEN
207  WRITE (iw, '(/,A,T20,A,T30,A,T36,A,T49,A,T71,A,/)') " Orbital energies", &
208  "State", "L", "Occupation", "Energy[a.u.]", "Energy[eV]"
209  DO l = 0, atom%state%maxl_calc
210  n = atom%state%maxn_calc(l)
211  DO i = 1, n
212  WRITE (iw, '(T23,I2,T30,I1,T36,F10.3,T46,F15.6,T66,F15.6)') &
213  i, l, atom%state%occupation(l, i), atom%orbitals%ener(i, l), atom%orbitals%ener(i, l)*evolt
214  END DO
215  IF (n > 0) WRITE (iw, *)
216  END DO
217  ELSE
218  WRITE (iw, '(/,A,T20,A,T30,A,T36,A,T42,A,T55,A,T71,A,/)') " Orbital energies", &
219  "State", "Spin", "L", "Occupation", "Energy[a.u.]", "Energy[eV]"
220  DO l = 0, atom%state%maxl_calc
221  n = atom%state%maxn_calc(l)
222  DO i = 1, n
223  WRITE (iw, '(T23,I2,T29,A,T36,I1,T42,F10.3,T52,F15.6,T68,F13.6)') &
224  i, "alpha", l, atom%state%occa(l, i), atom%orbitals%enera(i, l), atom%orbitals%enera(i, l)*evolt
225  END DO
226  DO i = 1, n
227  WRITE (iw, '(T23,I2,T29,A,T36,I1,T42,F10.3,T52,F15.6,T68,F13.6)') &
228  i, " beta", l, atom%state%occb(l, i), atom%orbitals%enerb(i, l), atom%orbitals%enerb(i, l)*evolt
229  END DO
230  IF (n > 0) WRITE (iw, *)
231  END DO
232  END IF
233 
234  CALL get_rho0(atom, drho)
235  WRITE (iw, '(/,A,T66,F15.6)') " Total Electron Density at R=0: ", drho
236 
237  END SUBROUTINE atom_print_energies
238 
239 ! **************************************************************************************************
240 !> \brief Printing of the atomic iterations when ZMP is active.
241 !> \param iter current iteration number
242 !> \param deps convergence
243 !> \param atom intormation about the atomic kind
244 !> \param iw output file unit
245 !> \author D. Varsano [daniele.varsano@nano.cnr.it]
246 ! **************************************************************************************************
247  SUBROUTINE atom_print_zmp_iteration(iter, deps, atom, iw)
248  INTEGER, INTENT(IN) :: iter
249  REAL(dp), INTENT(IN) :: deps
250  TYPE(atom_type), INTENT(IN) :: atom
251  INTEGER, INTENT(IN) :: iw
252 
253  IF (iter == 1) THEN
254  WRITE (iw, '(/," ",79("*"),/,T33,"Integral",T48,"Integral",/,T3,A,T16,A,T33,A,T46,A,T69,A/," ",79("*"))') &
255  "Iteration", "Convergence", "rho diff.", "rho*v_xc[au]", "Energy[au]"
256  END IF
257  WRITE (iw, '(T3,I9,T15,G13.6,T30,G13.6,T46,G13.6,T61,F20.12)') iter, deps, atom%rho_diff_integral, &
258  atom%energy%exc, atom%energy%etot
259 
260  END SUBROUTINE atom_print_zmp_iteration
261 
262 ! **************************************************************************************************
263 !> \brief Print convergence information.
264 !> \param iter current iteration number
265 !> \param deps convergency
266 !> \param etot total energy
267 !> \param iw output file unit
268 !> \par History
269 !> * 08.2008 created [Juerg Hutter]
270 ! **************************************************************************************************
271  SUBROUTINE atom_print_iteration(iter, deps, etot, iw)
272  INTEGER, INTENT(IN) :: iter
273  REAL(dp), INTENT(IN) :: deps, etot
274  INTEGER, INTENT(IN) :: iw
275 
276  IF (iter == 1) THEN
277  WRITE (iw, '(/," ",79("*"),/,T19,A,T38,A,T70,A,/," ",79("*"))') &
278  "Iteration", "Convergence", "Energy [au]"
279  END IF
280  WRITE (iw, '(T20,i8,T34,G14.6,T61,F20.12)') iter, deps, etot
281 
282  END SUBROUTINE atom_print_iteration
283 
284 ! **************************************************************************************************
285 !> \brief Print atomic basis set.
286 !> \param atom_basis atomic basis set
287 !> \param iw output file unit
288 !> \param title header to print on top of the basis set
289 !> \par History
290 !> * 09.2008 created [Juerg Hutter]
291 ! **************************************************************************************************
292  SUBROUTINE atom_print_basis(atom_basis, iw, title)
293  TYPE(atom_basis_type) :: atom_basis
294  INTEGER, INTENT(IN) :: iw
295  CHARACTER(len=*) :: title
296 
297  INTEGER :: i, j, l
298 
299  WRITE (iw, '(/,A)') trim(title)
300  SELECT CASE (atom_basis%basis_type)
301  CASE (gto_basis)
302  IF (atom_basis%geometrical) THEN
303  WRITE (iw, '(/," ",21("*"),A,22("*"))') " Geometrical Gaussian Type Orbitals "
304  WRITE (iw, '(A,F15.8,T41,A,F15.8)') " Initial exponent: ", atom_basis%aval, &
305  " Proportionality factor: ", atom_basis%cval
306  ELSE
307  WRITE (iw, '(/," ",21("*"),A,21("*"))') " Uncontracted Gaussian Type Orbitals "
308  END IF
309  DO l = 0, lmat
310  IF (atom_basis%nbas(l) > 0) THEN
311  SELECT CASE (l)
312  CASE DEFAULT
313  WRITE (iw, '(/,T2,A,(T30,I5,T51,F30.8))') &
314  "X Exponents: ", (i, atom_basis%am(i, l), i=1, atom_basis%nbas(l))
315  CASE (0)
316  WRITE (iw, '(/,T2,A,(T30,I5,T51,F30.8))') &
317  "s Exponents: ", (i, atom_basis%am(i, 0), i=1, atom_basis%nbas(0))
318  CASE (1)
319  WRITE (iw, '(/,T2,A,(T30,I5,T51,F30.8))') &
320  "p Exponents: ", (i, atom_basis%am(i, 1), i=1, atom_basis%nbas(1))
321  CASE (2)
322  WRITE (iw, '(/,T2,A,(T30,I5,T51,F30.8))') &
323  "d Exponents: ", (i, atom_basis%am(i, 2), i=1, atom_basis%nbas(2))
324  CASE (3)
325  WRITE (iw, '(/,T2,A,(T30,I5,T51,F30.8))') &
326  "f Exponents: ", (i, atom_basis%am(i, 3), i=1, atom_basis%nbas(3))
327  END SELECT
328  END IF
329  END DO
330  WRITE (iw, '(" ",79("*"))')
331  CASE (cgto_basis)
332  WRITE (iw, '(/," ",22("*"),A,22("*"))') " Contracted Gaussian Type Orbitals "
333  DO l = 0, lmat
334  IF (atom_basis%nbas(l) > 0) THEN
335  IF (l == 0) WRITE (iw, '(A)') " s Functions"
336  IF (l == 1) WRITE (iw, '(A)') " p Functions"
337  IF (l == 2) WRITE (iw, '(A)') " d Functions"
338  IF (l == 3) WRITE (iw, '(A)') " f Functions"
339  IF (l >= 3) WRITE (iw, '(A)') " x Functions"
340  DO i = 1, atom_basis%nprim(l)
341  WRITE (iw, '(F15.6,5(T21,6F10.6,/))') &
342  atom_basis%am(i, l), (atom_basis%cm(i, j, l), j=1, atom_basis%nbas(l))
343  END DO
344  END IF
345  END DO
346  WRITE (iw, '(" ",79("*"))')
347  CASE (sto_basis)
348  WRITE (iw, '(/," ",28("*"),A,29("*"))') " Slater Type Orbitals "
349  DO l = 0, lmat
350  DO i = 1, atom_basis%nbas(l)
351  SELECT CASE (l)
352  CASE DEFAULT
353  WRITE (iw, '(T10,I1,A,T40,F25.12)') atom_basis%ns(i, l), "X Exponent :", atom_basis%as(i, l)
354  CASE (0)
355  WRITE (iw, '(T10,I1,A,T40,F25.12)') atom_basis%ns(i, 0), "S Exponent :", atom_basis%as(i, 0)
356  CASE (1)
357  WRITE (iw, '(T10,I1,A,T40,F25.12)') atom_basis%ns(i, 1), "P Exponent :", atom_basis%as(i, 1)
358  CASE (2)
359  WRITE (iw, '(T10,I1,A,T40,F25.12)') atom_basis%ns(i, 2), "D Exponent :", atom_basis%as(i, 2)
360  CASE (3)
361  WRITE (iw, '(T10,I1,A,T40,F25.12)') atom_basis%ns(i, 3), "F Exponent :", atom_basis%as(i, 3)
362  END SELECT
363  END DO
364  END DO
365  WRITE (iw, '(" ",79("*"))')
366  CASE (num_basis)
367  cpabort("")
368  CASE DEFAULT
369  cpabort("")
370  END SELECT
371 
372  END SUBROUTINE atom_print_basis
373 
374 ! **************************************************************************************************
375 !> \brief Print the optimized atomic basis set into a file.
376 !> \param atom_basis atomic basis set
377 !> \param wfn ...
378 !> \par History
379 !> * 11.2016 revised output format [Matthias Krack]
380 !> * 11.2011 Slater basis functions [Juerg Hutter]
381 !> * 03.2011 created [Juerg Hutter]
382 !> \note The basis set is stored as the file 'OPT_BASIS' inside the current working directory.
383 !> It may be a good idea, however, to specify the name of this file via some input section.
384 ! **************************************************************************************************
385  SUBROUTINE atom_print_basis_file(atom_basis, wfn)
386  TYPE(atom_basis_type) :: atom_basis
387  REAL(kind=dp), DIMENSION(:, :, 0:), OPTIONAL :: wfn
388 
389  INTEGER :: i, im, iw, l
390  REAL(kind=dp) :: expzet, prefac, zeta
391 
392  CALL open_file(file_name="OPT_BASIS", file_status="UNKNOWN", file_action="WRITE", unit_number=iw)
393  SELECT CASE (atom_basis%basis_type)
394  CASE (gto_basis)
395  IF (atom_basis%geometrical) THEN
396  WRITE (iw, '(/," ",21("*"),A,22("*"))') " Geometrical Gaussian Type Orbitals "
397  WRITE (iw, '(A,F15.8,T41,A,F15.8)') " Initial exponent: ", atom_basis%aval, &
398  " Proportionality factor: ", atom_basis%cval
399  ELSE
400  WRITE (iw, '(T3,A)') "BASIS_TYPE GAUSSIAN"
401  END IF
402  DO l = 0, lmat
403  IF (atom_basis%nbas(l) > 0) THEN
404  SELECT CASE (l)
405  CASE DEFAULT
406  WRITE (iw, '(T3,A,(T15,F20.8,:," \"))') &
407  "X_EXPONENTS ", (atom_basis%am(i, l), i=1, atom_basis%nbas(l))
408  CASE (0)
409  WRITE (iw, '(T3,A,(T15,F20.8,:," \"))') &
410  "S_EXPONENTS ", (atom_basis%am(i, 0), i=1, atom_basis%nbas(0))
411  CASE (1)
412  WRITE (iw, '(T3,A,(T15,F20.8,:," \"))') &
413  "P_EXPONENTS ", (atom_basis%am(i, 1), i=1, atom_basis%nbas(1))
414  CASE (2)
415  WRITE (iw, '(T3,A,(T15,F20.8,:," \"))') &
416  "D_EXPONENTS ", (atom_basis%am(i, 2), i=1, atom_basis%nbas(2))
417  CASE (3)
418  WRITE (iw, '(T3,A,(T15,F20.8,:," \"))') &
419  "F_EXPONENTS ", (atom_basis%am(i, 3), i=1, atom_basis%nbas(3))
420  END SELECT
421  END IF
422  END DO
423  CASE (cgto_basis)
424  cpabort("")
425  CASE (sto_basis)
426  WRITE (iw, '(T3,A)') "BASIS_TYPE SLATER"
427  DO l = 0, lmat
428  IF (atom_basis%nbas(l) > 0) THEN
429  SELECT CASE (l)
430  CASE DEFAULT
431  WRITE (iw, '(T3,A,(T15,F20.8,:," \"))') &
432  "X_EXPONENTS ", (atom_basis%as(i, l), i=1, atom_basis%nbas(l))
433  WRITE (iw, '(T3,A,60I3)') &
434  "X_QUANTUM_NUMBERS ", (atom_basis%ns(i, l), i=1, atom_basis%nbas(l))
435  CASE (0)
436  WRITE (iw, '(T3,A,(T15,F20.8,:," \"))') &
437  "S_EXPONENTS ", (atom_basis%as(i, 0), i=1, atom_basis%nbas(0))
438  WRITE (iw, '(T3,A,60I3)') &
439  "S_QUANTUM_NUMBERS ", (atom_basis%ns(i, 0), i=1, atom_basis%nbas(0))
440  CASE (1)
441  WRITE (iw, '(T3,A,(T15,F20.8,:," \"))') &
442  "P_EXPONENTS ", (atom_basis%as(i, 1), i=1, atom_basis%nbas(1))
443  WRITE (iw, '(T3,A,60I3)') &
444  "P_QUANTUM_NUMBERS ", (atom_basis%ns(i, 1), i=1, atom_basis%nbas(1))
445  CASE (2)
446  WRITE (iw, '(T3,A,(T15,F20.8,:," \"))') &
447  "D_EXPONENTS ", (atom_basis%as(i, 2), i=1, atom_basis%nbas(2))
448  WRITE (iw, '(T3,A,60I3)') &
449  "D_QUANTUM_NUMBERS ", (atom_basis%ns(i, 2), i=1, atom_basis%nbas(2))
450  CASE (3)
451  WRITE (iw, '(T3,A,(T15,F20.8,:," \"))') &
452  "F_EXPONENTS ", (atom_basis%as(i, 3), i=1, atom_basis%nbas(3))
453  WRITE (iw, '(T3,A,60I3)') &
454  "F_QUANTUM_NUMBERS ", (atom_basis%ns(i, 3), i=1, atom_basis%nbas(3))
455  END SELECT
456  END IF
457  END DO
458  CASE (num_basis)
459  cpabort("")
460  CASE DEFAULT
461  cpabort("")
462  END SELECT
463 
464  IF (PRESENT(wfn)) THEN
465  SELECT CASE (atom_basis%basis_type)
466  CASE DEFAULT
467  CASE (gto_basis)
468  IF (.NOT. atom_basis%geometrical) THEN
469  WRITE (iw, '(/,T3,A)') "ORBITAL COEFFICENTS (Quickstep normalization)"
470  im = min(6, SIZE(wfn, 2))
471  DO l = 0, lmat
472  IF (atom_basis%nbas(l) > 0) THEN
473  WRITE (iw, '(T3,A,I3)') "L Quantum Number:", l
474  ! Quickstep normalization
475  expzet = 0.25_dp*real(2*l + 3, dp)
476  prefac = sqrt(rootpi/2._dp**(l + 2)*dfac(2*l + 1))
477  DO i = 1, atom_basis%nbas(l)
478  zeta = (2._dp*atom_basis%am(i, l))**expzet
479  WRITE (iw, '(T5,F14.8,2x,6F12.8)') atom_basis%am(i, l), wfn(i, 1:im, l)*prefac/zeta
480  END DO
481  END IF
482  END DO
483  END IF
484  END SELECT
485  END IF
486 
487  CALL close_file(unit_number=iw)
488 
489  END SUBROUTINE atom_print_basis_file
490 
491 ! **************************************************************************************************
492 !> \brief Print information about the electronic structure method in use.
493 !> \param atom information about the atomic kind
494 !> \param iw output file unit
495 !> \par History
496 !> * 09.2015 direct use of the LibXC Fortran interface [Andreas Gloess]
497 !> * 10.2012 LibXC interface [Fabien Tran]
498 !> * 02.2010 unrestricted KS and HF methods [Juerg Hutter]
499 !> * 04.2009 print geometrical Gaussian type orbitals [Juerg Hutter]
500 !> * 09.2008 new subroutine's prototype; print relativistic methods [Juerg Hutter]
501 !> * 09.2008 created [Juerg Hutter]
502 ! **************************************************************************************************
503  SUBROUTINE atom_print_method(atom, iw)
504  TYPE(atom_type) :: atom
505  INTEGER, INTENT(IN) :: iw
506 
507  CHARACTER(len=160) :: shortform
508  CHARACTER(LEN=20) :: tmpstr
509  CHARACTER(len=:), ALLOCATABLE :: reference
510  INTEGER :: ifun, il, meth, myfun, reltyp
511  LOGICAL :: lsd
512  TYPE(section_vals_type), POINTER :: xc_fun, xc_fun_section, xc_section
513 
514  NULLIFY (xc_fun, xc_fun_section, xc_section)
515 
516  meth = atom%method_type
517 
518  xc_section => atom%xc_section
519  xc_fun_section => section_vals_get_subs_vals(xc_section, "XC_FUNCTIONAL")
520  SELECT CASE (meth)
521  CASE DEFAULT
522  cpabort("")
523  CASE (do_rks_atom)
524  CALL section_vals_val_get(xc_fun_section, "_SECTION_PARAMETERS_", i_val=myfun)
525  CASE (do_uks_atom)
526  CALL section_vals_val_get(xc_fun_section, "_SECTION_PARAMETERS_", i_val=myfun)
527  CASE (do_rhf_atom)
528  myfun = xc_none
529  CASE (do_uhf_atom)
530  myfun = xc_none
531  CASE (do_rohf_atom)
532  myfun = xc_none
533  END SELECT
534 
535  SELECT CASE (meth)
536  CASE DEFAULT
537  cpabort("")
538  CASE (do_rks_atom)
539  IF (iw > 0) WRITE (iw, fmt="(/,' METHOD | Restricted Kohn-Sham Calculation')")
540  CASE (do_uks_atom)
541  IF (iw > 0) WRITE (iw, fmt="(/,' METHOD | Unrestricted Kohn-Sham Calculation')")
542  CASE (do_rhf_atom)
543  IF (iw > 0) WRITE (iw, fmt="(/,' METHOD | Restricted Hartree-Fock Calculation')")
544  CASE (do_uhf_atom)
545  IF (iw > 0) WRITE (iw, fmt="(/,' METHOD | Unrestricted Hartree-Fock Calculation')")
546  CASE (do_rohf_atom)
547  IF (iw > 0) WRITE (iw, fmt="(/,' METHOD | Restricted Open-Shell Kohn-Sham Calculation')")
548  END SELECT
549 
550  ! zmp
551  IF (atom%do_zmp) THEN
552  IF (iw > 0) WRITE (iw, fmt="(' ZMP | Method on atomic radial density')")
553  IF (iw > 0) WRITE (iw, fmt="(' ZMP | Lambda : ',F5.1)") atom%lambda
554  IF (iw > 0) WRITE (iw, fmt="(' ZMP | Reading external density : ',A20)") atom%ext_file
555  IF (atom%dm) THEN
556  IF (iw > 0) WRITE (iw, fmt="(' ZMP | The file is in the form of a density matrix')")
557  ELSE
558  IF (iw > 0) WRITE (iw, fmt="(' ZMP | The file is in the form of a linear density')")
559  END IF
560  IF (atom%doread) THEN
561  IF (iw > 0) WRITE (iw, fmt="(' ZMP | Restarting calculation from ',A20,' file if present')") atom%zmp_restart_file
562  END IF
563  ELSE IF (atom%read_vxc) THEN
564  IF (iw > 0) WRITE (iw, fmt="(' ZMP | Calculating density from external V_xc')")
565  IF (iw > 0) WRITE (iw, fmt="(' ZMP | Reading external v_xc file : ',A20)") atom%ext_vxc_file
566  END IF
567 
568  IF (atom%pp_calc) THEN
569  IF (iw > 0) WRITE (iw, fmt="(' METHOD | Nonrelativistic Calculation')")
570  ELSE
571  reltyp = atom%relativistic
572 
573  SELECT CASE (reltyp)
574  CASE DEFAULT
575  cpabort("")
576  CASE (do_nonrel_atom)
577  IF (iw > 0) WRITE (iw, fmt="(' METHOD | Nonrelativistic Calculation')")
578  CASE (do_zoramp_atom)
579  IF (iw > 0) WRITE (iw, fmt="(' METHOD | Relativistic Calculation using ZORA(MP)')")
580  CASE (do_sczoramp_atom)
581  IF (iw > 0) WRITE (iw, fmt="(' METHOD | Relativistic Calculation using scaled ZORA(MP)')")
582  CASE (do_dkh0_atom)
583  IF (iw > 0) WRITE (iw, fmt="(' METHOD | Relativistic Calculation using Douglas-Kroll 0th order')")
584  IF (iw > 0) WRITE (iw, fmt="(' METHOD | Relativistic Calculation using kietic energy scaling')")
585  CASE (do_dkh1_atom)
586  IF (iw > 0) WRITE (iw, fmt="(' METHOD | Relativistic Calculation using Douglas-Kroll 1st order')")
587  IF (iw > 0) WRITE (iw, fmt="(' METHOD | Relativistic Calculation using Foldy-Wouthuysen transformation')")
588  CASE (do_dkh2_atom)
589  IF (iw > 0) WRITE (iw, fmt="(' METHOD | Relativistic Calculation using Douglas-Kroll 2nd order')")
590  CASE (do_dkh3_atom)
591  IF (iw > 0) WRITE (iw, fmt="(' METHOD | Relativistic Calculation using Douglas-Kroll 3rd order')")
592  END SELECT
593  END IF
594 
595  lsd = (meth == do_uks_atom)
596 
597  IF (myfun /= xc_none) THEN
598  CALL section_vals_val_get(xc_section, "FUNCTIONAL_ROUTINE", c_val=tmpstr)
599  IF (iw > 0) WRITE (iw, fmt="(' FUNCTIONAL| ROUTINE=',a)") trim(tmpstr)
600  CALL xc_functionals_expand(xc_fun_section, xc_section)
601  IF (iw > 0) THEN
602  ifun = 0
603  DO
604  ifun = ifun + 1
605  xc_fun => section_vals_get_subs_vals2(xc_fun_section, i_section=ifun)
606  IF (.NOT. ASSOCIATED(xc_fun)) EXIT
607  IF (libxc_check_existence_in_libxc(xc_fun)) THEN
608  ALLOCATE (CHARACTER(LEN=libxc_get_reference_length(xc_fun, lsd)) :: reference)
609  ELSE
610  ALLOCATE (CHARACTER(LEN=20*default_string_length) :: reference)
611  END IF
612  CALL xc_functional_get_info(xc_fun, lsd=lsd, reference=reference, shortform=shortform)
613  WRITE (iw, fmt="(' FUNCTIONAL| ',a,':')") &
614  trim(xc_fun%section%name)
615  DO il = 1, len_trim(reference), 67
616  WRITE (iw, fmt="(' FUNCTIONAL| ',a67)") reference(il:)
617  END DO
618  DEALLOCATE (reference)
619  END DO
620  END IF
621  ELSE
622  IF (iw > 0) WRITE (iw, fmt="(' FUNCTIONAL| NO EXCHANGE-CORRELATION FUNCTIONAL USED.')")
623  END IF
624 
625  END SUBROUTINE atom_print_method
626 
627 ! **************************************************************************************************
628 !> \brief Print information about the pseudo-potential.
629 !> \param potential pseudo-potential
630 !> \param iw output file unit
631 !> \par History
632 !> * 05.2017 SGP pseudo-potentials [Juerg Hutter]
633 !> * 02.2016 pseudo-potential in Quantum Espresso UPF format [Juerg Hutter]
634 !> * 01.2016 new confinement potential form [Juerg Hutter]
635 !> * 03.2010 extension of GTH pseudo-potential definition [Juerg Hutter]
636 !> * 05.2009 GTH pseudo-potential [Juerg Hutter]
637 !> * 09.2008 created [Juerg Hutter]
638 ! **************************************************************************************************
639  SUBROUTINE atom_print_potential(potential, iw)
640  TYPE(atom_potential_type) :: potential
641  INTEGER, INTENT(IN) :: iw
642 
643  CHARACTER(len=60) :: pline
644  INTEGER :: i, j, k, l
645 
646  SELECT CASE (potential%ppot_type)
647  CASE (no_pseudo)
648  WRITE (iw, '(/," ",28("*"),A,27("*"))') " All Electron Potential "
649  CASE (gth_pseudo)
650  WRITE (iw, '(/," ",29("*"),A,29("*"))') " GTH Pseudopotential "
651  WRITE (iw, '(T10,A,T76,F5.1)') " Core Charge ", potential%gth_pot%zion
652  WRITE (iw, '(T10,A,T66,F15.6)') " Rc ", potential%gth_pot%rc
653  WRITE (pline, '(5F12.6)') (potential%gth_pot%cl(i), i=1, potential%gth_pot%ncl)
654  WRITE (iw, '(T10,A,T21,A60)') " C1 C2 ... ", adjustr(pline)
655  IF (potential%gth_pot%lpotextended) THEN
656  DO k = 1, potential%gth_pot%nexp_lpot
657  WRITE (iw, '(T10,A,F10.6,T38,A,4F10.6)') " LPot: rc=", potential%gth_pot%alpha_lpot(k), &
658  "CX=", (potential%gth_pot%cval_lpot(i, k), i=1, potential%gth_pot%nct_lpot(k))
659  END DO
660  END IF
661  IF (potential%gth_pot%nlcc) THEN
662  DO k = 1, potential%gth_pot%nexp_nlcc
663  WRITE (iw, '(T10,A,F10.6,T38,A,4F10.6)') " LSDPot: rc=", potential%gth_pot%alpha_nlcc(k), &
664  "CX=", (potential%gth_pot%cval_nlcc(i, k)*4.0_dp*pi, i=1, potential%gth_pot%nct_nlcc(k))
665  END DO
666  END IF
667  IF (potential%gth_pot%lsdpot) THEN
668  DO k = 1, potential%gth_pot%nexp_lsd
669  WRITE (iw, '(T10,A,F10.6,T38,A,4F10.6)') " LSDPot: rc=", potential%gth_pot%alpha_lsd(k), &
670  "CX=", (potential%gth_pot%cval_lsd(i, k), i=1, potential%gth_pot%nct_lsd(k))
671  END DO
672  END IF
673  DO l = 0, lmat
674  IF (potential%gth_pot%nl(l) > 0) THEN
675  WRITE (iw, '(T10,A,T76,I5)') " Angular momentum ", l
676  WRITE (iw, '(T10,A,T66,F15.6)') " Rcnl ", potential%gth_pot%rcnl(l)
677  WRITE (iw, '(T10,A,T76,I5)') " Nl ", potential%gth_pot%nl(l)
678  WRITE (pline, '(5F12.6)') (potential%gth_pot%hnl(1, j, l), j=1, potential%gth_pot%nl(l))
679  WRITE (iw, '(T10,A,T21,A60)') " Hnl ", adjustr(pline)
680  DO i = 2, potential%gth_pot%nl(l)
681  WRITE (pline, '(T21,5F12.6)') (potential%gth_pot%hnl(i, j, l), j=i, potential%gth_pot%nl(l))
682  WRITE (iw, '(T21,A60)') adjustr(pline)
683  END DO
684  END IF
685  END DO
686  CASE (upf_pseudo)
687  WRITE (iw, '(/," ",29("*"),A,29("*"))') " UPF Pseudopotential "
688  DO k = 1, potential%upf_pot%maxinfo
689  WRITE (iw, '(A80)') potential%upf_pot%info(k)
690  END DO
691  CASE (sgp_pseudo)
692  WRITE (iw, '(/," ",29("*"),A,29("*"))') " SGP Pseudopotential "
693  WRITE (iw, '(T10,A,T76,F5.1)') " Core Charge ", potential%sgp_pot%zion
694  CASE (ecp_pseudo)
695  WRITE (iw, '(/," ",26("*"),A,27("*"))') " Effective Core Potential "
696  WRITE (iw, '(T10,A,T76,F5.1)') " Core Charge ", potential%ecp_pot%zion
697  DO k = 1, potential%ecp_pot%nloc
698  IF (k == 1) THEN
699  WRITE (iw, '(T10,A,T40,I3,T49,2F16.8)') " Local Potential ", potential%ecp_pot%nrloc(k), &
700  potential%ecp_pot%bloc(k), potential%ecp_pot%aloc(k)
701  ELSE
702  WRITE (iw, '(T40,I3,T49,2F16.8)') potential%ecp_pot%nrloc(k), &
703  potential%ecp_pot%bloc(k), potential%ecp_pot%aloc(k)
704  END IF
705  END DO
706  DO l = 0, potential%ecp_pot%lmax
707  WRITE (iw, '(T10,A,I3)') " ECP l-value ", l
708  DO k = 1, potential%ecp_pot%npot(l)
709  WRITE (iw, '(T40,I3,T49,2F16.8)') potential%ecp_pot%nrpot(k, l), &
710  potential%ecp_pot%bpot(k, l), potential%ecp_pot%apot(k, l)
711  END DO
712  END DO
713  CASE DEFAULT
714  cpabort("")
715  END SELECT
716  IF (potential%confinement) THEN
717  IF (potential%conf_type == poly_conf) THEN
718  WRITE (iw, '(/,T10,A,T51,F12.6," * (R /",F6.2,")**",F6.2)') &
719  " Confinement Potential ", potential%acon, potential%rcon, potential%scon
720  ELSE IF (potential%conf_type == barrier_conf) THEN
721  WRITE (iw, '(/,T10,A)') " Confinement Potential s*F[(r-ron)/w] "
722  WRITE (iw, '(T57,A,F12.6,A)') "s =", potential%acon, " Ha"
723  WRITE (iw, '(T57,A,F12.6,A)') "w =", potential%rcon, " Bohr"
724  WRITE (iw, '(T57,A,F12.6,A)') "ron =", potential%scon, " Bohr"
725  ELSE
726  cpabort("")
727  END IF
728  ELSE
729  WRITE (iw, '(/,T10,A)') " No Confinement Potential is applied "
730  END IF
731  WRITE (iw, '(" ",79("*"))')
732 
733  END SUBROUTINE atom_print_potential
734 
735 ! **************************************************************************************************
736 !> \brief Print GTH pseudo-potential parameters.
737 !> \param gthpot pseudo-potential
738 !> \param iunit output file unit
739 !> \par History
740 !> * 09.2012 created [Juerg Hutter]
741 !> \note The pseudo-potential is written into the 'iunit' file unit or as the file 'GTH-PARAMETER'
742 !> inside the current working directory if the I/O unit is not given explicitly.
743 ! **************************************************************************************************
744  SUBROUTINE atom_write_pseudo_param(gthpot, iunit)
745  TYPE(atom_gthpot_type), INTENT(INOUT) :: gthpot
746  INTEGER, INTENT(IN), OPTIONAL :: iunit
747 
748  INTEGER :: i, iw, j, k, n
749 
750  IF (PRESENT(iunit)) THEN
751  iw = iunit
752  ELSE
753  CALL open_file(file_name="GTH-PARAMETER", file_status="UNKNOWN", file_action="WRITE", unit_number=iw)
754  END IF
755  WRITE (iw, '(A2,A)') gthpot%symbol, adjustl(trim(gthpot%pname))
756  WRITE (iw, '(4I5)') gthpot%econf(0:3)
757  WRITE (iw, '(F20.14,I8,5F20.14)') gthpot%rc, gthpot%ncl, (gthpot%cl(i), i=1, gthpot%ncl)
758  IF (gthpot%lpotextended) THEN
759  WRITE (iw, '(A,I5)') " LPOT", gthpot%nexp_lpot
760  DO i = 1, gthpot%nexp_lpot
761  WRITE (iw, '(F20.14,I8,5F20.14)') gthpot%alpha_lpot(i), gthpot%nct_lpot(i), &
762  (gthpot%cval_lpot(j, i), j=1, gthpot%nct_lpot(i))
763  END DO
764  END IF
765  IF (gthpot%lsdpot) THEN
766  WRITE (iw, '(A,I5)') " LSD ", gthpot%nexp_lsd
767  DO i = 1, gthpot%nexp_lsd
768  WRITE (iw, '(F20.14,I8,5F20.14)') gthpot%alpha_lsd(i), gthpot%nct_lsd(i), &
769  (gthpot%cval_lsd(j, i), j=1, gthpot%nct_lsd(i))
770  END DO
771  END IF
772  IF (gthpot%nlcc) THEN
773  WRITE (iw, '(A,I5)') " NLCC ", gthpot%nexp_nlcc
774  DO i = 1, gthpot%nexp_nlcc
775  WRITE (iw, '(F20.14,I8,5F20.14)') gthpot%alpha_nlcc(i), gthpot%nct_nlcc(i), &
776  (gthpot%cval_nlcc(j, i)*4.0_dp*pi, j=1, gthpot%nct_nlcc(i))
777  END DO
778  END IF
779  n = 0
780  DO i = lmat, 0, -1
781  IF (gthpot%nl(i) > 0) THEN
782  n = i + 1
783  EXIT
784  END IF
785  END DO
786  WRITE (iw, '(I8)') n
787  DO i = 0, n - 1
788  WRITE (iw, '(F20.14,I8,5F20.14)') gthpot%rcnl(i), gthpot%nl(i), (gthpot%hnl(1, k, i), k=1, gthpot%nl(i))
789  SELECT CASE (gthpot%nl(i))
790  CASE (2)
791  WRITE (iw, '(T49,F20.14)') gthpot%hnl(2, 2, i)
792  CASE (3)
793  WRITE (iw, '(T49,2F20.14)') gthpot%hnl(2, 2, i), gthpot%hnl(2, 3, i)
794  WRITE (iw, '(T69,F20.14)') gthpot%hnl(3, 3, i)
795  CASE DEFAULT
796  DO j = 2, gthpot%nl(i)
797  WRITE (iw, '(T29,5F20.14)') (gthpot%hnl(j, k, i), k=j, gthpot%nl(i))
798  END DO
799  END SELECT
800  END DO
801  IF (.NOT. PRESENT(iunit)) CALL close_file(unit_number=iw)
802 
803  END SUBROUTINE atom_write_pseudo_param
804 
805 ! **************************************************************************************************
806 !> \brief Print atomic orbitals.
807 !> \param atom information about the atomic kind
808 !> \param iw output file unit
809 !> \par History
810 !> * 04.2013 created [Juerg Hutter]
811 ! **************************************************************************************************
812  SUBROUTINE atom_print_orbitals(atom, iw)
813  TYPE(atom_type), POINTER :: atom
814  INTEGER, INTENT(IN) :: iw
815 
816  SELECT CASE (atom%method_type)
817  CASE DEFAULT
818  cpabort("")
819  CASE (do_rks_atom)
820  CALL atom_print_orbitals_helper(atom, atom%orbitals%wfn, "", iw)
821  CASE (do_uks_atom)
822  CALL atom_print_orbitals_helper(atom, atom%orbitals%wfna, "Alpha", iw)
823  CALL atom_print_orbitals_helper(atom, atom%orbitals%wfnb, "Beta", iw)
824  CASE (do_rhf_atom)
825  CALL atom_print_orbitals_helper(atom, atom%orbitals%wfn, "", iw)
826  CASE (do_uhf_atom)
827  CALL atom_print_orbitals_helper(atom, atom%orbitals%wfna, "Alpha", iw)
828  CALL atom_print_orbitals_helper(atom, atom%orbitals%wfnb, "Beta", iw)
829  CASE (do_rohf_atom)
830  cpabort("")
831  END SELECT
832 
833  END SUBROUTINE atom_print_orbitals
834 
835 ! **************************************************************************************************
836 !> \brief Print atomic orbitals of the given spin.
837 !> \param atom information about the atomic kind
838 !> \param wfn atomic orbitals
839 !> \param description description string
840 !> \param iw output file unit
841 !> \par History
842 !> * 04.2013 created [Juerg Hutter]
843 ! **************************************************************************************************
844  SUBROUTINE atom_print_orbitals_helper(atom, wfn, description, iw)
845  TYPE(atom_type), POINTER :: atom
846  REAL(kind=dp), DIMENSION(:, :, 0:), INTENT(INOUT) :: wfn
847  CHARACTER(len=*), INTENT(IN) :: description
848  INTEGER, INTENT(IN) :: iw
849 
850  INTEGER :: b, l, maxl, nb, nv, v
851 
852  WRITE (iw, '(/,A,A,A)') " Atomic orbital expansion coefficients [", description, "]"
853 
854  maxl = atom%state%maxl_calc
855  DO l = 0, maxl
856 
857  nb = atom%basis%nbas(l)
858  nv = atom%state%maxn_calc(l)
859  IF (nb > 0 .AND. nv > 0) THEN
860  nv = min(nv, SIZE(wfn, 2))
861  DO v = 1, nv
862  WRITE (iw, '(/," ORBITAL L = ",I1," State = ",I3)') l, v
863  DO b = 1, nb
864  WRITE (iw, '(" ",ES23.15)') wfn(b, v, l)
865  END DO
866  END DO
867  END IF
868  END DO
869  END SUBROUTINE atom_print_orbitals_helper
870 
871 END MODULE atom_output
Routines that print various information about an atomic kind.
Definition: atom_output.F:11
subroutine, public atom_print_energies(atom, iw)
Print energy components.
Definition: atom_output.F:175
subroutine, public atom_print_basis(atom_basis, iw, title)
Print atomic basis set.
Definition: atom_output.F:293
subroutine, public atom_print_iteration(iter, deps, etot, iw)
Print convergence information.
Definition: atom_output.F:272
subroutine, public atom_print_basis_file(atom_basis, wfn)
Print the optimized atomic basis set into a file.
Definition: atom_output.F:386
subroutine, public atom_write_pseudo_param(gthpot, iunit)
Print GTH pseudo-potential parameters.
Definition: atom_output.F:745
subroutine, public atom_print_state(state, iw)
Print information about electronic state.
Definition: atom_output.F:83
subroutine, public atom_print_zmp_iteration(iter, deps, atom, iw)
Printing of the atomic iterations when ZMP is active.
Definition: atom_output.F:248
subroutine, public atom_print_method(atom, iw)
Print information about the electronic structure method in use.
Definition: atom_output.F:504
subroutine, public atom_print_orbitals(atom, iw)
Print atomic orbitals.
Definition: atom_output.F:813
subroutine, public atom_print_potential(potential, iw)
Print information about the pseudo-potential.
Definition: atom_output.F:640
subroutine, public atom_print_info(zval, info, iw)
Print an information string related to the atomic kind.
Definition: atom_output.F:64
Define the atom type and its sub types.
Definition: atom_types.F:15
integer, parameter, public num_basis
Definition: atom_types.F:69
integer, parameter, public cgto_basis
Definition: atom_types.F:69
integer, parameter, public gto_basis
Definition: atom_types.F:69
integer, parameter, public sto_basis
Definition: atom_types.F:69
integer, parameter, public lmat
Definition: atom_types.F:67
Some basic routines for atomic calculations.
Definition: atom_utils.F:15
pure integer function, dimension(0:lmat), public get_maxn_occ(occupation)
Return the maximum principal quantum number of occupied orbitals.
Definition: atom_utils.F:302
subroutine, public get_rho0(atom, rho0)
Calculate the total electron density at R=0.
Definition: atom_utils.F:2219
pure integer function, public get_maxl_occ(occupation)
Return the maximum orbital quantum number of occupied orbitals.
Definition: atom_utils.F:282
Definition: atom.F:9
Utility routines to open and close files. Tracking of preconnections.
Definition: cp_files.F:16
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.
Definition: cp_files.F:308
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.
Definition: cp_files.F:119
collects all constants needed in input so that they can be used without circular dependencies
integer, parameter, public do_rhf_atom
integer, parameter, public do_rks_atom
integer, parameter, public do_dkh3_atom
integer, parameter, public do_nonrel_atom
integer, parameter, public do_dkh0_atom
integer, parameter, public do_uhf_atom
integer, parameter, public poly_conf
integer, parameter, public do_dkh2_atom
integer, parameter, public do_uks_atom
integer, parameter, public barrier_conf
integer, parameter, public do_zoramp_atom
integer, parameter, public do_dkh1_atom
integer, parameter, public do_rohf_atom
integer, parameter, public xc_none
integer, parameter, public do_sczoramp_atom
checks the input and perform some automatic "magic" on it
subroutine, public xc_functionals_expand(functionals, xc_section)
expand a shortcutted functional section
objects that represent the structure of input sections and the data contained in an input section
type(section_vals_type) function, pointer, public section_vals_get_subs_vals2(section_vals, i_section, i_rep_section)
returns the values of the n-th non default subsection (null if no such section exists (not so many no...
recursive type(section_vals_type) function, pointer, public section_vals_get_subs_vals(section_vals, subsection_name, i_rep_section, can_return_null)
returns the values of the requested subsection
subroutine, public section_vals_val_get(section_vals, keyword_name, i_rep_section, i_rep_val, n_rep_val, val, l_val, i_val, r_val, c_val, l_vals, i_vals, r_vals, c_vals, explicit)
returns the requested value
Defines the basic variable types.
Definition: kinds.F:23
integer, parameter, public dp
Definition: kinds.F:34
integer, parameter, public default_string_length
Definition: kinds.F:57
Definition of mathematical constants and functions.
Definition: mathconstants.F:16
real(kind=dp), parameter, public pi
real(kind=dp), dimension(-1:2 *maxfac+1), parameter, public dfac
Definition: mathconstants.F:61
real(kind=dp), parameter, public rootpi
Periodic Table related data definitions.
type(atom), dimension(0:nelem), public ptable
Definition of physical constants:
Definition: physcon.F:68
real(kind=dp), parameter, public evolt
Definition: physcon.F:183
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
Definition: xc_libxc.F:28
logical function, public libxc_check_existence_in_libxc(libxc_params)
This function checks whether a functional name belongs to LibXC.
Definition: xc_libxc.F:138
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.
Definition: xc_libxc.F:159