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