Octopus
ps.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-2012 M. Marques, A. Castro, A. Rubio, G. Bertsch, M. Oliveira
2!!
3!! This program is free software; you can redistribute it and/or modify
4!! it under the terms of the GNU General Public License as published by
5!! the Free Software Foundation; either version 2, or (at your option)
6!! any later version.
7!!
8!! This program is distributed in the hope that it will be useful,
9!! but WITHOUT ANY WARRANTY; without even the implied warranty of
10!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
11!! GNU General Public License for more details.
12!!
13!! You should have received a copy of the GNU General Public License
14!! along with this program; if not, write to the Free Software
15!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
16!! 02110-1301, USA.
17!!
18
19#include "global.h"
20
21module ps_oct_m
22 use atomic_oct_m
23 use debug_oct_m
24 use global_oct_m
25 use io_oct_m
26 use, intrinsic :: iso_fortran_env
30 use parser_oct_m
31 use logrid_oct_m
32 use math_oct_m
36 use ps_cpi_oct_m
37 use ps_fhi_oct_m
38 use ps_hgh_oct_m
39 use ps_xml_oct_m
41#ifdef HAVE_PSPIO
42 use fpspio_m
43#endif
44 use ps_psf_oct_m
45 use pseudo_oct_m
48 use string_oct_m
49
50 implicit none
51
52 private
53 public :: &
54 ps_t, &
55 ps_init, &
58 ps_filter, &
61 ps_debug, &
62 ps_niwfs, &
64 ps_end, &
70
71 integer, parameter, public :: &
72 PS_FILTER_NONE = 0, &
73 ps_filter_ts = 2, &
75
76 integer, public, parameter :: &
77 PROJ_NONE = 0, &
78 proj_hgh = 1, &
79 proj_kb = 2, &
80 proj_rkb = 3
81
82 integer, public, parameter :: &
83 PROJ_J_INDEPENDENT = 0, & !< Non-relativisitic or scalar-relativistic pseudopotentials
84 proj_j_dependent = 1, &
86
87 integer, parameter, public :: INVALID_L = 333
88
89 character(len=4), parameter :: ps_name(PSEUDO_FORMAT_UPF1:PSEUDO_FORMAT_PSP8) = &
90 (/"upf1", "upf2", "qso ", "psml", "psf ", "cpi ", "fhi ", "hgh ", "psp8"/)
91
93 type ps_t
94 ! Components are public by default
95 integer :: projector_type
96 integer :: relativistic_treatment
98
99 character(len=10), private :: label
100
101 integer, private :: ispin
102
103 real(real64), private :: z
104 real(real64) :: z_val
105 type(valconf_t) :: conf
106
107 type(logrid_t), private :: g
108
109 type(spline_t), allocatable :: ur(:, :)
110 logical, allocatable :: bound(:, :)
111
112 ! Kleinman-Bylander projectors stuff
113 integer :: lmax
114 integer :: llocal
115
116 type(spline_t) :: vl
117 logical :: no_vl = .false.
118
119 real(real64) :: projectors_sphere_threshold
126 real(real64) :: rc_max
127
128 integer :: kbc
129 integer :: projectors_per_l(5)
130 real(real64), allocatable :: h(:,:,:), k(:, :, :)
131 type(spline_t), allocatable :: kb(:, :)
132 type(spline_t), allocatable :: dkb(:, :)
133
134 logical :: nlcc
135 type(spline_t) :: core
136 type(spline_t) :: core_der
137
138
139 !LONG-RANGE PART OF THE LOCAL POTENTIAL
140
141 logical, private :: has_long_range
142 logical, private :: is_separated
143
144 type(spline_t), private :: vlr
145 type(spline_t) :: nlr
146
147 real(real64) :: sigma_erf
148
149 logical, private :: has_density
150 type(spline_t), allocatable :: density(:)
151 type(spline_t), allocatable :: density_der(:)
152
153 logical :: local
154 integer :: file_format
155 integer, private :: pseudo_type
156 integer :: exchange_functional
157 integer :: correlation_functional
158 end type ps_t
159
160 real(real64), parameter :: eps = 1.0e-8_real64
161
162contains
163
164
165 ! ---------------------------------------------------------
166 subroutine ps_init(ps, namespace, label, z, user_lmax, user_llocal, ispin, filename)
167 type(ps_t), intent(out) :: ps
168 type(namespace_t), intent(in) :: namespace
169 character(len=10), intent(in) :: label
170 integer, intent(in) :: user_lmax
171 integer, intent(in) :: user_llocal
172 integer, intent(in) :: ispin
173 real(real64), intent(in) :: z
174 character(len=*), intent(in) :: filename
175
176 integer :: l, ii, ll, is, ierr
177 type(ps_psf_t) :: ps_psf
178 type(ps_cpi_t) :: ps_cpi
179 type(ps_fhi_t) :: ps_fhi
180 type(ps_hgh_t) :: ps_hgh
181 type(ps_xml_t) :: ps_xml
182 real(real64), allocatable :: eigen(:, :)
183
184 push_sub(ps_init)
185
186 ps%exchange_functional = pseudo_exchange_unknown
187 ps%correlation_functional = pseudo_correlation_unknown
189 ! Fix the threshold to calculate the radius of the projector-function localization spheres:
191 call messages_obsolete_variable(namespace, 'SpecieProjectorSphereThreshold', 'SpeciesProjectorSphereThreshold')
192
193 !%Variable SpeciesProjectorSphereThreshold
194 !%Type float
195 !%Default 0.001
196 !%Section System::Species
197 !%Description
198 !% The pseudopotentials may be composed of a local part, and a linear combination of nonlocal
199 !% operators. These nonlocal projectors have "projector" form, <math> \left| v \right> \left< v \right| </math>
200 !% (or, more generally speaking, <math> \left| u \right> \left< v \right| </math>).
201 !% These projectors are localized in real space -- that is, the function <math>v</math>
202 !% has a finite support around the nucleus. This region where the projectors are localized should
203 !% be small or else the computation time required to operate with them will be very large.
204 !%
205 !% In practice, this localization is fixed by requiring the definition of the projectors to be
206 !% contained in a sphere of a certain radius. This radius is computed by making sure that the
207 !% absolute value of the projector functions, at points outside the localization sphere, is
208 !% below a certain threshold. This threshold is set by <tt>SpeciesProjectorSphereThreshold</tt>.
209 !%End
210 call parse_variable(namespace, 'SpeciesProjectorSphereThreshold', 0.001_real64, ps%projectors_sphere_threshold)
211 if (ps%projectors_sphere_threshold <= m_zero) call messages_input_error(namespace, 'SpeciesProjectorSphereThreshold')
213 ps%file_format = pseudo_detect_format(string_f_to_c(filename))
215 if (ps%file_format == pseudo_format_file_not_found) then
216 call messages_write("Cannot open pseudopotential file '"//trim(filename)//"'.")
217 call messages_fatal(namespace=namespace)
218 end if
219
220 if (ps%file_format == pseudo_format_unknown) then
221 call messages_write("Cannot determine the pseudopotential type for species '"//trim(label)//"' from", &
222 new_line = .true.)
223 call messages_write("file '"//trim(filename)//"'.")
224 call messages_fatal(namespace=namespace)
225 end if
227 ps%label = label
228 ps%ispin = ispin
229 ps%relativistic_treatment = proj_j_independent
230 ps%projector_type = proj_kb
231 ps%sigma_erf = 0.625_real64 ! This is hard-coded to a reasonable value
232
233 ps%projectors_per_l = 0
234
235 select case (ps%file_format)
237 ps%has_density = .true.
238 case default
239 ps%has_density = .false.
240 end select
241
242 select case (ps%file_format)
243 case (pseudo_format_psf)
244 ps%pseudo_type = pseudo_type_semilocal
246 call ps_psf_init(ps_psf, ispin, filename, namespace)
247
248 call valconf_copy(ps%conf, ps_psf%conf)
249 ps%z = z
250 ps%conf%z = nint(z) ! atomic number
251 ps%kbc = 1 ! only one projector per angular momentum
253 ps%lmax = ps_psf%ps_grid%no_l_channels - 1
254
255 if (user_lmax /= invalid_l) then
256 ps%lmax = min(ps%lmax, user_lmax) ! Maybe the file does not have enough components.
257 if (user_lmax /= ps%lmax) then
258 message(1) = "lmax in Species block for " // trim(label) // &
259 " is larger than number available in pseudopotential."
260 call messages_fatal(1, namespace=namespace)
261 end if
262 end if
263
264 ps%conf%p = ps_psf%ps_grid%no_l_channels
265 if (ps%lmax == 0) ps%llocal = 0 ! Vanderbilt is not acceptable if ps%lmax == 0.
266
267 ! the local part of the pseudo
268 if (user_llocal == invalid_l) then
269 ps%llocal = 0
270 else
271 ps%llocal = user_llocal
272 end if
273
274 ps%projectors_per_l(1:ps%lmax+1) = 1
275 if (ps%llocal > -1) ps%projectors_per_l(ps%llocal+1) = 0
276
277 call ps_psf_process(ps_psf, namespace, ps%lmax, ps%llocal)
278 call logrid_copy(ps_psf%ps_grid%g, ps%g)
279
281 ps%pseudo_type = pseudo_type_semilocal
282
283 if (ps%file_format == pseudo_format_cpi) then
284 call ps_cpi_init(ps_cpi, trim(filename), namespace)
285 ps%conf%p = ps_cpi%ps_grid%no_l_channels
286 else
287 call ps_fhi_init(ps_fhi, trim(filename), namespace)
288 ps%conf%p = ps_fhi%ps_grid%no_l_channels
289 end if
290
291 ps%conf%z = nint(z)
292 ps%conf%symbol = label(1:2)
293 ps%conf%type = 1
294 do l = 1, ps%conf%p
295 ps%conf%l(l) = l - 1
296 end do
297
298 ps%z = z
299 ps%kbc = 1 ! only one projector per angular momentum
300
301 ps%lmax = ps%conf%p - 1
302
303 if (user_lmax /= invalid_l) then
304 ps%lmax = min(ps%lmax, user_lmax) ! Maybe the file does not have enough components.
305 if (user_lmax /= ps%lmax) then
306 message(1) = "lmax in Species block for " // trim(label) // &
307 " is larger than number available in pseudopotential."
308 call messages_fatal(1, namespace=namespace)
309 end if
310 end if
311
312 if (ps%lmax == 0) ps%llocal = 0 ! Vanderbilt is not acceptable if ps%lmax == 0.
313
314 ! the local part of the pseudo
315 if (user_llocal == invalid_l) then
316 ps%llocal = 0
317 else
318 ps%llocal = user_llocal
319 end if
320
321 ps%projectors_per_l(1:ps%lmax+1) = 1
322 if (ps%llocal > -1) ps%projectors_per_l(ps%llocal+1) = 0
323
324 if (ps%file_format == pseudo_format_cpi) then
325 call ps_cpi_process(ps_cpi, ps%llocal, namespace)
326 call logrid_copy(ps_cpi%ps_grid%g, ps%g)
327 else
328 call ps_fhi_process(ps_fhi, ps%lmax, ps%llocal, namespace)
329 call logrid_copy(ps_fhi%ps_grid%g, ps%g)
330 end if
331
332 case (pseudo_format_hgh)
333 ps%pseudo_type = pseudo_type_kleinman_bylander
334 ps%projector_type = proj_hgh
335
336 call hgh_init(ps_hgh, trim(filename), namespace)
337 call valconf_copy(ps%conf, ps_hgh%conf)
338
339 ps%z = z
340 ps%z_val = ps_hgh%z_val
341 ps%kbc = 3
342 ps%llocal = -1
343 ps%lmax = ps_hgh%l_max
344 ps%conf%symbol = label(1:2)
345 ps%sigma_erf = ps_hgh%rlocal ! We use the correct value
346
347 ps%projectors_per_l(1:ps%lmax+1) = 1
348
349 ! Get the occupations from the valence charge of the atom
350 ps%conf%occ = m_zero
351 ! We impose here non-spin-polarized occupations, to preserve the behavior of the code
352 ! We might want to change this to get a better LCAO guess
353 call ps_guess_atomic_occupations(namespace, ps%z, ps%z_val, 1, ps%conf)
354 ! We need the information to solve the Schrodinder equation
355 call valconf_copy(ps_hgh%conf, ps%conf)
356
357 call hgh_process(ps_hgh, namespace)
358 call logrid_copy(ps_hgh%g, ps%g)
359
360 ! In case of spin-polarized calculations, we properly distribute the electrons
361 if (ispin == 2) then
363 end if
364
365 case (pseudo_format_qso, pseudo_format_upf1, pseudo_format_upf2, pseudo_format_psml, pseudo_format_psp8)
366
367 call ps_xml_init(ps_xml, namespace, trim(filename), ps%file_format, ierr)
368
369 ps%pseudo_type = pseudo_type(ps_xml%pseudo)
370 ps%exchange_functional = pseudo_exchange(ps_xml%pseudo)
371 ps%correlation_functional = pseudo_correlation(ps_xml%pseudo)
372
373 ps%z = z
374 ps%conf%z = nint(z)
375
376 if (ps_xml%kleinman_bylander) then
377 ps%conf%p = ps_xml%nwavefunctions
378 else
379 ps%conf%p = ps_xml%lmax + 1
380 end if
381
382 do ll = 0, ps_xml%lmax
383 ps%conf%l(ll + 1) = ll
384 end do
385
386 ps%kbc = ps_xml%nchannels
387 ps%lmax = ps_xml%lmax
388
389 ps%projectors_per_l = 0
390 do ll = 0, ps_xml%lmax
391 ps%projectors_per_l(ll+1) = pseudo_nprojectors_per_l(ps_xml%pseudo, ll)
392 end do
393
394 if (ps_xml%kleinman_bylander) then
395 ps%llocal = ps_xml%llocal
396 else
397 ! we have several options
398 ps%llocal = 0 ! the default
399 if (ps_xml%llocal >= 0) ps%llocal = ps_xml%llocal ! the one given in the pseudopotential file
400 if (user_llocal /= invalid_l) ps%llocal = user_llocal ! user supplied local component
401 assert(ps%llocal >= 0)
402 assert(ps%llocal <= ps%lmax)
403 end if
404
405 ps%g%nrval = ps_xml%grid_size
406
407 safe_allocate(ps%g%rofi(1:ps%g%nrval))
408 safe_allocate(ps%g%r2ofi(1:ps%g%nrval))
409
410 do ii = 1, ps%g%nrval
411 ps%g%rofi(ii) = ps_xml%grid(ii)
412 ps%g%r2ofi(ii) = ps_xml%grid(ii)**2
413 end do
414
415 end select
416
417 ps%local = (ps%lmax == 0 .and. ps%llocal == 0 ) .or. (ps%lmax == -1 .and. ps%llocal == -1)
418
419 ! We allocate all the stuff
420 safe_allocate(ps%kb (0:ps%lmax, 1:ps%kbc))
421 safe_allocate(ps%dkb (0:ps%lmax, 1:ps%kbc))
422 safe_allocate(ps%ur (1:ps%conf%p, 1:ps%ispin))
423 safe_allocate(ps%bound(1:ps%conf%p, 1:ps%ispin))
424 safe_allocate(ps%h (0:ps%lmax, 1:ps%kbc, 1:ps%kbc))
425 safe_allocate(ps%density(1:ps%ispin))
426 safe_allocate(ps%density_der(1:ps%ispin))
427
428 call spline_init(ps%kb)
429 call spline_init(ps%dkb)
430 call spline_init(ps%vl)
431 call spline_init(ps%core)
432 call spline_init(ps%core_der)
433 call spline_init(ps%density)
434 call spline_init(ps%density_der)
435
436 safe_allocate(eigen(1:ps%conf%p, 1:ps%ispin))
437 eigen = m_zero
438
439 ! Now we load the necessary information.
440 select case (ps%file_format)
441 case (pseudo_format_psf)
442 call ps_psf_get_eigen(ps_psf, eigen)
443 call ps_grid_load(ps, ps_psf%ps_grid)
444 call ps_psf_end(ps_psf)
445 case (pseudo_format_cpi)
446 call ps_grid_load(ps, ps_cpi%ps_grid)
447 call ps_cpi_end(ps_cpi)
448 case (pseudo_format_fhi)
449 call ps_grid_load(ps, ps_fhi%ps_grid)
450 call ps_fhi_end(ps_fhi)
451 case (pseudo_format_hgh)
452 call hgh_get_eigen(ps_hgh, eigen)
453 safe_allocate(ps%k (0:ps%lmax, 1:ps%kbc, 1:ps%kbc))
454 call hgh_load(ps, ps_hgh)
455 call hgh_end(ps_hgh)
456 case (pseudo_format_qso, pseudo_format_upf1, pseudo_format_upf2, pseudo_format_psml, pseudo_format_psp8)
457 call ps_xml_load(ps, ps_xml, namespace)
458 call ps_xml_end(ps_xml)
459 end select
460
461 if (ps_has_density(ps)) then
462 do is = 1, ps%ispin
463 call spline_der(ps%density(is), ps%density_der(is), ps%projectors_sphere_threshold)
464 end do
465 end if
466
467 if (ps_has_nlcc(ps)) then
468 call spline_der(ps%core, ps%core_der, ps%projectors_sphere_threshold)
469 end if
470
471 call ps_check_bound(ps, eigen)
472
473 ps%has_long_range = .true.
474 ps%is_separated = .false.
475
476 call ps_info(ps, filename, namespace)
477
478 safe_deallocate_a(eigen)
479
480 pop_sub(ps_init)
481 end subroutine ps_init
482
483 !------------------------------------------------------------------------
484
485 subroutine ps_info(ps, filename, namespace)
486 type(ps_t), intent(in) :: ps
487 character(len=*), intent(in) :: filename
488 type(namespace_t), intent(in) :: namespace
489
490 call messages_write(" Species '"//trim(ps%label)//"'", new_line = .true.)
491 call messages_write(" type : pseudopotential", new_line = .true.)
492 call messages_write(" file : '"//trim(filename)//"'")
493 call messages_info(namespace=namespace)
494
495 call messages_write(" file format :")
496 select case (ps%file_format)
497 case (pseudo_format_upf1)
498 call messages_write(" UPF1")
499 case (pseudo_format_upf2)
500 call messages_write(" UPF2")
501 case (pseudo_format_qso)
502 call messages_write(" QSO")
503 case (pseudo_format_psml)
504 call messages_write(" PSML")
505 case (pseudo_format_psp8)
506 call messages_write(" PSP8")
507 case (pseudo_format_psf)
508 call messages_write(" PSF")
509 case (pseudo_format_cpi)
510 call messages_write(" CPI")
511 case (pseudo_format_fhi)
512 call messages_write(" FHI")
513 case (pseudo_format_hgh)
514 call messages_write(" HGH")
515 end select
516 call messages_new_line()
517
518 call messages_write(" valence charge :")
519 call messages_write(ps%z_val, align_left = .true., fmt = '(f4.1)')
520 call messages_info(namespace=namespace)
521
522 call messages_write(" atomic number :")
523 call messages_write(nint(ps%z), fmt = '(i4)')
524 call messages_info(namespace=namespace)
525
526 call messages_write(" form on file :")
527 select case (ps%pseudo_type)
529 call messages_write(" ultrasoft")
531 call messages_write(" semilocal")
533 call messages_write(" kleinman-bylander")
534 case (pseudo_type_paw)
535 call messages_write(" paw")
536 end select
537 call messages_info(namespace=namespace)
538
539 if (ps%pseudo_type == pseudo_type_semilocal) then
540 call messages_write(" orbital origin :")
541 select case (ps%file_format)
542 case (pseudo_format_psf)
543 call messages_write(" calculated")
544 case default
545 call messages_write(" from file")
546 end select
547 call messages_info(namespace=namespace)
548 end if
549
550 call messages_write(" lmax :")
551 call messages_write(ps%lmax, fmt = '(i2)')
552 call messages_info(namespace=namespace)
553
554 call messages_write(" llocal :")
555 if (ps%llocal >= 0) then
556 call messages_write(ps%llocal, fmt = '(i2)')
557 else
558 call messages_write(ps%llocal, fmt = '(i3)')
559 end if
560 call messages_info(namespace=namespace)
561
562 call messages_write(" projectors per l :")
563 call messages_write(ps%kbc, fmt = '(i2)')
564 call messages_info(namespace=namespace)
565
566 call messages_write(" total projectors :")
567 if (ps%llocal < 0) then
568 call messages_write(ps%kbc*(ps%lmax + 1), fmt = '(i2)')
569 else
570 call messages_write(ps%kbc*ps%lmax, fmt = '(i2)')
571 end if
572 call messages_info(namespace=namespace)
573
574 if (ps%local) then
575 call messages_write(" application form : local")
576 else
577 call messages_write(" application form : kleinman-bylander")
578 end if
579 call messages_info(namespace=namespace)
581 call messages_write(" orbitals :")
582 call messages_write(ps_niwfs(ps), fmt='(i3)')
583 call messages_info(namespace=namespace)
584 call messages_write(" bound orbitals :")
585 call messages_write(ps_bound_niwfs(ps), fmt='(i3)')
586 call messages_info(namespace=namespace)
587
588 call messages_info(namespace=namespace)
589
590 end subroutine ps_info
591
592
593 ! ---------------------------------------------------------
596 subroutine ps_separate(ps)
597 type(ps_t), intent(inout) :: ps
598
599 real(real64), allocatable :: vsr(:), vlr(:), nlr(:)
600 real(real64) :: r, exp_arg, max_val_vsr
601 integer :: ii
602
603 push_sub(ps_separate)
604
605 assert(ps%g%nrval > 0)
606
607 safe_allocate(vsr(1:ps%g%nrval))
608 safe_allocate(vlr(1:ps%g%nrval))
609 safe_allocate(nlr(1:ps%g%nrval))
610
611 ps%no_vl = .false.
612
613 max_val_vsr = m_zero
614
615 do ii = 1, ps%g%nrval
616 r = ps%g%rofi(ii)
617
618 vlr(ii) = long_range_potential(r, ps%sigma_erf, ps%z_val)
619
620 vsr(ii) = spline_eval(ps%vl, r) - vlr(ii)
621 if(abs(vsr(ii)) < 1.0e-14_real64) vsr(ii) = m_zero
622 max_val_vsr = max(max_val_vsr, abs(vsr(ii)))
623
624 exp_arg = -m_half*r**2/ps%sigma_erf**2
625 if (exp_arg > m_min_exp_arg) then
626 nlr(ii) = -ps%z_val/(ps%sigma_erf*sqrt(m_two*m_pi))**3*exp(exp_arg)
627 else
628 nlr(ii) = m_zero
629 end if
630 end do
631
632 call spline_init(ps%vlr)
633 call spline_fit(ps%g%nrval, ps%g%rofi, vlr, ps%vlr, ps%projectors_sphere_threshold)
634
635 call spline_init(ps%nlr)
636 call spline_fit(ps%g%nrval, ps%g%rofi, nlr, ps%nlr, ps%projectors_sphere_threshold)
637
638 !overwrite vl
639 call spline_end(ps%vl)
640 call spline_init(ps%vl)
641 call spline_fit(ps%g%nrval, ps%g%rofi, vsr, ps%vl, ps%projectors_sphere_threshold)
642 if (max_val_vsr < 1.0e-12_real64) ps%no_vl = .true.
643
644 safe_deallocate_a(vsr)
645 safe_deallocate_a(vlr)
646 safe_deallocate_a(nlr)
647
648 ps%is_separated = .true.
649
650 pop_sub(ps_separate)
651 end subroutine ps_separate
652
654 real(real64) pure function long_range_potential(r, sigma, z_val) result(vlr)
655 real(real64), intent(in) :: r, sigma, z_val
656
657 real(real64) :: arg
658
659 arg = r/(sigma*sqrt(m_two))
660 ! Note that the threshold is taken such that the Taylor expansion produces an error smaller than \epsilon
661 ! given a third order Taylor expansion
662 if(arg > 1e-3_real64) then
663 vlr = -z_val * erf(arg) / r
664 else
665 vlr = -z_val * m_two / (sigma*sqrt(m_two*m_pi)) * (m_one - arg**2 / m_three)
666 end if
667 end function long_range_potential
668
669 ! ---------------------------------------------------------
670 subroutine ps_getradius(ps)
671 type(ps_t), intent(inout) :: ps
672 integer :: l, j
673
674 push_sub(ps_getradius)
675
676 ps%rc_max = m_zero
677
678 do l = 0, ps%lmax
679 if (l == ps%llocal) cycle
680 do j = 1, ps%kbc
681 ps%rc_max = max(ps%rc_max, ps%kb(l, j)%x_threshold)
682 end do
683 end do
684
685 pop_sub(ps_getradius)
686 end subroutine ps_getradius
687
688
689 ! ---------------------------------------------------------
690 subroutine ps_derivatives(ps)
691 type(ps_t), intent(inout) :: ps
692 integer :: l, j
693
694 push_sub(ps_derivatives)
695
696 do l = 0, ps%lmax
697 do j = 1, ps%kbc
698 call spline_der(ps%kb(l, j), ps%dkb(l, j), ps%projectors_sphere_threshold)
699 end do
700 end do
701
702 pop_sub(ps_derivatives)
703 end subroutine ps_derivatives
704
705
706 ! ---------------------------------------------------------
707 subroutine ps_filter(ps, filter, gmax)
708 type(ps_t), intent(inout) :: ps
709 integer, intent(in) :: filter
710 real(real64), intent(in) :: gmax
711
712 integer :: l, k, ispin
713
714 real(real64) :: alpha, beta_fs, rmax, rcut, gamma, beta_rs
715
716 push_sub(ps_filter)
717 call profiling_in("PS_FILTER")
718
719 select case (filter)
720 case (ps_filter_none)
721
722 case (ps_filter_ts)
723 alpha = 1.1_real64
724 gamma = m_two
725
726 if(.not. ps%no_vl) then
727 rmax = ps%vl%x_threshold
728 ! For the special case l=-1, we use l=0 when doing the filtering, but it is not clear this is correct
729 call spline_filter_mask(ps%vl, max(0, ps%llocal), rmax, gmax, alpha, gamma, ps%projectors_sphere_threshold)
730 end if
731
732 do l = 0, ps%lmax
733 if (l == ps%llocal) cycle
734 do k = 1, ps%kbc
735 call spline_filter_mask(ps%kb(l, k), l, ps%rc_max, gmax, alpha, gamma, ps%projectors_sphere_threshold)
736 end do
737 end do
738
739 if (ps_has_nlcc(ps)) then
740 rmax = ps%core%x_threshold
741 call spline_filter_mask(ps%core, 0, rmax, gmax, alpha, gamma, ps%projectors_sphere_threshold)
742 end if
743
744 if (ps_has_density(ps)) then
745 do ispin = 1, ps%ispin
746 if (abs(spline_integral(ps%density(ispin))) > 1.0e-12_real64) then
747 rmax = ps%density(ispin)%x_threshold
748 call spline_filter_mask(ps%density(ispin), 0, rmax, gmax, alpha, gamma, ps%projectors_sphere_threshold)
749 call spline_force_pos(ps%density(ispin), ps%projectors_sphere_threshold)
750 end if
751
752 if (abs(spline_integral(ps%density_der(ispin))) > 1.0e-12_real64) then
753 rmax = ps%density_der(ispin)%x_threshold
754 call spline_filter_mask(ps%density_der(ispin), 0, rmax, gmax, alpha, gamma, ps%projectors_sphere_threshold)
755 end if
756 end do
757 end if
758
759 case (ps_filter_bsb)
760 alpha = 0.7_real64 ! The original was M_FOUR/7.0_real64
761 beta_fs = 18.0_real64
762 rcut = 2.5_real64
763 beta_rs = 0.4_real64
764
765 ! For the special case l=-1, we use l=0 when doing the filtering, but it is not clear this is correct
766 call spline_filter_bessel(ps%vl, max(0, ps%llocal), gmax, alpha, beta_fs, rcut, beta_rs, ps%projectors_sphere_threshold)
767 do l = 0, ps%lmax
768 if (l == ps%llocal) cycle
769 do k = 1, ps%kbc
770 call spline_filter_bessel(ps%kb(l, k), l, gmax, alpha, beta_fs, rcut, beta_rs, ps%projectors_sphere_threshold)
771 end do
772 end do
773
774 if (ps_has_nlcc(ps)) then
775 call spline_filter_bessel(ps%core, 0, gmax, alpha, beta_fs, rcut, beta_rs, ps%projectors_sphere_threshold)
776 end if
777
778 if (ps_has_density(ps)) then
779 do ispin = 1, ps%ispin
780 call spline_filter_bessel(ps%density(ispin), 0, gmax, alpha, beta_fs, rcut, beta_rs, ps%projectors_sphere_threshold)
781 call spline_force_pos(ps%density(ispin), ps%projectors_sphere_threshold)
782 call spline_filter_bessel(ps%density_der(ispin), 0, gmax, alpha, beta_fs, rcut, beta_rs, ps%projectors_sphere_threshold)
783 end do
784 end if
786 end select
787
788 call profiling_out("PS_FILTER")
789 pop_sub(ps_filter)
790 end subroutine ps_filter
791
792 ! ---------------------------------------------------------
793 subroutine ps_check_bound(ps, eigen)
794 type(ps_t), intent(inout) :: ps
795 real(real64), intent(in) :: eigen(:,:)
796
797 integer :: i, is, ir
798 real(real64) :: ur1, ur2
799
800 push_sub(ps_check_bound)
801
802 ! Unbound states have positive eigenvalues
803 where(eigen > m_zero)
804 ps%bound = .false.
805 elsewhere
806 ps%bound = .true.
807 end where
808
809 ! We might not have information about the eigenvalues, so we need to check the wavefunctions
810 do i = 1, ps%conf%p
811 do is = 1, ps%ispin
812 if (.not. ps%bound(i, is)) cycle
813
814 do ir = ps%g%nrval, 3, -1
815 ! First we look for the outmost value that is not zero
816 if (abs(spline_eval(ps%ur(i, is), ps%g%rofi(ir))*ps%g%rofi(ir)) > m_zero) then
817 ! Usually bound states have exponentially decaying wavefunctions,
818 ! while unbound states have exponentially diverging
819 ! wavefunctions. Therefore we check if the wavefunctions
820 ! value is increasing with increasing radius. The fact
821 ! that we do not use the wavefunctions outmost value that
822 ! is not zero is on purpose, as some pseudopotential
823 ! generators do funny things with that point.
824 ur1 = spline_eval(ps%ur(i, is), ps%g%rofi(ir-2))*ps%g%rofi(ir-2)
825 ur2 = spline_eval(ps%ur(i, is), ps%g%rofi(ir-1))*ps%g%rofi(ir-1)
826 if ((ur1*ur2 > m_zero) .and. (abs(ur2) > abs(ur1))) ps%bound(i, is) = .false.
827 exit
828 end if
829 end do
830 end do
831 end do
832
833 pop_sub(ps_check_bound)
834 end subroutine ps_check_bound
835
836
837 ! ---------------------------------------------------------
838 subroutine ps_debug(ps, dir, namespace, gmax)
839 type(ps_t), intent(in) :: ps
840 character(len=*), intent(in) :: dir
841 type(namespace_t), intent(in) :: namespace
842 real(real64), intent(in) :: gmax
843
844 ! We will plot also some Fourier transforms.
845 type(spline_t), allocatable :: fw(:, :)
846
847 integer :: iunit
848 integer :: j, k, l
849
850 push_sub(ps_debug)
851
852 ! A text file with some basic data.
853 iunit = io_open(trim(dir)//'/pseudo-info', namespace, action='write')
854 write(iunit,'(a,/)') ps%label
855 write(iunit,'(a,a,/)') 'Format : ', ps_name(ps%file_format)
856 write(iunit,'(a,f6.3)') 'z : ', ps%z
857 write(iunit,'(a,f6.3,/)') 'zval : ', ps%z_val
858 write(iunit,'(a,i4)') 'lmax : ', ps%lmax
859 write(iunit,'(a,i4)') 'lloc : ', ps%llocal
860 write(iunit,'(a,i4,/)') 'kbc : ', ps%kbc
861 write(iunit,'(a,f9.5,/)') 'rcmax : ', ps%rc_max
862 write(iunit,'(a,/)') 'h matrix:'
863 do l = 0, ps%lmax
864 do k = 1, ps%kbc
865 write(iunit,'(10f9.5)') (ps%h(l, k, j), j = 1, ps%kbc)
866 end do
867 end do
868 if (allocated(ps%k)) then
869 write(iunit,'(/,a,/)') 'k matrix:'
870 do l = 0, ps%lmax
871 do k = 1, ps%kbc
872 write(iunit,'(10f9.5)') (ps%k(l, k, j), j = 1, ps%kbc)
873 end do
874 end do
875 end if
876
877 write(iunit,'(/,a)') 'orbitals:'
878 do j = 1, ps%conf%p
879 if (ps%ispin == 2) then
880 write(iunit,'(1x,a,i2,3x,a,i2,3x,a,f5.1,3x,a,l1,3x,a,f5.1,3x,a,f5.1)') 'n = ', ps%conf%n(j), 'l = ', ps%conf%l(j), &
881 'j = ', ps%conf%j(j), 'bound = ', all(ps%bound(j,:)), &
882 'occ(1) = ', ps%conf%occ(j, 1), 'occ(2) = ', ps%conf%occ(j, 2)
883 else
884 write(iunit,'(1x,a,i2,3x,a,i2,3x,a,f5.1,3x,a,l1,3x,a,f5.1)') 'n = ', ps%conf%n(j), 'l = ', ps%conf%l(j), &
885 'j = ', ps%conf%j(j), 'bound = ', all(ps%bound(j,:)), &
886 'occ = ', ps%conf%occ(j, 1)
887 end if
888 end do
889
890
891 call io_close(iunit)
892
893 ! Local part of the pseudopotential
894 iunit = io_open(trim(dir)//'/local', namespace, action='write')
895 call spline_print(ps%vl, iunit)
896 call io_close(iunit)
897
898 ! Local part of the pseudopotential
899 iunit = io_open(trim(dir)//'/local_long_range', namespace, action='write')
900 call spline_print(ps%vlr, iunit)
901 call io_close(iunit)
902
903 ! Local part of the pseudopotential
904 iunit = io_open(trim(dir)//'/local_long_range_density', namespace, action='write')
905 call spline_print(ps%nlr, iunit)
906 call io_close(iunit)
907
908 ! Fourier transform of the local part
909 iunit = io_open(trim(dir)//'/local_ft', namespace, action='write')
910 safe_allocate(fw(1, 1))
911 call spline_init(fw(1, 1))
912 call spline_3dft(ps%vl, fw(1, 1), ps%projectors_sphere_threshold, gmax = gmax)
913 call spline_print(fw(1, 1), iunit)
914 call spline_end(fw(1, 1))
915 safe_deallocate_a(fw)
916 call io_close(iunit)
917
918 ! Kleinman-Bylander projectors
919 iunit = io_open(trim(dir)//'/nonlocal', namespace, action='write')
920 call spline_print(ps%kb, iunit)
921 call io_close(iunit)
922
923 iunit = io_open(trim(dir)//'/nonlocal_derivative', namespace, action='write')
924 call spline_print(ps%dkb, iunit)
925 call io_close(iunit)
926
927 iunit = io_open(trim(dir)//'/nonlocal_ft', namespace, action='write')
928 safe_allocate(fw(0:ps%lmax, 1:ps%kbc))
929 call spline_init(fw)
930 do k = 0, ps%lmax
931 do j = 1, ps%kbc
932 call spline_besselft(ps%kb(k, j), fw(k, j), k, threshold=ps%projectors_sphere_threshold, &
933 gmax=m_four*gmax)
934 end do
935 end do
936 call spline_print(fw, iunit)
937 call spline_end(fw)
938 safe_deallocate_a(fw)
939 call io_close(iunit)
940
941 ! Pseudo-wavefunctions
942 iunit = io_open(trim(dir)//'/wavefunctions', namespace, action='write')
943 call spline_print(ps%ur, iunit)
944 call io_close(iunit)
945
946 ! Density
947 if (ps%has_density) then
948 iunit = io_open(trim(dir)//'/density', namespace, action='write')
949 call spline_print(ps%density, iunit)
950 call io_close(iunit)
951
952 iunit = io_open(trim(dir)//'/density_derivative', namespace, action='write')
953 call spline_print(ps%density_der, iunit)
954 call io_close(iunit)
955 end if
956
957 ! Non-linear core-corrections
958 if (ps_has_nlcc(ps)) then
959 iunit = io_open(trim(dir)//'/nlcc', namespace, action='write')
960 call spline_print(ps%core, iunit)
961 call io_close(iunit)
962 end if
963
964 pop_sub(ps_debug)
965 end subroutine ps_debug
966
967
968 ! ---------------------------------------------------------
969 subroutine ps_end(ps)
970 type(ps_t), intent(inout) :: ps
971
972 if (.not. allocated(ps%kb)) return
973
974 push_sub(ps_end)
975
976 if (ps%is_separated) then
977 call spline_end(ps%vlr)
978 call spline_end(ps%nlr)
979 end if
980
981 call spline_end(ps%kb)
982 call spline_end(ps%dkb)
983 call spline_end(ps%ur)
984
985 call spline_end(ps%vl)
986 call spline_end(ps%core)
987 call spline_end(ps%core_der)
988
989 if (allocated(ps%density)) call spline_end(ps%density)
990 if (allocated(ps%density_der)) call spline_end(ps%density_der)
991
992 call logrid_end(ps%g)
993
994 safe_deallocate_a(ps%kb)
995 safe_deallocate_a(ps%dkb)
996 safe_deallocate_a(ps%ur)
997 safe_deallocate_a(ps%bound)
998 safe_deallocate_a(ps%h)
999 safe_deallocate_a(ps%k)
1000 safe_deallocate_a(ps%density)
1001 safe_deallocate_a(ps%density_der)
1002
1003 pop_sub(ps_end)
1004 end subroutine ps_end
1005
1006
1007 ! ---------------------------------------------------------
1008 subroutine hgh_load(ps, ps_hgh)
1009 type(ps_t), intent(inout) :: ps
1010 type(ps_hgh_t), intent(inout) :: ps_hgh
1011
1012 push_sub(hgh_load)
1013
1014 ! Fixes some components of ps
1015 ps%nlcc = .false.
1016 if (ps%lmax >= 0) then
1017 ps%rc_max = 1.1_real64 * maxval(ps_hgh%kbr(0:ps%lmax)) ! Increase a little.
1018 else
1019 ps%rc_max = m_zero
1020 end if
1021 ps%h(0:ps%lmax, 1:ps%kbc, 1:ps%kbc) = ps_hgh%h(0:ps%lmax, 1:ps%kbc, 1:ps%kbc)
1022 ps%k(0:ps%lmax, 1:ps%kbc, 1:ps%kbc) = ps_hgh%k(0:ps%lmax, 1:ps%kbc, 1:ps%kbc)
1023
1024 ! now we fit the splines
1025 call get_splines()
1026
1027 pop_sub(hgh_load)
1028
1029 contains
1030
1031 ! ---------------------------------------------------------
1032 subroutine get_splines()
1033 integer :: l, is, j, ip
1034 real(real64), allocatable :: hato(:), dens(:)
1035
1036 push_sub(hgh_load.get_splines)
1037
1038 safe_allocate(hato(1:ps_hgh%g%nrval))
1039 safe_allocate(dens(1:ps_hgh%g%nrval))
1040
1041 ! Interpolate the KB-projection functions
1042 do l = 0, ps_hgh%l_max
1043 do j = 1, 3
1044 hato(:) = ps_hgh%kb(:, l, j)
1045 call spline_fit(ps_hgh%g%nrval, ps_hgh%g%rofi, hato, ps%kb(l, j), ps%projectors_sphere_threshold)
1046 end do
1047 end do
1048
1049 ! Now the part corresponding to the local pseudopotential
1050 ! where the asymptotic part is subtracted
1051 call spline_fit(ps_hgh%g%nrval, ps_hgh%g%rofi, ps_hgh%vlocal, ps%vl, ps%projectors_sphere_threshold)
1052
1053 ! Define the table for the pseudo-wavefunction components (using splines)
1054 ! with a correct normalization function
1055 do is = 1, ps%ispin
1056 dens = m_zero
1057 do l = 1, ps%conf%p
1058 hato(2:ps_hgh%g%nrval) = ps_hgh%rphi(2:ps_hgh%g%nrval, l)/ps_hgh%g%rofi(2:ps_hgh%g%nrval)
1059 hato(1) = hato(2)
1060
1061 do ip = 1, ps_hgh%g%nrval
1062 dens(ip) = dens(ip) + ps%conf%occ(l, is)*hato(ip)**2/(m_four*m_pi)
1063 end do
1065 call spline_fit(ps_hgh%g%nrval, ps_hgh%g%rofi, hato, ps%ur(l, is), ps%projectors_sphere_threshold)
1066 end do
1067 call spline_fit(ps_hgh%g%nrval, ps_hgh%g%rofi, dens, ps%density(is), ps%projectors_sphere_threshold)
1068 end do
1069
1070 safe_deallocate_a(hato)
1071 safe_deallocate_a(dens)
1072
1073 pop_sub(hgh_load.get_splines)
1074 end subroutine get_splines
1075 end subroutine hgh_load
1076
1077
1078 ! ---------------------------------------------------------
1079 subroutine ps_grid_load(ps, ps_grid)
1080 type(ps_t), intent(inout) :: ps
1081 type(ps_in_grid_t), intent(in) :: ps_grid
1082
1083 push_sub(ps_grid_load)
1084
1085 ! Fixes some components of ps, read in ps_grid
1086 ps%z_val = ps_grid%zval
1087
1088 ps%nlcc = ps_grid%core_corrections
1089
1090 ps%h(0:ps%lmax, 1, 1) = ps_grid%dkbcos(1:ps%lmax+1)
1091
1092 ! Increasing radius a little, just in case.
1093 ! I have hard-coded a larger increase of the cutoff for the filtering.
1094 ps%rc_max = maxval(ps_grid%kb_radius(1:ps%lmax+1)) * 1.5_real64
1095
1096 ! now we fit the splines
1097 call get_splines(ps_grid%g)
1098
1099 ! Passes from Rydbergs to Hartrees.
1100 ps%h(0:ps%lmax,:,:) = ps%h(0:ps%lmax,:,:) / m_two
1101
1102 pop_sub(ps_grid_load)
1104 contains
1105
1106 subroutine get_splines(g)
1107 type(logrid_t), intent(in) :: g
1108
1109 real(real64), allocatable :: hato(:), dens(:)
1110 integer :: is, l, ir, nrc, ip
1111
1112 push_sub(ps_grid_load.get_splines)
1113
1114 safe_allocate(hato(1:g%nrval))
1115 safe_allocate(dens(1:g%nrval))
1116
1117 ! the wavefunctions
1118 do is = 1, ps%ispin
1119
1120 dens = m_zero
1121
1122 do l = 1, ps_grid%no_l_channels
1123 hato(2:) = ps_grid%rphi(2:, l, 1+is)/g%rofi(2:)
1124 hato(1) = first_point_extrapolate(g%rofi, hato)
1125
1126 if(ps%conf%occ(l, is) > m_epsilon) then
1127 do ip = 1, g%nrval
1128 dens(ip) = dens(ip) + ps%conf%occ(l, is)*hato(ip)**2/(m_four*m_pi)
1129 end do
1130 end if
1131
1132 call spline_fit(g%nrval, g%rofi, hato, ps%ur(l, is), ps%projectors_sphere_threshold)
1133
1134 end do
1135 call spline_fit(g%nrval, g%rofi, dens, ps%density(is), ps%projectors_sphere_threshold)
1136 end do
1137
1138
1139 ! the Kleinman-Bylander projectors
1140 do l = 1, ps%lmax+1
1141 nrc = logrid_index(g, ps_grid%kb_radius(l)) + 1
1142 hato(1:nrc) = ps_grid%KB(1:nrc, l)
1143 hato(nrc+1:g%nrval) = m_zero
1144
1145 call spline_fit(g%nrval, g%rofi, hato, ps%kb(l-1, 1), ps%projectors_sphere_threshold)
1146 end do
1147
1148 ! Now the part corresponding to the local pseudopotential
1149 ! where the asymptotic part is subtracted
1150 hato(:) = ps_grid%vlocal(:)/m_two
1151 call spline_fit(g%nrval, g%rofi, hato, ps%vl, ps%projectors_sphere_threshold)
1152
1153 if (ps_grid%core_corrections) then
1154 ! find cutoff radius
1155 hato(2:) = ps_grid%chcore(2:)/(m_four*m_pi*g%r2ofi(2:))
1156
1157 do ir = g%nrval-1, 2, -1
1158 if (hato(ir) > eps) then
1159 nrc = ir + 1
1160 exit
1161 end if
1162 end do
1163
1164 hato(nrc:g%nrval) = m_zero
1165 hato(1) = first_point_extrapolate(g%rofi, hato)
1166
1167 call spline_fit(g%nrval, g%rofi, hato, ps%core, ps%projectors_sphere_threshold)
1168 end if
1169
1170 safe_deallocate_a(hato)
1171 safe_deallocate_a(dens)
1172
1173 pop_sub(ps_grid_load.get_splines)
1174 end subroutine get_splines
1175 end subroutine ps_grid_load
1176
1177 ! ---------------------------------------------------------
1181 subroutine ps_xml_load(ps, ps_xml, namespace)
1182 type(ps_t), intent(inout) :: ps
1183 type(ps_xml_t), intent(in) :: ps_xml
1184 type(namespace_t), intent(in) :: namespace
1185
1186 integer :: ll, ip, is, ic, jc, ir, nrc, ii
1187 real(real64) :: rr, kbcos, kbnorm, dnrm, avgv, volume_element
1188 real(real64), allocatable :: vlocal(:), kbprojector(:), wavefunction(:), nlcc_density(:), dens(:,:)
1189 integer, allocatable :: cmap(:, :)
1190 real(real64), allocatable :: matrix(:, :), eigenvalues(:)
1191 logical :: density_is_known
1192 real(real64), parameter :: tol_diagonal=1.0e-10_real64 ! tolerance for taking d_ij as a diagonal matrix
1193
1194 push_sub(ps_xml_load)
1195
1196 ! Vanderbilt-Kleinman-Bylander, i.e. more than one projector per l. Hamann uses 2 for ONCVPSP.
1197 if (pseudo_has_total_angular_momentum(ps_xml%pseudo)) then
1198 if (ps%file_format == pseudo_format_psp8) then
1199 ps%relativistic_treatment = proj_j_average
1200 message(1) = "SOC from PSP8 is not currently supported."
1201 message(2) = "Only scalar relativistic effects will be considered."
1202 call messages_warning(2, namespace=namespace)
1203 else
1204 ps%relativistic_treatment = proj_j_dependent
1205 end if
1206 end if
1207
1208 ps%nlcc = ps_xml%nlcc
1209
1210 ps%z_val = ps_xml%valence_charge
1211
1212 ps%has_density = ps_xml%has_density
1213
1214 ! We start with the local potential
1215
1216 ! Outside of the grid, we use the analytical long-range tail
1217 ! TODO: check usefulness: ps%g%nrval is initialized as ps_xml%grid_size
1218 ! This makes no sense, ps_xml%grid can only go up to ps_xml%grid_size
1219 ! We should just fit ps_xml%potential directly
1220 ! Also, if ip > ps_xml%grid_size, most likely ps_xml%grid is not defined
1221 safe_allocate(vlocal(1:ps%g%nrval))
1222 do ip = 1, ps_xml%grid_size
1223 rr = ps_xml%grid(ip)
1224 if (ip <= ps_xml%grid_size) then
1225 vlocal(ip) = ps_xml%potential(ip, ps%llocal)
1226 else
1227 vlocal(ip) = -ps_xml%valence_charge/rr
1228 end if
1229 end do
1230
1231 ! We then fit this by a spline
1232 call spline_fit(ps%g%nrval, ps%g%rofi, vlocal, ps%vl, ps%projectors_sphere_threshold)
1233
1234 safe_deallocate_a(vlocal)
1235
1236 safe_allocate(kbprojector(1:ps%g%nrval))
1237 safe_allocate(wavefunction(1:ps%g%nrval))
1238
1239 kbprojector = m_zero
1240 wavefunction = m_zero
1241
1242 density_is_known = .false.
1243
1244 ! We then proceed with the nonlocal projectors and the orbitals
1245 if (ps_xml%kleinman_bylander) then
1246
1247 safe_allocate(cmap(0:ps_xml%lmax, 1:ps_xml%nchannels))
1248
1249 ! We order the different projectors and create some mappings
1250 ! the order of the channels is determined by spin orbit and the j value
1251 do ll = 0, ps_xml%lmax
1252 do ic = 1, ps_xml%nchannels
1253 cmap(ll, ic) = ic
1254
1255 if (ll == 0) cycle
1256 if (ll == ps_xml%llocal) cycle
1257 if (ps%relativistic_treatment /= proj_j_dependent) cycle
1258 if (ic > pseudo_nprojectors_per_l(ps_xml%pseudo, ll)) cycle ! This occurs for O and F for pseudodojo FR PBE for instance
1259
1260 ! This is Octopus convention:
1261 ! We treat l+1/2 first and l-1/2 after, so we order the projectors this way
1262
1263 assert(mod(ps_xml%nchannels, 2)==0)
1264 if (pseudo_projector_2j(ps_xml%pseudo, ll, ic) == 2*ll - 1) then
1265 cmap(ll, ic) = int((ic-1)/2)*2 + 2
1266 else
1267 assert(pseudo_projector_2j(ps_xml%pseudo, ll, ic) == 2*ll + 1)
1268 cmap(ll, ic) = int((ic-1)/2)*2 + 1
1269 end if
1270
1271 end do
1272
1273 ! check that all numbers are present for each l
1274 assert(sum(cmap(ll, 1:ps_xml%nchannels)) == (ps_xml%nchannels + 1)*ps_xml%nchannels/2)
1275 end do
1277 assert(all(cmap >= 0 .and. cmap <= ps_xml%nchannels))
1278
1279
1280 ! Weight of the projectors
1281 ! This is a diagonal matrix after the following treatment
1282 ps%h = m_zero
1283
1284 ! We now take the matrix dij and we diagonalize it
1285 ! The nonlocal projectors are changed accordingly, and then fitted by splines
1286 if (pseudo_nprojectors(ps_xml%pseudo) > 0) then
1287 safe_allocate(matrix(1:ps_xml%nchannels, 1:ps_xml%nchannels))
1288 safe_allocate(eigenvalues(1:ps_xml%nchannels))
1289
1290 do ll = 0, ps_xml%lmax
1291
1292 if (is_diagonal(ps_xml%nchannels, ps_xml%dij(ll, :, :), tol_diagonal) .or. &
1293 pseudo_has_total_angular_momentum(ps_xml%pseudo)) then
1294 matrix = m_zero
1295 do ic = 1, ps_xml%nchannels
1296 eigenvalues(ic) = ps_xml%dij(ll, ic, ic)
1297 matrix(ic, ic) = m_one
1298 end do
1299 else
1300 ! diagonalize the coefficient matrix
1301 matrix(1:ps_xml%nchannels, 1:ps_xml%nchannels) = ps_xml%dij(ll, 1:ps_xml%nchannels, 1:ps_xml%nchannels)
1302 call lalg_eigensolve(ps_xml%nchannels, matrix, eigenvalues)
1303 end if
1304
1305 do ic = 1, ps_xml%nchannels
1306 kbprojector = m_zero
1307 do jc = 1, ps_xml%nchannels
1308 call lalg_axpy(ps_xml%grid_size, matrix(jc, ic), ps_xml%projector(:, ll, jc), kbprojector)
1309 end do
1310
1311 call spline_fit(ps%g%nrval, ps%g%rofi, kbprojector, ps%kb(ll, cmap(ll, ic)), ps%projectors_sphere_threshold)
1312
1313 ps%h(ll, cmap(ll, ic), cmap(ll, ic)) = eigenvalues(ic)
1314
1315 end do
1316 end do
1317
1318 safe_deallocate_a(matrix)
1319 safe_deallocate_a(eigenvalues)
1320 end if
1321
1322 ps%conf%p = ps_xml%nwavefunctions
1323
1324 ! If we do not have a density but we have wavefunctions, we compute the
1325 ! pseudo-atomic density, from the pseudo-atomic wavefunctions
1326 ! In case we are doing a spin-polarized calculation, it is better to use the
1327 ! wavefunctions with spin-dependent occupations, instead the spin unresolved density
1328 ! given in the pseudopotential
1329 if ((.not. ps_has_density(ps) .or. ps%ispin == 2) .and. ps_xml%nwavefunctions > 0) then
1330 safe_allocate(dens(1:ps%g%nrval, 1:ps%ispin))
1331 dens = m_zero
1332 end if
1333
1334 do ii = 1, ps_xml%nwavefunctions
1335
1336 ps%conf%n(ii) = ps_xml%wf_n(ii)
1337 ps%conf%l(ii) = ps_xml%wf_l(ii)
1338
1339 if (ps%ispin == 2) then
1340 ps%conf%occ(ii, 1) = min(ps_xml%wf_occ(ii), m_two*ps_xml%wf_l(ii) + m_one)
1341 ps%conf%occ(ii, 2) = ps_xml%wf_occ(ii) - ps%conf%occ(ii, 1)
1342 else
1343 ps%conf%occ(ii, 1) = ps_xml%wf_occ(ii)
1344 end if
1345
1346 ps%conf%j(ii) = m_zero
1347 if (pseudo_has_total_angular_momentum(ps_xml%pseudo)) then
1348 ps%conf%j(ii) = m_half*pseudo_wavefunction_2j(ps_xml%pseudo, ii)
1349 end if
1350
1351 do ip = 1, ps%g%nrval
1352 if (ip <= ps_xml%grid_size) then
1353 wavefunction(ip) = ps_xml%wavefunction(ip, ii)
1354 else
1355 wavefunction(ip) = m_zero
1356 end if
1357 end do
1358
1359 do is = 1, ps%ispin
1360 call spline_fit(ps%g%nrval, ps%g%rofi, wavefunction, ps%ur(ii, is), ps%projectors_sphere_threshold)
1361 end do
1362
1363 if (.not. ps_has_density(ps) .or. ps%ispin == 2) then
1364 do is = 1, ps%ispin
1365 do ip = 1, ps_xml%grid_size
1366 dens(ip, is) = dens(ip, is) + ps%conf%occ(ii, is)*wavefunction(ip)**2/(m_four*m_pi)
1367 end do
1368 end do
1369 end if
1370
1371 end do
1372
1373 !If we assigned all the valence electrons, we can compute the (spin-resolved) atomic density
1374 if ((.not. ps_has_density(ps) .or. ps%ispin == 2) .and. ps_xml%nwavefunctions > 0) then
1375 do is = 1, ps%ispin
1376 call spline_fit(ps%g%nrval, ps%g%rofi, dens(:,is), ps%density(is), ps%projectors_sphere_threshold)
1377 end do
1378 safe_deallocate_a(dens)
1379 density_is_known = .true.
1380 ps%has_density = .true.
1381 end if
1382
1383 safe_deallocate_a(cmap)
1384
1385 else !Not ps_xml%kleinman_bylander
1386
1387 !Get the occupations from the valence charge of the atom
1388 ps%conf%occ = m_zero
1389 ps%conf%symbol = ps%label(1:2)
1390 call ps_guess_atomic_occupations(namespace, ps%z, ps%z_val, ps%ispin, ps%conf)
1391
1392 ! In order to work in the following, we need to sort the occupations by angular momentum
1393 call valconf_sort_by_l(ps%conf, ps_xml%lmax)
1394
1395 !If we assigned all the valence electrons, we can compute the (spin-resolved) atomic density
1396 if (abs(sum(ps%conf%occ) - ps%z_val ) < m_epsilon) then
1397 safe_allocate(dens(1:ps%g%nrval, 1:ps%ispin))
1398 dens = m_zero
1399 end if
1400
1401 do ll = 0, ps_xml%lmax
1402 ! we need to build the KB projectors
1403 ! the procedure was copied from ps_in_grid.F90 (r12967)
1404 dnrm = m_zero
1405 avgv = m_zero
1406 do ip = 1, ps_xml%grid_size
1407 rr = ps_xml%grid(ip)
1408 volume_element = rr**2*ps_xml%weights(ip)
1409 kbprojector(ip) = (ps_xml%potential(ip, ll) - ps_xml%potential(ip, ps%llocal))*ps_xml%wavefunction(ip, ll)
1410 dnrm = dnrm + kbprojector(ip)**2*volume_element
1411 avgv = avgv + kbprojector(ip)*ps_xml%wavefunction(ip, ll)*volume_element
1412 end do
1413
1414 kbcos = dnrm/(safe_tol(avgv,1.0e-20_real64))
1415 kbnorm = m_one/(safe_tol(sqrt(dnrm),1.0e-20_real64))
1416
1417 if (ll /= ps%llocal) then
1418 ps%h(ll, 1, 1) = kbcos
1419 kbprojector = kbprojector*kbnorm
1420 else
1421 ps%h(ll, 1, 1) = m_zero
1422 end if
1423
1424 call spline_fit(ps%g%nrval, ps%g%rofi, kbprojector, ps%kb(ll, 1), ps%projectors_sphere_threshold)
1425
1426 ! wavefunctions, for the moment we pad them with zero
1427 call lalg_copy(ps_xml%grid_size, ps_xml%wavefunction(:, ll), wavefunction)
1428 wavefunction(ps_xml%grid_size+1:ps%g%nrval) = m_zero
1429
1430 do is = 1, ps%ispin
1431 call spline_fit(ps%g%nrval, ps%g%rofi, wavefunction, ps%ur(ll + 1, is), ps%projectors_sphere_threshold)
1432 end do
1433
1434 !If we assigned all the valence electrons, we can compute the (spin-resolved) atomic density
1435 if (abs(sum(ps%conf%occ) - ps%z_val) < m_epsilon) then
1436 do is = 1, ps%ispin
1437 do ip = 1, ps_xml%grid_size
1438 dens(ip, is) = dens(ip, is) + ps%conf%occ(ll+1, is)*wavefunction(ip)**2/(m_four*m_pi)
1439 end do
1440 end do
1441 end if
1442 end do
1443
1444 !If we assigned all the valence electrons, we can compute the (spin-resolved) atomic density
1445 if (abs(sum(ps%conf%occ) - ps%z_val) < m_epsilon) then
1446 do is = 1, ps%ispin
1447 call spline_fit(ps%g%nrval, ps%g%rofi, dens(:,is), ps%density(is), ps%projectors_sphere_threshold)
1448 end do
1449 safe_deallocate_a(dens)
1450 density_is_known = .true.
1451 ps%has_density = .true.
1452 end if
1453 end if ! ps_xml%kleinman_bylander
1454
1455
1456 if (ps_has_density(ps) .and. .not. density_is_known) then
1457
1458 safe_allocate(dens(1:ps%g%nrval, 1))
1459
1460 dens(1:ps_xml%grid_size, 1) = ps_xml%density(1:ps_xml%grid_size)/ps%ispin
1461 dens(ps_xml%grid_size + 1:ps%g%nrval, 1) = m_zero
1462
1463 do is = 1, ps%ispin
1464 call spline_fit(ps%g%nrval, ps%g%rofi, dens(:, 1), ps%density(is), ps%projectors_sphere_threshold)
1465 end do
1466
1467 safe_deallocate_a(dens)
1468 end if
1469
1470 ! Non-linear core-corrections
1471 ! We truncate the NLCC when below eps
1472 if (ps_xml%nlcc) then
1473
1474 safe_allocate(nlcc_density(1:ps%g%nrval))
1475
1476 call lalg_copy(ps_xml%grid_size, ps_xml%nlcc_density, nlcc_density)
1477
1478 ! find cutoff radius
1479 do ir = ps_xml%grid_size - 1, 1, -1
1480 if (nlcc_density(ir) > eps) then
1481 nrc = ir + 1
1482 exit
1483 end if
1484 end do
1485
1486 nlcc_density(nrc:ps%g%nrval) = m_zero
1487
1488 call spline_fit(ps%g%nrval, ps%g%rofi, nlcc_density, ps%core, ps%projectors_sphere_threshold)
1489
1490 safe_deallocate_a(nlcc_density)
1491 end if
1492
1493 call ps_getradius(ps)
1494
1495 !To be consistent with the other pseudopotentials, we are increasing the radius here
1496 ps%rc_max = ps%rc_max*1.05_real64
1497
1498 safe_deallocate_a(kbprojector)
1499 safe_deallocate_a(wavefunction)
1500
1501 pop_sub(ps_xml_load)
1502 end subroutine ps_xml_load
1503
1504
1505 ! ---------------------------------------------------------
1507 pure integer function ps_niwfs(ps)
1508 type(ps_t), intent(in) :: ps
1509
1510 integer :: i, l
1511
1512 ps_niwfs = 0
1513 do i = 1, ps%conf%p
1514 l = ps%conf%l(i)
1515 ps_niwfs = ps_niwfs + (2*l+1)
1516 end do
1517
1518 end function ps_niwfs
1519
1520 ! ---------------------------------------------------------
1522 pure integer function ps_bound_niwfs(ps)
1523 type(ps_t), intent(in) :: ps
1524
1525 integer :: i, l
1526
1527 ps_bound_niwfs = 0
1528 do i = 1, ps%conf%p
1529 l = ps%conf%l(i)
1530 if (any(.not. ps%bound(i,:))) cycle
1531 ps_bound_niwfs = ps_bound_niwfs + (2*l+1)
1532 end do
1533
1534 end function ps_bound_niwfs
1535
1536 !---------------------------------------
1537
1538 pure logical function ps_has_density(ps) result(has_density)
1539 type(ps_t), intent(in) :: ps
1540
1541 has_density = ps%has_density
1542
1543 end function ps_has_density
1544
1545 !---------------------------------------
1546
1547 pure logical function ps_has_nlcc(ps) result(has_nlcc)
1548 type(ps_t), intent(in) :: ps
1549
1550 has_nlcc = ps%nlcc
1551
1552 end function ps_has_nlcc
1553
1554 !---------------------------------------
1555 real(real64) function ps_density_volume(ps, namespace) result(volume)
1556 type(ps_t), intent(in) :: ps
1557 type(namespace_t), intent(in) :: namespace
1558
1559 integer :: ip, ispin
1560 real(real64) :: rr
1561 real(real64), allocatable ::vol(:)
1562 type(spline_t) :: volspl
1563
1564 push_sub(ps_density_volume)
1565
1566 if (.not. ps_has_density(ps)) then
1567 message(1) = "The pseudopotential does not contain an atomic density"
1568 call messages_fatal(1, namespace=namespace)
1569 end if
1570
1571 safe_allocate(vol(1:ps%g%nrval))
1572
1573 do ip = 1, ps%g%nrval
1574 rr = ps%g%rofi(ip)
1575 vol(ip) = m_zero
1576 do ispin = 1, ps%ispin
1577 vol(ip) = vol(ip) + spline_eval(ps%density(ispin), rr)*m_four*m_pi*rr**5
1578 end do
1579 end do
1580
1581 call spline_init(volspl)
1582 call spline_fit(ps%g%nrval, ps%g%rofi, vol, volspl, ps%projectors_sphere_threshold)
1583 volume = spline_integral(volspl)
1584 call spline_end(volspl)
1585
1586 safe_deallocate_a(vol)
1587
1588 pop_sub(ps_density_volume)
1589 end function ps_density_volume
1590
1594 subroutine ps_guess_atomic_occupations(namespace, zz, valcharge, ispin, conf)
1595 type(namespace_t), intent(in) :: namespace
1596 real(real64), intent(in) :: zz
1597 real(real64), intent(in) :: valcharge
1598 integer, intent(in) :: ispin
1599 type(valconf_t), intent(inout) :: conf
1600
1601 real(real64) :: val
1604
1605 val = valcharge
1606 conf%p = 0
1607 conf%l = 0
1608 conf%j = m_zero
1609 conf%n = 0
1610
1611 write(message(1), '(a,a,a)') 'Debug: Guessing the atomic occupations for ', trim(conf%symbol), "."
1612 call messages_info(1, namespace=namespace, debug_only=.true.)
1613
1614 assert(valcharge <= zz)
1615
1616 ! Here we populate the core states
1617 ! 1s state for all atoms after He
1618 if(int(zz) > 2 .and. val > zz - 2) then
1619 call fill_s_orbs(val, 2, 1)
1620 end if
1621 ! 2s state for all atoms after Be
1622 if(int(zz) > 4 .and. val > zz - 4) then
1623 call fill_s_orbs(val, 2, 2)
1624 end if
1625 ! 2p state for all atoms after Ne
1626 ! For pseudopotentials Al-Ar, we fill the 2s but not the 2p
1627 if(int(zz) > 18 .and. val > zz - 10) then
1628 call fill_p_orbs(val, 6, 2)
1629 end if
1630 print *, val, zz - 12, int(zz) > 12 .and. val > zz - 12
1631 ! 3s state for all atoms after Mg
1632 if(int(zz) > 12 .and. val > zz - 12) then
1633 print *, 'Filling 3s2'
1634 call fill_s_orbs(val, 2, 3)
1635 end if
1636 ! 3p state for all atoms after Ar
1637 if(int(zz) > 18 .and. val > zz - 18) then
1638 print *, 'Filling 3p6'
1639 call fill_p_orbs(val, 6, 3)
1640 end if
1641 ! 3d states for all atoms after Ni
1642 if(int(zz) > 28 .and. val > zz - 28) then
1643 call fill_d_orbs(val, 10, 3)
1644 end if
1645 ! 4s states for all atoms after Zn
1646 if(int(zz) > 30 .and. val > zz - 30) then
1647 call fill_s_orbs(val, 2, 4)
1648 end if
1649 ! 4p states for all atoms after Kr
1650 if(int(zz) > 36 .and. val > zz - 36) then
1651 call fill_p_orbs(val, 6, 4)
1652 end if
1653 ! 4d states for all atoms after Pd
1654 if(int(zz) > 46 .and. val > zz - 46) then
1655 call fill_d_orbs(val, 10, 4)
1656 end if
1657 ! 5s states for all atoms after Cd
1658 ! For Z=71 ot Z=80, the 4f is filled before the 5s/5p
1659 if((int(zz) > 48 .and. val > zz - 48) .or. &
1660 (int(zz) > 70 .and. int(zz) <= 81 .and. val > zz - 62)) then
1661 call fill_s_orbs(val, 2, 5)
1662 end if
1663 ! 5p states for all atoms after Xe
1664 ! For Z=71 ot Z=80, the 4f is filled before the 5s/5p
1665 if((int(zz) > 54 .and. val > zz - 54) .or. &
1666 (int(zz) > 70 .and. int(zz) <= 81 .and. val > zz - 68) ) then
1667 call fill_p_orbs(val, 6, 5)
1668 end if
1669 ! 4f states for all atoms after Yb
1670 ! Only after Z=80 for pseudopotentials
1671 if(int(zz) > 80 .and. val > zz - 68) then
1672 call fill_f_orbs(val, 14, 4)
1673 end if
1674
1675
1676 ! We leave here the valence states
1677 select case (int(zz))
1678 case (1)
1679 call fill_s_orbs(val, 1, 1) ! H 1s^1
1680 case (2)
1681 call fill_s_orbs(val, 2, 1) ! He 1s^2
1682 case (3)
1683 call fill_s_orbs(val, 1, 2) ! Li 1s^2 2s^1
1684 case (4)
1685 call fill_s_orbs(val, 2, 2) ! Be 1s^2 ; 2s^2
1686 case (5)
1687 call fill_p_orbs(val, 1, 2) ! B 1s^2 ; 2s^2 2p^1
1688 case (6)
1689 call fill_p_orbs(val, 2, 2) ! C 1s^2 ; 2s^2 2p^2
1690 case (7)
1691 call fill_p_orbs(val, 3, 2) ! N 1s^2 ; 2s^2 2
1692 case (8)
1693 call fill_p_orbs(val, 4, 2) ! O 1s^2 ; 2s^2 2p^4
1694 case (9)
1695 call fill_p_orbs(val, 5, 2) ! F 1s^2 ; 2s^2 2p^5
1696 case (10)
1697 call fill_p_orbs(val, 6, 2) ! Ne 1s^2 ; 2s^2 2p^6
1698 case (11)
1699 if(val > 6) call fill_p_orbs(val, 6, 2)
1700 call fill_s_orbs(val, 1, 3) ! Na 1s^2 2s^2 2p^6 ; 3s^1
1701 case (12)
1702 if(val > 6) call fill_p_orbs(val, 6, 2)
1703 call fill_s_orbs(val, 2, 3) ! Mg 1s^2 2s^2 2p^6 ; 3s^2
1704 case (13)
1705 if(val > 6) call fill_p_orbs(val, 6, 2)
1706 call fill_p_orbs(val, 1, 3) ! Al 1s^2 2s^2 2p^6 ; 3s^2 3p^1
1707 case (14)
1708 if(val > 6) call fill_p_orbs(val, 6, 2)
1709 call fill_p_orbs(val, 2, 3) ! Si 1s^2 2s^2 2p^6 ; 3s^2 3p^2
1710 case (15)
1711 if(val > 6) call fill_p_orbs(val, 6, 2)
1712 call fill_p_orbs(val, 3, 3) ! P 1s^2 2s^2 2p^6 ; 3s^2 3p^3
1713 case (16)
1714 if(val > 6) call fill_p_orbs(val, 6, 2)
1715 call fill_p_orbs(val, 4, 3) ! S 1s^2 2s^2 2p^6 ; 3s^2 3p^4
1716 case (17)
1717 if(val > 6) call fill_p_orbs(val, 6, 2)
1718 call fill_p_orbs(val, 5, 3) ! Cl 1s^2 2s^2 2p^6 ; 3s^2 3p^5
1719 case (18)
1720 if(val > 6) call fill_p_orbs(val, 6, 2)
1721 call fill_p_orbs(val, 6, 3) ! Ar 1s^2 2s^2 2p^6 ; 3s^2 3p^6
1722 case (19)
1723 call fill_s_orbs(val, 1, 4) ! K 1s^2 2s^2 2p^6 3s^2 3p^6 ; 4s^1
1724 case (20)
1725 call fill_s_orbs(val, 2, 4) ! Ca 1s^2 2s^2 2p^6 3s^2 3p^6 ; 4s^2
1726 case (21)
1727 if (val > 1) call fill_s_orbs(val, 2, 4)
1728 call fill_d_orbs(val, 1, 3) ! Sc 1s^2 2s^2 2p^6 3s^2 3p^6 ; 4s^2 3d^1
1729 case (22)
1730 if (val > 2) call fill_s_orbs(val, 2, 4)
1731 call fill_d_orbs(val, 2, 3) ! Ti 1s^2 2s^2 2p^6 3s^2 3p^6 ; 4s^2 3d^2
1732 case (23)
1733 print *, 'V', val
1734 if (val > 3) call fill_s_orbs(val, 2, 4)
1735 call fill_d_orbs(val, 3, 3) ! V 1s^2 2s^2 2p^6 3s^2 3p^6 ; 4s^2 3d^3
1736 case (24)
1737 if (val > 4) call fill_s_orbs(val, 2, 4)
1738 call fill_d_orbs(val, 4, 3) ! Cr 1s^2 2s^2 2p^6 3s^2 3p^6 ; 4s^2 3d^4
1739 case (25)
1740 if (val > 5) call fill_s_orbs(val, 2, 4)
1741 call fill_d_orbs(val, 5, 3) ! Mn 1s^2 2s^2 2p^6 3s^2 3p^6 ; 4s^2 3d^5
1742 case (26)
1743 if (val > 6) call fill_s_orbs(val, 2, 4)
1744 call fill_d_orbs(val, 6, 3) ! Fe 1s^2 2s^2 2p^6 3s^2 3p^6 ; 4s^2 3d^6
1745 case (27)
1746 if (val > 7) call fill_s_orbs(val, 2, 4)
1747 call fill_d_orbs(val, 7, 3) ! Co 1s^2 2s^2 2p^6 3s^2 3p^6 ; 4s^2 3d^7
1748 case (28)
1749 if (val > 8 ) call fill_s_orbs(val, 2, 4)
1750 call fill_d_orbs(val, 8, 3) ! Ni 1s^2 2s^2 2p^6 3s^2 3p^6 ; 4s^2 3d^8
1751 case (29)
1752 call fill_s_orbs(val, 1, 4) ! Cu 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 ; 4s^1
1753 case (30)
1754 call fill_s_orbs(val, 2, 4) ! Zn 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 ; 4s^2
1755 case (31)
1756 call fill_p_orbs(val, 1, 4) ! Ga 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 ; 4s^2 4p^1
1757 case (32)
1758 call fill_p_orbs(val, 2, 4) ! Ge 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 ; 4s^2 4p^2
1759 case (33)
1760 call fill_p_orbs(val, 3, 4) ! As 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 ; 4s^2 4p^3
1761 case (34)
1762 call fill_p_orbs(val, 4, 4) ! Se 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 ; 4s^2 4p^4
1763 case (35)
1764 call fill_p_orbs(val, 5, 4) ! Br 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 ; 4s^2 4p^5
1765 case (36)
1766 call fill_p_orbs(val, 6, 4) ! Kr 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 ; 4s^2 4p^6
1767 case (37)
1768 call fill_s_orbs(val, 1, 5) ! Rb 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 4s^2 4p^6 ; 5s^1
1769 case (38)
1770 call fill_s_orbs(val, 2, 5) ! Sr 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 4s^2 4p^6 ; 5s^2
1771 case (39)
1772 if (val > 2) call fill_d_orbs(val, 1, 4)
1773 call fill_s_orbs(val, 2, 5) ! Y 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 4s^2 4p^6 ; 4d^1 5s^2
1774 case (40)
1775 if (val > 2) call fill_d_orbs(val, 2, 4)
1776 call fill_s_orbs(val, 2, 5) ! Zr 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 4s^2 4p^6 ; 4d^2 5s^2
1777 case (41)
1778 if (val > 1) call fill_d_orbs(val, 4, 4)
1779 call fill_s_orbs(val, 1, 5) ! Nb 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 4s^2 4p^6 ; 4d^4 5s^1
1780 case (42)
1781 if (val > 1) call fill_d_orbs(val, 5, 4)
1782 call fill_s_orbs(val, 1, 5) ! Mo 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 4s^2 4p^6 ; 4d^5 5s^1
1783 case (43)
1784 if (val > 2) call fill_d_orbs(val, 5, 4)
1785 call fill_s_orbs(val, 2, 5) ! Tc 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 4s^2 4p^6 ; 4d^5 5s^2
1786 case (44)
1787 if (val > 1) call fill_d_orbs(val, 7, 4)
1788 call fill_s_orbs(val, 1, 5) ! Ru 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 4s^2 4p^6 ; 4d^7 5s^1
1789 case (45)
1790 if (val > 1) call fill_d_orbs(val, 8, 4)
1791 call fill_s_orbs(val, 1, 5) ! Rh 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 4s^2 4p^6 ; 4d^8 5s^1
1792 case (46)
1793 call fill_d_orbs(val, 10, 4) ! Pd 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 4s^2 4p^6 ; 4d^10
1794 case (47)
1795 call fill_s_orbs(val, 1, 5) ! Ag 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 4s^2 4p^6 4d^10 ; 5s^1
1796 case (48)
1797 call fill_s_orbs(val, 2, 5) ! Cd 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 4s^2 4p^6 4d^10 ; 5s^2
1798 case (49)
1799 call fill_p_orbs(val, 1, 5) ! In 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 4s^2 4p^6 4d^10 ; 5s^2 5p^1
1800 case (50)
1801 call fill_p_orbs(val, 2, 5) ! Sn 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 4s^2 4p^6 4d^10 ; 5s^2 5p^2
1802 case (51)
1803 call fill_p_orbs(val, 3, 5) ! Sb 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 4s^2 4p^6 4d^10 ; 5s^2 5p^3
1804 case (52)
1805 call fill_p_orbs(val, 4, 5) ! Te 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 4s^2 4p^6 4d^10 ; 5s^2 5p^4
1806 case (53)
1807 call fill_p_orbs(val, 5, 5) ! I 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 4s^2 4p^6 4d^10 ; 5s^2 5p^5
1808 case (54)
1809 call fill_p_orbs(val, 6, 5) ! Xe 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 4s^2 4p^6 4d^10 ; 5s^2 5p^6
1810 case (55)
1811 call fill_s_orbs(val, 1, 6) ! Cs 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 4s^2 4p^6 4d^10 5s^2 5p^6 ; 6s^1
1812 case (56)
1813 call fill_s_orbs(val, 2, 6) ! Ba 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 4s^2 4p^6 4d^10 5s^2 5p^6 ; 6s^2
1814 case (57)
1815 if (val > 2) call fill_d_orbs(val, 1, 5)
1816 call fill_s_orbs(val, 2, 6) ! La 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 4s^2 4p^6 4d^10 5s^2 5p^6 ; 5d^1 6s^2
1817 case (58)
1818 if (val > 3) call fill_f_orbs(val, 1, 4)
1819 if (val > 2) call fill_d_orbs(val, 1, 5)
1820 call fill_s_orbs(val, 2, 6) ! Ce 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 4s^2 4p^6 4d^10 5s^2 5p^6 ; 4f^1 5d^1 6s^2
1821 case (59)
1822 if (val > 2) call fill_f_orbs(val, 3, 4)
1823 call fill_s_orbs(val, 2, 6) ! Pr 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 4s^2 4p^6 4d^10 5s^2 5p^6 ; 4f^3 6s^2
1824 case (60)
1825 if (val > 2) call fill_f_orbs(val, 4, 4)
1826 call fill_s_orbs(val, 2, 6) ! Nd 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 4s^2 4p^6 4d^10 5s^2 5p^6 ; 4f^4 6s^2
1827 case (61)
1828 if (val > 2) call fill_f_orbs(val, 5, 4)
1829 call fill_s_orbs(val, 2, 6) ! Pm 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 4s^2 4p^6 4d^10 5s^2 5p^6 ; 4f^5 6s^2
1830 case (62)
1831 if (val > 2) call fill_f_orbs(val, 6, 4)
1832 call fill_s_orbs(val, 2, 6) ! Sm 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 4s^2 4p^6 4d^10 5s^2 5p^6 ; 4f^6 6s^2
1833 case (63)
1834 if (val > 2) call fill_f_orbs(val, 7, 4)
1835 call fill_s_orbs(val, 2, 6) ! Eu 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 4s^2 4p^6 4d^10 5s^2 5p^6 ; 4f^7 6s^2
1836 case (64)
1837 if (val > 3) call fill_f_orbs(val, 7, 4)
1838 if (val > 2) call fill_d_orbs(val, 1, 5)
1839 call fill_s_orbs(val, 2, 6) ! Gd 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 4s^2 4p^6 4d^10 5s^2 5p^6 ; 4f^7 5d^1 6s^2
1840 case (65)
1841 if (val > 2) call fill_f_orbs(val, 9, 4)
1842 call fill_s_orbs(val, 2, 6) ! Tb 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 4s^2 4p^6 4d^10 5s^2 5p^6 ; 4f^9 6s^2
1843 case (66)
1844 if (val > 2) call fill_f_orbs(val, 10, 4)
1845 call fill_s_orbs(val, 2, 6) ! Dy 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 4s^2 4p^6 4d^10 5s^2 5p^6 ; 4f^10 6s^2
1846 case (67)
1847 if (val > 2) call fill_f_orbs(val, 11, 4)
1848 call fill_s_orbs(val, 2, 6) ! Ho 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 4s^2 4p^6 4d^10 5s^2 5p^6 ; 4f^11 6s^2
1849 case (68)
1850 if (val > 2) call fill_f_orbs(val, 12, 4)
1851 call fill_s_orbs(val, 2, 6) ! Er 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 4s^2 4p^6 4d^10 5s^2 5p^6 ; 4f^12 6s^2
1852 case (69)
1853 if (val > 2) call fill_f_orbs(val, 13, 4)
1854 call fill_s_orbs(val, 2, 6) ! Tm 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 4s^2 4p^6 4d^10 5s^2 5p^6 ; 4f^13 6s^2
1855 case (70)
1856 if (val > 2) call fill_f_orbs(val, 14, 4)
1857 call fill_s_orbs(val, 2, 6) ! Yb 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 4s^2 4p^6 4d^10 5s^2 5p^6 ; 4f^14 6s^2
1858 case (71)
1859 if (val > 3) call fill_f_orbs(val, 14, 4)
1860 if (val > 2) call fill_d_orbs(val, 1, 5)
1861 call fill_s_orbs(val, 2, 6) ! Lu 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 4s^2 4p^6 4d^10 5s^2 5p^6 4f^14 ; 5d^1 6s^2
1862 case (72)
1863 if (val > 4) call fill_f_orbs(val, 14, 4)
1864 if (val > 2) call fill_d_orbs(val, 2, 5)
1865 call fill_s_orbs(val, 2, 6) ! Hf 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 4s^2 4p^6 4d^10 5s^2 5p^6 4f^14 ; 5d^2 6s^2
1866 case (73)
1867 if (val > 5) call fill_f_orbs(val, 14, 4)
1868 if (val > 2) call fill_d_orbs(val, 3, 5)
1869 call fill_s_orbs(val, 2, 6) ! Ta 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 4s^2 4p^6 4d^10 5s^2 5p^6 4f^14 ; 5d^3 6s^2
1870 case (74)
1871 if (val > 6) call fill_f_orbs(val, 14, 4)
1872 if (val > 2) call fill_d_orbs(val, 4, 5)
1873 call fill_s_orbs(val, 2, 6) ! W 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 4s^2 4p^6 4d^10 5s^2 5p^6 4f^14 ; 5d^4 6s^2
1874 case (75)
1875 if (val > 7) call fill_f_orbs(val, 14, 4)
1876 if (val > 2) call fill_d_orbs(val, 5, 5)
1877 call fill_s_orbs(val, 2, 6) ! Re 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 4s^2 4p^6 4d^10 5s^2 5p^6 4f^14 ; 5d^5 6s^2
1878 case (76)
1879 if (val > 8) call fill_f_orbs(val, 14, 4)
1880 if (val > 2) call fill_d_orbs(val, 6, 5)
1881 call fill_s_orbs(val, 2, 6) ! Os 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 4s^2 4p^6 4d^10 5s^2 5p^6 4f^14 ; 5d^6 6s^2
1882 case (77)
1883 if (val > 9) call fill_f_orbs(val, 14, 4)
1884 if (val > 2) call fill_d_orbs(val, 7, 5)
1885 call fill_s_orbs(val, 2, 6) ! Ir 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 4s^2 4p^6 4d^10 5s^2 5p^6 4f^14 ; 5d^7 6s^2
1886 case (78)
1887 if (val > 10) call fill_f_orbs(val, 14, 4)
1888 if (val > 1) call fill_d_orbs(val, 9, 5)
1889 call fill_s_orbs(val, 1, 6) ! Pt 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 4s^2 4p^6 4d^10 5s^2 5p^6 4f^14 ; 5d^9 6s^1
1890 case (79)
1891 if (val > 21) call fill_f_orbs(val, 14, 4)
1892 if (val > 1) call fill_d_orbs(val, 10, 5)
1893 call fill_s_orbs(val, 1, 6) ! Au 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 4s^2 4p^6 4d^10 5s^2 5p^6 4f^14 5d^10 ; 6s^1
1894 case (80)
1895 if (val > 12) call fill_f_orbs(val, 14, 4)
1896 if (val > 2) call fill_d_orbs(val, 10, 5)
1897 call fill_s_orbs(val, 2, 6) ! Hg 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 4s^2 4p^6 4d^10 5s^2 5p^6 4f^14 5d^10 ; 6s^2
1898 case (81)
1899 if (val > 10) call fill_d_orbs(val, 10, 5)
1900 if (val > 1) call fill_s_orbs(val, 2, 6)
1901 call fill_p_orbs(val, 1, 6) ! Tl 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 4s^2 4p^6 4d^10 5s^2 5p^6 4f^14 5d^10 ; 6s^2 6p^1
1902 case (82)
1903 if (val > 10) call fill_d_orbs(val, 10, 5)
1904 if (val > 2) call fill_s_orbs(val, 2, 6)
1905 call fill_p_orbs(val, 2, 6) ! Pb 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 4s^2 4p^6 4d^10 5s^2 5p^6 4f^14 5d^10 ; 6s^2 6p^2
1906 case (83)
1907 if (val > 10) call fill_d_orbs(val, 10, 5)
1908 if (val > 3) call fill_s_orbs(val, 2, 6)
1909 call fill_p_orbs(val, 3, 6) ! Bi 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 4s^2 4p^6 4d^10 5s^2 5p^6 4f^14 5d^10 ; 6s^2 6p^3
1910 case (84)
1911 if (val > 10) call fill_d_orbs(val, 10, 5)
1912 if (val > 4) call fill_s_orbs(val, 2, 6)
1913 call fill_p_orbs(val, 4, 6) ! Po 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 4s^2 4p^6 4d^10 5s^2 5p^6 4f^14 5d^10 ; 6s^2 6p^4
1914 case (85)
1915 if (val > 10) call fill_d_orbs(val, 10, 5)
1916 if (val > 5) call fill_s_orbs(val, 2, 6)
1917 call fill_p_orbs(val, 5, 6) ! At 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 4s^2 4p^6 4d^10 5s^2 5p^6 4f^14 5d^10 ; 6s^2 6p^5
1918 case (86)
1919 if (val > 10) call fill_d_orbs(val, 10, 5)
1920 if (val > 6) call fill_s_orbs(val, 2, 6)
1921 call fill_p_orbs(val, 6, 6) ! Rn 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 4s^2 4p^6 4d^10 5s^2 5p^6 4f^14 5d^10 ; 6s^2 6p^6
1922 case (87)
1923 if (val > 1) call fill_p_orbs(val, 6, 6)
1924 call fill_s_orbs(val, 1, 7) ! Fr 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 4s^2 4p^6 4d^10 5s^2 5p^6 4f^14 5d^10 6s^2 6p^6 ; 7s1
1925 case (88)
1926 if (val > 2) call fill_p_orbs(val, 6, 6)
1927 call fill_s_orbs(val, 2, 7) ! Ra 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 4s^2 4p^6 4d^10 5s^2 5p^6 4f^14 5d^10 6s^2 6p^6 ; 7s1
1928 case (89)
1929 if (val > 3) call fill_p_orbs(val, 6, 6)
1930 if (val > 2) call fill_d_orbs(val, 1, 6)
1931 call fill_s_orbs(val, 2, 7) ! Ac 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 4s^2 4p^6 4d^10 5s^2 5p^6 4f^14 5d^10 6s^2 6p^6 ; 6d1 7s2
1932 case (90)
1933 if (val > 4) call fill_p_orbs(val, 6, 6)
1934 if (val > 2) call fill_d_orbs(val, 2, 6)
1935 call fill_s_orbs(val, 2, 7) ! Th 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 4s^2 4p^6 4d^10 5s^2 5p^6 4f^14 5d^10 6s^2 6p^6 ; 6d2 7s2
1936 case (91)
1937 if (val > 5) call fill_p_orbs(val, 6, 6)
1938 if (val > 3) call fill_f_orbs(val, 2, 5)
1939 if (val > 2) call fill_d_orbs(val, 1, 6)
1940 call fill_s_orbs(val, 2, 7) ! Pa 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 4s^2 4p^6 4d^10 5s^2 5p^6 4f^14 5d^10 6s^2 6p^6 ; 5f2 6d1 7s2
1941 case (92)
1942 if (val > 6) call fill_p_orbs(val, 6, 6)
1943 if (val > 3) call fill_f_orbs(val, 3, 5)
1944 if (val > 2) call fill_d_orbs(val, 1, 6)
1945 call fill_s_orbs(val, 2, 7) ! U 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 4s^2 4p^6 4d^10 5s^2 5p^6 4f^14 5d^10 6s^2 6p^6 ; 5f3 6d1 7s2
1946 case (93)
1947 if (val > 7) call fill_p_orbs(val, 6, 6)
1948 if (val > 3) call fill_f_orbs(val, 4, 5)
1949 if (val > 2) call fill_d_orbs(val, 1, 6)
1950 call fill_s_orbs(val, 2, 7) ! Np 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 4s^2 4p^6 4d^10 5s^2 5p^6 4f^14 5d^10 6s^2 6p^6 ; 5f4 6d1 7s2
1951 case (94)
1952 if (val > 8) call fill_p_orbs(val, 6, 6)
1953 if (val > 2) call fill_f_orbs(val, 6, 5)
1954 call fill_s_orbs(val, 2, 7) ! Pu 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 4s^2 4p^6 4d^10 5s^2 5p^6 4f^14 5d^10 6s^2 6p^6 ; 5f6 7s2
1955 case (95)
1956 if (val > 9) call fill_p_orbs(val, 6, 6)
1957 if (val > 2) call fill_f_orbs(val, 7, 5)
1958 call fill_s_orbs(val, 2, 7) ! Am 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 4s^2 4p^6 4d^10 5s^2 5p^6 4f^14 5d^10 6s^2 6p^6 ; 5f7 7s2
1959 case (96)
1960 if (val > 10) call fill_p_orbs(val, 6, 6)
1961 if (val > 3) call fill_f_orbs(val, 7, 5)
1962 if (val > 2) call fill_d_orbs(val, 1, 6)
1963 call fill_s_orbs(val, 2, 7) ! Cm 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 4s^2 4p^6 4d^10 5s^2 5p^6 4f^14 5d^10 6s^2 6p^6 ; 5f7 6d1 7s2
1964 case (97)
1965 if (val > 11) call fill_p_orbs(val, 6, 6)
1966 if (val > 2) call fill_f_orbs(val, 9, 5)
1967 call fill_s_orbs(val, 2, 7) ! Bk 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 4s^2 4p^6 4d^10 5s^2 5p^6 4f^14 5d^10 6s^2 6p^6 ; 5f9 7s2
1968 case (98)
1969 if (val > 12) call fill_p_orbs(val, 6, 6)
1970 if (val > 2) call fill_f_orbs(val, 10, 5)
1971 call fill_s_orbs(val, 2, 7) ! Cf 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 4s^2 4p^6 4d^10 5s^2 5p^6 4f^14 5d^10 6s^2 6p^6 ; 5f10 7s2
1972 case (99)
1973 if (val > 13) call fill_p_orbs(val, 6, 6)
1974 if (val > 2) call fill_f_orbs(val, 11, 5)
1975 call fill_s_orbs(val, 2, 7) ! Es 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 4s^2 4p^6 4d^10 5s^2 5p^6 4f^14 5d^10 6s^2 6p^6 ; 5f11 7s2
1976 case (100)
1977 if (val > 14) call fill_p_orbs(val, 6, 6)
1978 if (val > 2) call fill_f_orbs(val, 12, 5)
1979 call fill_s_orbs(val, 2, 7) ! Fm 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 4s^2 4p^6 4d^10 5s^2 5p^6 4f^14 5d^10 6s^2 6p^6 ; 5f12 7s2
1980 case (101)
1981 if (val > 15) call fill_p_orbs(val, 6, 6)
1982 if (val > 2) call fill_f_orbs(val, 13, 5)
1983 call fill_s_orbs(val, 2, 7) ! Md 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 4s^2 4p^6 4d^10 5s^2 5p^6 4f^14 5d^10 6s^2 6p^6 ; 5f13 7s2
1984 case (102)
1985 if (val > 16) call fill_p_orbs(val, 6, 6)
1986 if (val > 2) call fill_f_orbs(val, 14, 5)
1987 call fill_s_orbs(val, 2, 7) ! No 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 4s^2 4p^6 4d^10 5s^2 5p^6 4f^14 5d^10 6s^2 6p^6 ; 5f14 7s2
1988 case (103)
1989 if (val > 3) call fill_f_orbs(val, 14, 5)
1990 if (val > 1) call fill_s_orbs(val, 2, 7)
1991 call fill_p_orbs(val, 1, 7) ! Lr 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 4s^2 4p^6 4d^10 5s^2 5p^6 4f^14 5d^10 6s^2 6p^6 ; 5f14 7s2 7p1
1992 end select
1993
1994 !If we attributed all the electrons, everything went fine
1995 if (val < m_epsilon) then
1996 !In case of spin-polarized calculations, we properly distribute the electrons
1997 if (ispin == 2) then
1998 call valconf_unpolarized_to_polarized(conf)
1999 end if
2000 else
2001 conf%occ = m_zero
2002 message(1) = "Error in attributing atomic occupations"
2003 call messages_warning(1, namespace=namespace)
2004 end if
2005
2007
2008 contains
2009 subroutine fill_s_orbs(val, max_occ, nn)
2010 real(real64), intent(inout) :: val
2011 integer, intent(in) :: max_occ, nn
2012
2013 conf%p = conf%p + 1
2014 conf%occ(conf%p,1) = min(val, real(max_occ, real64) )
2015 val = val - conf%occ(conf%p,1)
2016 conf%l(conf%p) = 0
2017 conf%n(conf%p) = nn
2018 print *, 'Filling ', conf%n(conf%p), conf%l(conf%p), conf%occ(conf%p,1)
2019 end subroutine fill_s_orbs
2020
2021 subroutine fill_p_orbs(val, max_occ, nn)
2022 real(real64), intent(inout) :: val
2023 integer, intent(in) :: max_occ, nn
2024
2025 conf%p = conf%p + 1
2026 conf%occ(conf%p,1) = min(val, real(max_occ, real64) )
2027 val = val - conf%occ(conf%p,1)
2028 conf%l(conf%p) = 1
2029 conf%n(conf%p) = nn
2030 end subroutine fill_p_orbs
2031
2032 subroutine fill_d_orbs(val, max_occ, nn)
2033 real(real64), intent(inout) :: val
2034 integer, intent(in) :: max_occ, nn
2035
2036 conf%p = conf%p + 1
2037 conf%occ(conf%p,1) = min(val, real(max_occ, real64) )
2038 val = val - conf%occ(conf%p,1)
2039 conf%l(conf%p) = 2
2040 conf%n(conf%p) = nn
2041 end subroutine fill_d_orbs
2042
2043 subroutine fill_f_orbs(val, max_occ, nn)
2044 real(real64), intent(inout) :: val
2045 integer, intent(in) :: max_occ, nn
2046
2047 conf%p = conf%p + 1
2048 conf%occ(conf%p,1) = min(val, real(max_occ, real64) )
2049 val = val - conf%occ(conf%p,1)
2050 conf%l(conf%p) = 3
2051 conf%n(conf%p) = nn
2052 end subroutine fill_f_orbs
2053
2054 end subroutine ps_guess_atomic_occupations
2055
2056#include "ps_pspio_inc.F90"
2057
2058end module ps_oct_m
2059
2060!! Local Variables:
2061!! mode: f90
2062!! coding: utf-8
2063!! End:
Some operations may be done for one spline-function, or for an array of them.
Definition: splines.F90:181
double exp(double __x) __attribute__((__nothrow__
double sqrt(double __x) __attribute__((__nothrow__
subroutine, public valconf_unpolarized_to_polarized(conf)
Definition: atomic.F90:264
subroutine, public valconf_copy(cout, cin)
Definition: atomic.F90:170
real(real64), parameter, public m_two
Definition: global.F90:193
real(real64), parameter, public m_zero
Definition: global.F90:191
real(real64), parameter, public m_pi
some mathematical constants
Definition: global.F90:189
real(real64), parameter, public m_min_exp_arg
Definition: global.F90:210
real(real64), parameter, public m_half
Definition: global.F90:197
real(real64), parameter, public m_one
Definition: global.F90:192
real(real64), parameter, public m_three
Definition: global.F90:194
Definition: io.F90:116
subroutine, public logrid_copy(grid_in, grid_out)
Definition: logrid.F90:230
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:117
subroutine, public messages_obsolete_variable(namespace, name, rep)
Definition: messages.F90:1023
subroutine, public messages_new_line()
Definition: messages.F90:1112
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:162
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:410
subroutine, public messages_input_error(namespace, var, details, row, column)
Definition: messages.F90:691
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:594
subroutine, public ps_cpi_end(ps_cpi)
Definition: ps_cpi.F90:185
subroutine, public ps_cpi_init(ps_cpi, filename, namespace)
Definition: ps_cpi.F90:150
subroutine, public ps_cpi_process(ps_cpi, lloc, namespace)
Definition: ps_cpi.F90:223
subroutine, public ps_fhi_end(ps_fhi)
Definition: ps_fhi.F90:185
subroutine, public ps_fhi_init(ps_fhi, filename, namespace)
Definition: ps_fhi.F90:151
subroutine, public ps_fhi_process(ps_fhi, lmax, lloc, namespace)
Definition: ps_fhi.F90:204
subroutine, public hgh_end(psp)
Definition: ps_hgh.F90:230
subroutine, public hgh_get_eigen(psp, eigen)
Definition: ps_hgh.F90:291
subroutine, public hgh_process(psp, namespace)
Definition: ps_hgh.F90:249
subroutine, public hgh_init(psp, filename, namespace)
Definition: ps_hgh.F90:185
Definition: ps.F90:116
pure logical function, public ps_has_density(ps)
Definition: ps.F90:1634
integer, parameter, public ps_filter_ts
Definition: ps.F90:166
subroutine, public ps_end(ps)
Definition: ps.F90:1065
subroutine hgh_load(ps, ps_hgh)
Definition: ps.F90:1104
subroutine, public ps_init(ps, namespace, label, z, user_lmax, user_llocal, ispin, filename)
Definition: ps.F90:262
subroutine ps_info(ps, filename, namespace)
Definition: ps.F90:581
subroutine, public ps_separate(ps)
Separate the local potential into (soft) long-ranged and (hard) short-ranged parts.
Definition: ps.F90:692
subroutine ps_check_bound(ps, eigen)
Definition: ps.F90:889
subroutine ps_grid_load(ps, ps_grid)
Definition: ps.F90:1175
subroutine, public ps_debug(ps, dir, namespace, gmax)
Definition: ps.F90:934
real(real64) pure function, public long_range_potential(r, sigma, z_val)
Evaluate the long-range potential at a given distance.
Definition: ps.F90:750
integer, parameter, public proj_hgh
Definition: ps.F90:171
pure integer function, public ps_bound_niwfs(ps)
Returns the number of bound atomic orbitals taking into account then m quantum number multiplicity.
Definition: ps.F90:1618
real(real64), parameter eps
Cutoff for determining the radius of the NLCC.
Definition: ps.F90:255
subroutine, public ps_getradius(ps)
Definition: ps.F90:766
subroutine, public ps_guess_atomic_occupations(namespace, zz, valcharge, ispin, conf)
This routines provides, given Z and the number of valence electron the occupations of the orbitals....
Definition: ps.F90:1690
pure logical function, public ps_has_nlcc(ps)
Definition: ps.F90:1643
pure integer function, public ps_niwfs(ps)
Returns the number of atomic orbitals taking into account then m quantum number multiplicity.
Definition: ps.F90:1603
subroutine, public ps_pspio_init(ps, namespace, label, z, lmax, ispin, filename)
Definition: ps.F90:2172
integer, parameter, public ps_filter_none
Definition: ps.F90:166
integer, parameter, public proj_j_average
Fully-relativistic pseudopotentials with separate j-average and SOC terms.
Definition: ps.F90:177
integer, parameter, public proj_rkb
Definition: ps.F90:171
integer, parameter, public proj_j_dependent
Fully-relativistic j-dependent pseudopotentials.
Definition: ps.F90:177
subroutine, public ps_derivatives(ps)
Definition: ps.F90:786
subroutine, public ps_filter(ps, filter, gmax)
Definition: ps.F90:803
real(real64) function, public ps_density_volume(ps, namespace)
Definition: ps.F90:1651
character(len=4), dimension(pseudo_format_upf1:pseudo_format_psp8), parameter ps_name
Definition: ps.F90:184
integer, parameter, public ps_filter_bsb
Definition: ps.F90:166
integer, parameter, public proj_kb
Definition: ps.F90:171
subroutine ps_xml_load(ps, ps_xml, namespace)
Loads XML files for QSO, UPF1, UPF2, PSML, and PSP8 formats.
Definition: ps.F90:1277
subroutine, public ps_psf_init(pstm, ispin, filename, namespace)
Definition: ps_psf.F90:153
subroutine, public ps_psf_process(ps_psf, namespace, lmax, lloc)
Definition: ps_psf.F90:303
subroutine, public ps_psf_get_eigen(ps_psf, eigen)
Definition: ps_psf.F90:337
subroutine, public ps_psf_end(ps_psf)
Definition: ps_psf.F90:204
subroutine, public ps_xml_end(this)
Definition: ps_xml.F90:340
subroutine, public ps_xml_init(this, namespace, filename, fmt, ierr)
Definition: ps_xml.F90:166
integer, parameter, public pseudo_format_file_not_found
Definition: pseudo.F90:175
integer, parameter, public pseudo_type_kleinman_bylander
Definition: pseudo.F90:160
integer, parameter, public pseudo_format_qso
Definition: pseudo.F90:175
integer, parameter, public pseudo_format_hgh
Definition: pseudo.F90:175
integer, parameter, public pseudo_format_upf2
Definition: pseudo.F90:175
integer, parameter, public pseudo_type_paw
Definition: pseudo.F90:160
integer, parameter, public pseudo_format_cpi
Definition: pseudo.F90:175
integer, parameter, public pseudo_exchange_unknown
Definition: pseudo.F90:190
integer, parameter, public pseudo_type_semilocal
Definition: pseudo.F90:160
integer, parameter, public pseudo_type_ultrasoft
Definition: pseudo.F90:160
integer, parameter, public pseudo_correlation_unknown
Definition: pseudo.F90:194
integer, parameter, public pseudo_format_psf
Definition: pseudo.F90:175
integer, parameter, public pseudo_format_fhi
Definition: pseudo.F90:175
integer, parameter, public pseudo_format_unknown
Definition: pseudo.F90:175
integer, parameter, public pseudo_format_psml
Definition: pseudo.F90:175
subroutine, public spline_der(spl, dspl, threshold)
Definition: splines.F90:923
subroutine, public spline_fit(nrc, rofi, ffit, spl, threshold)
Definition: splines.F90:415
real(real64) function, public spline_eval(spl, x)
Definition: splines.F90:443
character(kind=c_char, len=1) function, dimension(:), allocatable, public string_f_to_c(f_string)
convert a Fortran string to a C string
Definition: string.F90:273
subroutine fill_s_orbs(val, max_occ, nn)
Definition: ps.F90:2105
subroutine fill_p_orbs(val, max_occ, nn)
Definition: ps.F90:2117
subroutine fill_d_orbs(val, max_occ, nn)
Definition: ps.F90:2128
subroutine get_splines()
Definition: ps.F90:1128
subroutine fill_f_orbs(val, max_occ, nn)
Definition: ps.F90:2139
remember that the FHI format is basically the CPI format with a header
Definition: ps_fhi.F90:139
The following data type contains: (a) the pseudopotential parameters, as read from a *....
Definition: ps_hgh.F90:146
A type storing the information and data about a pseudopotential.
Definition: ps.F90:188
int true(void)