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 implicit none
50
51 private
52 public :: &
53 ps_t, &
54 ps_init, &
57 ps_filter, &
60 ps_debug, &
61 ps_niwfs, &
63 ps_end, &
68
69 integer, parameter, public :: &
70 PS_FILTER_NONE = 0, &
71 ps_filter_ts = 2, &
73
74 integer, public, parameter :: &
75 PROJ_NONE = 0, &
76 proj_hgh = 1, &
77 proj_kb = 2, &
78 proj_rkb = 3
79
80 integer, public, parameter :: &
81 PROJ_J_INDEPENDENT = 0, & !< Non-relativisitic or scalar-relativistic pseudopotentials
82 proj_j_dependent = 1, &
84
85 integer, parameter, public :: INVALID_L = 333
86
87 character(len=4), parameter :: ps_name(PSEUDO_FORMAT_UPF1:PSEUDO_FORMAT_PSP8) = &
88 (/"upf1", "upf2", "qso ", "psml", "psf ", "cpi ", "fhi ", "hgh ", "psp8"/)
89
91 type ps_t
92 ! Components are public by default
93 integer :: projector_type
94 integer :: relativistic_treatment
96
97 character(len=10), private :: label
98
99 integer, private :: ispin
100
101 real(real64), private :: z
102 real(real64) :: z_val
103 type(valconf_t) :: conf
104
105 type(logrid_t), private :: g
106
107 type(spline_t), allocatable :: ur(:, :)
108 logical, allocatable :: bound(:, :)
109
110 ! Kleinman-Bylander projectors stuff
111 integer :: lmax
112 integer :: llocal
113
114 type(spline_t) :: vl
115 logical :: no_vl = .false.
116
117 real(real64) :: projectors_sphere_threshold
124 real(real64) :: rc_max
125
126 integer :: kbc
127 integer :: projectors_per_l(5)
128 real(real64), allocatable :: h(:,:,:), k(:, :, :)
129 type(spline_t), allocatable :: kb(:, :)
130 type(spline_t), allocatable :: dkb(:, :)
131
132 logical :: nlcc
133 type(spline_t) :: core
134 type(spline_t) :: core_der
135
136
137 !LONG-RANGE PART OF THE LOCAL POTENTIAL
138
139 logical, private :: has_long_range
140 logical, private :: is_separated
141
142 type(spline_t), private :: vlr
143 type(spline_t) :: vlr_sq
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
188
189 ! Fix the threshold to calculate the radius of the projector-function localization spheres:
191 call messages_obsolete_variable(namespace, 'SpecieProjectorSphereThreshold', 'SpeciesProjectorSphereThreshold')
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')
212
213 ps%file_format = pseudo_detect_format(filename)
214
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
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
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)
244 ps%pseudo_type = pseudo_type_semilocal
245
246 call ps_psf_init(ps_psf, ispin, filename, namespace)
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
252
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)
580
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, 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 arg = r/(ps%sigma_erf*sqrt(m_two))
618 ! Note that the threshold is taken such that the Taylor expansion produces an error smaller than \epsilon
619 ! given a third order Taylor expansion
620 if(arg > 1e-3_real64) then
621 vlr(ii) = -ps%z_val*loct_erf(arg)/r
622 else
623 vlr(ii) = -ps%z_val*m_two/(ps%sigma_erf*sqrt(m_two*m_pi)) * (m_one - arg**2 / m_three)
624 end if
625
626 vsr(ii) = spline_eval(ps%vl, r) - vlr(ii)
627 if(abs(vsr(ii)) < 1.0e-14_real64) vsr(ii) = m_zero
628 max_val_vsr = max(max_val_vsr, abs(vsr(ii)))
629
630 exp_arg = -m_half*r**2/ps%sigma_erf**2
631 if (exp_arg > m_min_exp_arg) then
632 nlr(ii) = -ps%z_val/(ps%sigma_erf*sqrt(m_two*m_pi))**3*exp(exp_arg)
633 else
634 nlr(ii) = m_zero
635 end if
636 end do
637
638 call spline_init(ps%vlr)
639 call spline_fit(ps%g%nrval, ps%g%rofi, vlr, ps%vlr, ps%projectors_sphere_threshold)
640
641 call spline_init(ps%vlr_sq)
642 call spline_fit(ps%g%nrval, ps%g%r2ofi, vlr, ps%vlr_sq, ps%projectors_sphere_threshold)
643
644 call spline_init(ps%nlr)
645 call spline_fit(ps%g%nrval, ps%g%rofi, nlr, ps%nlr, ps%projectors_sphere_threshold)
646
647 !overwrite vl
648 call spline_end(ps%vl)
649 call spline_init(ps%vl)
650 call spline_fit(ps%g%nrval, ps%g%rofi, vsr, ps%vl, ps%projectors_sphere_threshold)
651 if (max_val_vsr < 1.0e-12_real64) ps%no_vl = .true.
652
653 safe_deallocate_a(vsr)
654 safe_deallocate_a(vlr)
655 safe_deallocate_a(nlr)
656
657 ps%is_separated = .true.
658
659 pop_sub(ps_separate)
660 end subroutine ps_separate
661
662
663 ! ---------------------------------------------------------
664 subroutine ps_getradius(ps)
665 type(ps_t), intent(inout) :: ps
666 integer :: l, j
667
668 push_sub(ps_getradius)
669
670 ps%rc_max = m_zero
671
672 do l = 0, ps%lmax
673 if (l == ps%llocal) cycle
674 do j = 1, ps%kbc
675 ps%rc_max = max(ps%rc_max, ps%kb(l, j)%x_threshold)
676 end do
677 end do
678
679 pop_sub(ps_getradius)
680 end subroutine ps_getradius
681
682
683 ! ---------------------------------------------------------
684 subroutine ps_derivatives(ps)
685 type(ps_t), intent(inout) :: ps
686 integer :: l, j
687
688 push_sub(ps_derivatives)
690 do l = 0, ps%lmax
691 do j = 1, ps%kbc
692 call spline_der(ps%kb(l, j), ps%dkb(l, j), ps%projectors_sphere_threshold)
693 end do
694 end do
695
696 pop_sub(ps_derivatives)
697 end subroutine ps_derivatives
698
699
700 ! ---------------------------------------------------------
701 subroutine ps_filter(ps, filter, gmax)
702 type(ps_t), intent(inout) :: ps
703 integer, intent(in) :: filter
704 real(real64), intent(in) :: gmax
705
706 integer :: l, k, ispin
707
708 real(real64) :: alpha, beta_fs, rmax, rcut, gamma, beta_rs
709
710 push_sub(ps_filter)
711 call profiling_in("PS_FILTER")
712
713 select case (filter)
714 case (ps_filter_none)
715
716 case (ps_filter_ts)
717 alpha = 1.1_real64
718 gamma = m_two
719
720 if(.not. ps%no_vl) then
721 rmax = ps%vl%x_threshold
722 ! For the special case l=-1, we use l=0 when doing the filtering, but it is not clear this is correct
723 call spline_filter_mask(ps%vl, max(0, ps%llocal), rmax, gmax, alpha, gamma, ps%projectors_sphere_threshold)
724 end if
725
726 do l = 0, ps%lmax
727 if (l == ps%llocal) cycle
728 do k = 1, ps%kbc
729 call spline_filter_mask(ps%kb(l, k), l, ps%rc_max, gmax, alpha, gamma, ps%projectors_sphere_threshold)
730 end do
731 end do
732
733 if (ps_has_nlcc(ps)) then
734 rmax = ps%core%x_threshold
735 call spline_filter_mask(ps%core, 0, rmax, gmax, alpha, gamma, ps%projectors_sphere_threshold)
736 end if
737
738 if (ps_has_density(ps)) then
739 do ispin = 1, ps%ispin
740 if (abs(spline_integral(ps%density(ispin))) > 1.0e-12_real64) then
741 rmax = ps%density(ispin)%x_threshold
742 call spline_filter_mask(ps%density(ispin), 0, rmax, gmax, alpha, gamma, ps%projectors_sphere_threshold)
743 call spline_force_pos(ps%density(ispin), ps%projectors_sphere_threshold)
744 end if
745
746 if (abs(spline_integral(ps%density_der(ispin))) > 1.0e-12_real64) then
747 rmax = ps%density_der(ispin)%x_threshold
748 call spline_filter_mask(ps%density_der(ispin), 0, rmax, gmax, alpha, gamma, ps%projectors_sphere_threshold)
749 end if
750 end do
751 end if
752
753 case (ps_filter_bsb)
754 alpha = 0.7_real64 ! The original was M_FOUR/7.0_real64
755 beta_fs = 18.0_real64
756 rcut = 2.5_real64
757 beta_rs = 0.4_real64
758
759 ! For the special case l=-1, we use l=0 when doing the filtering, but it is not clear this is correct
760 call spline_filter_bessel(ps%vl, max(0, ps%llocal), gmax, alpha, beta_fs, rcut, beta_rs, ps%projectors_sphere_threshold)
761 do l = 0, ps%lmax
762 if (l == ps%llocal) cycle
763 do k = 1, ps%kbc
764 call spline_filter_bessel(ps%kb(l, k), l, gmax, alpha, beta_fs, rcut, beta_rs, ps%projectors_sphere_threshold)
765 end do
766 end do
767
768 if (ps_has_nlcc(ps)) then
769 call spline_filter_bessel(ps%core, 0, gmax, alpha, beta_fs, rcut, beta_rs, ps%projectors_sphere_threshold)
770 end if
771
772 if (ps_has_density(ps)) then
773 do ispin = 1, ps%ispin
774 call spline_filter_bessel(ps%density(ispin), 0, gmax, alpha, beta_fs, rcut, beta_rs, ps%projectors_sphere_threshold)
775 call spline_force_pos(ps%density(ispin), ps%projectors_sphere_threshold)
776 call spline_filter_bessel(ps%density_der(ispin), 0, gmax, alpha, beta_fs, rcut, beta_rs, ps%projectors_sphere_threshold)
777 end do
778 end if
779
780 end select
781
782 call profiling_out("PS_FILTER")
783 pop_sub(ps_filter)
784 end subroutine ps_filter
785
786 ! ---------------------------------------------------------
787 subroutine ps_check_bound(ps, eigen)
788 type(ps_t), intent(inout) :: ps
789 real(real64), intent(in) :: eigen(:,:)
790
791 integer :: i, is, ir
792 real(real64) :: ur1, ur2
793
795
796 ! Unbound states have positive eigenvalues
797 where(eigen > m_zero)
798 ps%bound = .false.
799 elsewhere
800 ps%bound = .true.
801 end where
802
803 ! We might not have information about the eigenvalues, so we need to check the wavefunctions
804 do i = 1, ps%conf%p
805 do is = 1, ps%ispin
806 if (.not. ps%bound(i, is)) cycle
807
808 do ir = ps%g%nrval, 3, -1
809 ! First we look for the outmost value that is not zero
810 if (abs(spline_eval(ps%ur(i, is), ps%g%rofi(ir))*ps%g%rofi(ir)) > m_zero) then
811 ! Usually bound states have exponentially decaying wavefunctions,
812 ! while unbound states have exponentially diverging
813 ! wavefunctions. Therefore we check if the wavefunctions
814 ! value is increasing with increasing radius. The fact
815 ! that we do not use the wavefunctions outmost value that
816 ! is not zero is on purpose, as some pseudopotential
817 ! generators do funny things with that point.
818 ur1 = spline_eval(ps%ur(i, is), ps%g%rofi(ir-2))*ps%g%rofi(ir-2)
819 ur2 = spline_eval(ps%ur(i, is), ps%g%rofi(ir-1))*ps%g%rofi(ir-1)
820 if ((ur1*ur2 > m_zero) .and. (abs(ur2) > abs(ur1))) ps%bound(i, is) = .false.
821 exit
822 end if
823 end do
824 end do
825 end do
826
827 pop_sub(ps_check_bound)
828 end subroutine ps_check_bound
829
830
831 ! ---------------------------------------------------------
832 subroutine ps_debug(ps, dir, namespace, gmax)
833 type(ps_t), intent(in) :: ps
834 character(len=*), intent(in) :: dir
835 type(namespace_t), intent(in) :: namespace
836 real(real64), intent(in) :: gmax
837
838 ! We will plot also some Fourier transforms.
839 type(spline_t), allocatable :: fw(:, :)
840
841 integer :: iunit
842 integer :: j, k, l
843
844 push_sub(ps_debug)
845
846 ! A text file with some basic data.
847 iunit = io_open(trim(dir)//'/pseudo-info', namespace, action='write')
848 write(iunit,'(a,/)') ps%label
849 write(iunit,'(a,a,/)') 'Format : ', ps_name(ps%file_format)
850 write(iunit,'(a,f6.3)') 'z : ', ps%z
851 write(iunit,'(a,f6.3,/)') 'zval : ', ps%z_val
852 write(iunit,'(a,i4)') 'lmax : ', ps%lmax
853 write(iunit,'(a,i4)') 'lloc : ', ps%llocal
854 write(iunit,'(a,i4,/)') 'kbc : ', ps%kbc
855 write(iunit,'(a,f9.5,/)') 'rcmax : ', ps%rc_max
856 write(iunit,'(a,/)') 'h matrix:'
857 do l = 0, ps%lmax
858 do k = 1, ps%kbc
859 write(iunit,'(10f9.5)') (ps%h(l, k, j), j = 1, ps%kbc)
860 end do
861 end do
862 if (allocated(ps%k)) then
863 write(iunit,'(/,a,/)') 'k matrix:'
864 do l = 0, ps%lmax
865 do k = 1, ps%kbc
866 write(iunit,'(10f9.5)') (ps%k(l, k, j), j = 1, ps%kbc)
867 end do
868 end do
869 end if
870
871 write(iunit,'(/,a)') 'orbitals:'
872 do j = 1, ps%conf%p
873 if (ps%ispin == 2) then
874 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), &
875 'j = ', ps%conf%j(j), 'bound = ', all(ps%bound(j,:)), &
876 'occ(1) = ', ps%conf%occ(j, 1), 'occ(2) = ', ps%conf%occ(j, 2)
877 else
878 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), &
879 'j = ', ps%conf%j(j), 'bound = ', all(ps%bound(j,:)), &
880 'occ = ', ps%conf%occ(j, 1)
881 end if
882 end do
883
884
885 call io_close(iunit)
886
887 ! Local part of the pseudopotential
888 iunit = io_open(trim(dir)//'/local', namespace, action='write')
889 call spline_print(ps%vl, iunit)
890 call io_close(iunit)
891
892 ! Local part of the pseudopotential
893 iunit = io_open(trim(dir)//'/local_long_range', namespace, action='write')
894 call spline_print(ps%vlr, iunit)
895 call io_close(iunit)
896
897 ! Local part of the pseudopotential
898 iunit = io_open(trim(dir)//'/local_long_range_density', namespace, action='write')
899 call spline_print(ps%nlr, iunit)
900 call io_close(iunit)
901
902 ! Fourier transform of the local part
903 iunit = io_open(trim(dir)//'/local_ft', namespace, action='write')
904 safe_allocate(fw(1, 1))
905 call spline_init(fw(1, 1))
906 call spline_3dft(ps%vl, fw(1, 1), ps%projectors_sphere_threshold, gmax = gmax)
907 call spline_print(fw(1, 1), iunit)
908 call spline_end(fw(1, 1))
909 safe_deallocate_a(fw)
910 call io_close(iunit)
911
912 ! Kleinman-Bylander projectors
913 iunit = io_open(trim(dir)//'/nonlocal', namespace, action='write')
914 call spline_print(ps%kb, iunit)
915 call io_close(iunit)
916
917 iunit = io_open(trim(dir)//'/nonlocal_derivative', namespace, action='write')
918 call spline_print(ps%dkb, iunit)
919 call io_close(iunit)
920
921 iunit = io_open(trim(dir)//'/nonlocal_ft', namespace, action='write')
922 safe_allocate(fw(0:ps%lmax, 1:ps%kbc))
923 call spline_init(fw)
924 do k = 0, ps%lmax
925 do j = 1, ps%kbc
926 call spline_besselft(ps%kb(k, j), fw(k, j), k, threshold=ps%projectors_sphere_threshold, &
927 gmax=m_four*gmax)
928 end do
929 end do
930 call spline_print(fw, iunit)
931 call spline_end(fw)
932 safe_deallocate_a(fw)
933 call io_close(iunit)
934
935 ! Pseudo-wavefunctions
936 iunit = io_open(trim(dir)//'/wavefunctions', namespace, action='write')
937 call spline_print(ps%ur, iunit)
938 call io_close(iunit)
939
940 ! Density
941 if (ps%has_density) then
942 iunit = io_open(trim(dir)//'/density', namespace, action='write')
943 call spline_print(ps%density, iunit)
944 call io_close(iunit)
945
946 iunit = io_open(trim(dir)//'/density_derivative', namespace, action='write')
947 call spline_print(ps%density_der, iunit)
948 call io_close(iunit)
949 end if
950
951 ! Non-linear core-corrections
952 if (ps_has_nlcc(ps)) then
953 iunit = io_open(trim(dir)//'/nlcc', namespace, action='write')
954 call spline_print(ps%core, iunit)
955 call io_close(iunit)
956 end if
957
958 pop_sub(ps_debug)
959 end subroutine ps_debug
960
961
962 ! ---------------------------------------------------------
963 subroutine ps_end(ps)
964 type(ps_t), intent(inout) :: ps
965
966 if (.not. allocated(ps%kb)) return
967
968 push_sub(ps_end)
969
970 if (ps%is_separated) then
971 call spline_end(ps%vlr)
972 call spline_end(ps%vlr_sq)
973 call spline_end(ps%nlr)
974 end if
975
976 call spline_end(ps%kb)
977 call spline_end(ps%dkb)
978 call spline_end(ps%ur)
979
980 call spline_end(ps%vl)
981 call spline_end(ps%core)
982 call spline_end(ps%core_der)
983
984 if (allocated(ps%density)) call spline_end(ps%density)
985 if (allocated(ps%density_der)) call spline_end(ps%density_der)
986
987 call logrid_end(ps%g)
988
989 safe_deallocate_a(ps%kb)
990 safe_deallocate_a(ps%dkb)
991 safe_deallocate_a(ps%ur)
992 safe_deallocate_a(ps%bound)
993 safe_deallocate_a(ps%h)
994 safe_deallocate_a(ps%k)
995 safe_deallocate_a(ps%density)
996 safe_deallocate_a(ps%density_der)
997
998 pop_sub(ps_end)
999 end subroutine ps_end
1000
1001
1002 ! ---------------------------------------------------------
1003 subroutine hgh_load(ps, ps_hgh)
1004 type(ps_t), intent(inout) :: ps
1005 type(ps_hgh_t), intent(inout) :: ps_hgh
1006
1007 push_sub(hgh_load)
1008
1009 ! Fixes some components of ps
1010 ps%nlcc = .false.
1011 if (ps%lmax >= 0) then
1012 ps%rc_max = 1.1_real64 * maxval(ps_hgh%kbr(0:ps%lmax)) ! Increase a little.
1013 else
1014 ps%rc_max = m_zero
1015 end if
1016 ps%h(0:ps%lmax, 1:ps%kbc, 1:ps%kbc) = ps_hgh%h(0:ps%lmax, 1:ps%kbc, 1:ps%kbc)
1017 ps%k(0:ps%lmax, 1:ps%kbc, 1:ps%kbc) = ps_hgh%k(0:ps%lmax, 1:ps%kbc, 1:ps%kbc)
1018
1019 ! now we fit the splines
1020 call get_splines()
1021
1022 pop_sub(hgh_load)
1023
1024 contains
1025
1026 ! ---------------------------------------------------------
1027 subroutine get_splines()
1028 integer :: l, is, j, ip
1029 real(real64), allocatable :: hato(:), dens(:)
1030
1031 push_sub(hgh_load.get_splines)
1032
1033 safe_allocate(hato(1:ps_hgh%g%nrval))
1034 safe_allocate(dens(1:ps_hgh%g%nrval))
1035
1036 ! Interpolate the KB-projection functions
1037 do l = 0, ps_hgh%l_max
1038 do j = 1, 3
1039 hato(:) = ps_hgh%kb(:, l, j)
1040 call spline_fit(ps_hgh%g%nrval, ps_hgh%g%rofi, hato, ps%kb(l, j), ps%projectors_sphere_threshold)
1041 end do
1042 end do
1043
1044 ! Now the part corresponding to the local pseudopotential
1045 ! where the asymptotic part is subtracted
1046 call spline_fit(ps_hgh%g%nrval, ps_hgh%g%rofi, ps_hgh%vlocal, ps%vl, ps%projectors_sphere_threshold)
1047
1048 ! Define the table for the pseudo-wavefunction components (using splines)
1049 ! with a correct normalization function
1050 do is = 1, ps%ispin
1051 dens = m_zero
1052 do l = 1, ps%conf%p
1053 hato(2:ps_hgh%g%nrval) = ps_hgh%rphi(2:ps_hgh%g%nrval, l)/ps_hgh%g%rofi(2:ps_hgh%g%nrval)
1054 hato(1) = hato(2)
1055
1056 do ip = 1, ps_hgh%g%nrval
1057 dens(ip) = dens(ip) + ps%conf%occ(l, is)*hato(ip)**2/(m_four*m_pi)
1058 end do
1059
1060 call spline_fit(ps_hgh%g%nrval, ps_hgh%g%rofi, hato, ps%ur(l, is), ps%projectors_sphere_threshold)
1061 end do
1062 call spline_fit(ps_hgh%g%nrval, ps_hgh%g%rofi, dens, ps%density(is), ps%projectors_sphere_threshold)
1063 end do
1064
1065 safe_deallocate_a(hato)
1066 safe_deallocate_a(dens)
1067
1068 pop_sub(hgh_load.get_splines)
1069 end subroutine get_splines
1070 end subroutine hgh_load
1071
1072
1073 ! ---------------------------------------------------------
1074 subroutine ps_grid_load(ps, ps_grid)
1075 type(ps_t), intent(inout) :: ps
1076 type(ps_in_grid_t), intent(in) :: ps_grid
1077
1078 push_sub(ps_grid_load)
1079
1080 ! Fixes some components of ps, read in ps_grid
1081 ps%z_val = ps_grid%zval
1082
1083 ps%nlcc = ps_grid%core_corrections
1084
1085 ps%h(0:ps%lmax, 1, 1) = ps_grid%dkbcos(1:ps%lmax+1)
1086
1087 ! Increasing radius a little, just in case.
1088 ! I have hard-coded a larger increase of the cutoff for the filtering.
1089 ps%rc_max = maxval(ps_grid%kb_radius(1:ps%lmax+1)) * 1.5_real64
1090
1091 ! now we fit the splines
1092 call get_splines(ps_grid%g)
1093
1094 ! Passes from Rydbergs to Hartrees.
1095 ps%h(0:ps%lmax,:,:) = ps%h(0:ps%lmax,:,:) / m_two
1097 pop_sub(ps_grid_load)
1098
1099 contains
1100
1101 subroutine get_splines(g)
1102 type(logrid_t), intent(in) :: g
1103
1104 real(real64), allocatable :: hato(:), dens(:)
1105 integer :: is, l, ir, nrc, ip
1106
1107 push_sub(ps_grid_load.get_splines)
1108
1109 safe_allocate(hato(1:g%nrval))
1110 safe_allocate(dens(1:g%nrval))
1111
1112 ! the wavefunctions
1113 do is = 1, ps%ispin
1114
1115 dens = m_zero
1116
1117 do l = 1, ps_grid%no_l_channels
1118 hato(2:) = ps_grid%rphi(2:, l, 1+is)/g%rofi(2:)
1119 hato(1) = first_point_extrapolate(g%rofi, hato)
1121 if(ps%conf%occ(l, is) > m_epsilon) then
1122 do ip = 1, g%nrval
1123 dens(ip) = dens(ip) + ps%conf%occ(l, is)*hato(ip)**2/(m_four*m_pi)
1124 end do
1125 end if
1126
1127 call spline_fit(g%nrval, g%rofi, hato, ps%ur(l, is), ps%projectors_sphere_threshold)
1128
1129 end do
1130 call spline_fit(g%nrval, g%rofi, dens, ps%density(is), ps%projectors_sphere_threshold)
1131 end do
1132
1133
1134 ! the Kleinman-Bylander projectors
1135 do l = 1, ps%lmax+1
1136 nrc = logrid_index(g, ps_grid%kb_radius(l)) + 1
1137 hato(1:nrc) = ps_grid%KB(1:nrc, l)
1138 hato(nrc+1:g%nrval) = m_zero
1139
1140 call spline_fit(g%nrval, g%rofi, hato, ps%kb(l-1, 1), ps%projectors_sphere_threshold)
1141 end do
1142
1143 ! Now the part corresponding to the local pseudopotential
1144 ! where the asymptotic part is subtracted
1145 hato(:) = ps_grid%vlocal(:)/m_two
1146 call spline_fit(g%nrval, g%rofi, hato, ps%vl, ps%projectors_sphere_threshold)
1147
1148 if (ps_grid%core_corrections) then
1149 ! find cutoff radius
1150 hato(2:) = ps_grid%chcore(2:)/(m_four*m_pi*g%r2ofi(2:))
1151
1152 do ir = g%nrval-1, 2, -1
1153 if (hato(ir) > eps) then
1154 nrc = ir + 1
1155 exit
1156 end if
1157 end do
1158
1159 hato(nrc:g%nrval) = m_zero
1160 hato(1) = first_point_extrapolate(g%rofi, hato)
1161
1162 call spline_fit(g%nrval, g%rofi, hato, ps%core, ps%projectors_sphere_threshold)
1163 end if
1164
1165 safe_deallocate_a(hato)
1166 safe_deallocate_a(dens)
1168 pop_sub(ps_grid_load.get_splines)
1169 end subroutine get_splines
1170 end subroutine ps_grid_load
1171
1172 ! ---------------------------------------------------------
1176 subroutine ps_xml_load(ps, ps_xml, namespace)
1177 type(ps_t), intent(inout) :: ps
1178 type(ps_xml_t), intent(in) :: ps_xml
1179 type(namespace_t), intent(in) :: namespace
1180
1181 integer :: ll, ip, is, ic, jc, ir, nrc, ii
1182 real(real64) :: rr, kbcos, kbnorm, dnrm, avgv, volume_element
1183 real(real64), allocatable :: vlocal(:), kbprojector(:), wavefunction(:), nlcc_density(:), dens(:,:)
1184 integer, allocatable :: cmap(:, :)
1185 real(real64), allocatable :: matrix(:, :), eigenvalues(:)
1186 logical :: density_is_known
1187 integer, allocatable :: order(:)
1188 real(real64), allocatable :: occ_tmp(:,:)
1189 real(real64), parameter :: tol_diagonal=1.0e-10_real64 ! tolerance for taking d_ij as a diagonal matrix
1190
1191 push_sub(ps_xml_load)
1192
1193 ! Vanderbilt-Kleinman-Bylander, i.e. more than one projector per l. Hamann uses 2 for ONCVPSP.
1194 if (pseudo_has_total_angular_momentum(ps_xml%pseudo)) then
1195 if (ps%file_format == pseudo_format_psp8) then
1196 ps%relativistic_treatment = proj_j_average
1197 message(1) = "SOC from PSP8 is not currently supported."
1198 message(2) = "Only scalar relativistic effects will be considered."
1199 call messages_warning(2, namespace=namespace)
1200 else
1201 ps%relativistic_treatment = proj_j_dependent
1202 end if
1203 end if
1204
1205 ps%nlcc = ps_xml%nlcc
1206
1207 ps%z_val = ps_xml%valence_charge
1208
1209 ps%has_density = ps_xml%has_density
1210
1211 ! We start with the local potential
1212
1213 ! Outside of the grid, we use the analytical long-range tail
1214 ! TODO: check usefulness: ps%g%nrval is initialized as ps_xml%grid_size
1215 ! This makes no sense, ps_xml%grid can only go up to ps_xml%grid_size
1216 ! We should just fit ps_xml%potential directly
1217 ! Also, if ip > ps_xml%grid_size, most likely ps_xml%grid is not defined
1218 safe_allocate(vlocal(1:ps%g%nrval))
1219 do ip = 1, ps_xml%grid_size
1220 rr = ps_xml%grid(ip)
1221 if (ip <= ps_xml%grid_size) then
1222 vlocal(ip) = ps_xml%potential(ip, ps%llocal)
1223 else
1224 vlocal(ip) = -ps_xml%valence_charge/rr
1225 end if
1226 end do
1227
1228 ! We then fit this by a spline
1229 call spline_fit(ps%g%nrval, ps%g%rofi, vlocal, ps%vl, ps%projectors_sphere_threshold)
1230
1231 safe_deallocate_a(vlocal)
1232
1233 safe_allocate(kbprojector(1:ps%g%nrval))
1234 safe_allocate(wavefunction(1:ps%g%nrval))
1235
1236 kbprojector = m_zero
1237 wavefunction = m_zero
1238
1239 density_is_known = .false.
1240
1241 ! We then proceed with the nonlocal projectors and the orbitals
1242 if (ps_xml%kleinman_bylander) then
1243
1244 safe_allocate(cmap(0:ps_xml%lmax, 1:ps_xml%nchannels))
1245
1246 ! We order the different projectors and create some mappings
1247 ! the order of the channels is determined by spin orbit and the j value
1248 do ll = 0, ps_xml%lmax
1249 do ic = 1, ps_xml%nchannels
1250 cmap(ll, ic) = ic
1251
1252 if (ll == 0) cycle
1253 if (ll == ps_xml%llocal) cycle
1254 if (ps%relativistic_treatment /= proj_j_dependent) cycle
1255 if (ic > pseudo_nprojectors_per_l(ps_xml%pseudo, ll)) cycle ! This occurs for O and F for pseudodojo FR PBE for instance
1256
1257 ! This is Octopus convention:
1258 ! We treat l+1/2 first and l-1/2 after, so we order the projectors this way
1259
1260 assert(mod(ps_xml%nchannels, 2)==0)
1261 if (pseudo_projector_2j(ps_xml%pseudo, ll, ic) == 2*ll - 1) then
1262 cmap(ll, ic) = int((ic-1)/2)*2 + 2
1263 else
1264 assert(pseudo_projector_2j(ps_xml%pseudo, ll, ic) == 2*ll + 1)
1265 cmap(ll, ic) = int((ic-1)/2)*2 + 1
1266 end if
1267
1268 end do
1270 ! check that all numbers are present for each l
1271 assert(sum(cmap(ll, 1:ps_xml%nchannels)) == (ps_xml%nchannels + 1)*ps_xml%nchannels/2)
1272 end do
1273
1274 assert(all(cmap >= 0 .and. cmap <= ps_xml%nchannels))
1275
1276
1277 ! Weight of the projectors
1278 ! This is a diagonal matrix after the following treatment
1279 ps%h = m_zero
1280
1281 ! We now take the matrix dij and we diagonalize it
1282 ! The nonlocal projectors are changed accordingly, and then fitted by splines
1283 if (pseudo_nprojectors(ps_xml%pseudo) > 0) then
1284 safe_allocate(matrix(1:ps_xml%nchannels, 1:ps_xml%nchannels))
1285 safe_allocate(eigenvalues(1:ps_xml%nchannels))
1286
1287 do ll = 0, ps_xml%lmax
1288
1289 if (is_diagonal(ps_xml%nchannels, ps_xml%dij(ll, :, :), tol_diagonal) .or. &
1290 pseudo_has_total_angular_momentum(ps_xml%pseudo)) then
1291 matrix = m_zero
1292 do ic = 1, ps_xml%nchannels
1293 eigenvalues(ic) = ps_xml%dij(ll, ic, ic)
1294 matrix(ic, ic) = m_one
1295 end do
1296 else
1297 ! diagonalize the coefficient matrix
1298 matrix(1:ps_xml%nchannels, 1:ps_xml%nchannels) = ps_xml%dij(ll, 1:ps_xml%nchannels, 1:ps_xml%nchannels)
1299 call lalg_eigensolve(ps_xml%nchannels, matrix, eigenvalues)
1300 end if
1301
1302 do ic = 1, ps_xml%nchannels
1303 kbprojector = m_zero
1304 do jc = 1, ps_xml%nchannels
1305 call lalg_axpy(ps_xml%grid_size, matrix(jc, ic), ps_xml%projector(:, ll, jc), kbprojector)
1306 end do
1307
1308 call spline_fit(ps%g%nrval, ps%g%rofi, kbprojector, ps%kb(ll, cmap(ll, ic)), ps%projectors_sphere_threshold)
1309
1310 ps%h(ll, cmap(ll, ic), cmap(ll, ic)) = eigenvalues(ic)
1311
1312 end do
1313 end do
1314
1315 safe_deallocate_a(matrix)
1316 safe_deallocate_a(eigenvalues)
1317 end if
1318
1319 ps%conf%p = ps_xml%nwavefunctions
1320
1321 ! If we do not have a density but we have wavefunctions, we compute the
1322 ! pseudo-atomic density, from the pseudo-atomic wavefunctions
1323 ! In case we are doing a spin-polarized calculation, it is better to use the
1324 ! wavefunctions with spin-dependent occupations, instead the spin unresolved density
1325 ! given in the pseudopotential
1326 if ((.not. ps_has_density(ps) .or. ps%ispin == 2) .and. ps_xml%nwavefunctions > 0) then
1327 safe_allocate(dens(1:ps%g%nrval, 1:ps%ispin))
1328 dens = m_zero
1329 end if
1330
1331 do ii = 1, ps_xml%nwavefunctions
1332
1333 ps%conf%n(ii) = ps_xml%wf_n(ii)
1334 ps%conf%l(ii) = ps_xml%wf_l(ii)
1335
1336 if (ps%ispin == 2) then
1337 ps%conf%occ(ii, 1) = min(ps_xml%wf_occ(ii), m_two*ps_xml%wf_l(ii) + m_one)
1338 ps%conf%occ(ii, 2) = ps_xml%wf_occ(ii) - ps%conf%occ(ii, 1)
1339 else
1340 ps%conf%occ(ii, 1) = ps_xml%wf_occ(ii)
1341 end if
1342
1343 ps%conf%j(ii) = m_zero
1344 if (pseudo_has_total_angular_momentum(ps_xml%pseudo)) then
1345 ps%conf%j(ii) = m_half*pseudo_wavefunction_2j(ps_xml%pseudo, ii)
1346 end if
1347
1348 do ip = 1, ps%g%nrval
1349 if (ip <= ps_xml%grid_size) then
1350 wavefunction(ip) = ps_xml%wavefunction(ip, ii)
1351 else
1352 wavefunction(ip) = m_zero
1353 end if
1354 end do
1355
1356 do is = 1, ps%ispin
1357 call spline_fit(ps%g%nrval, ps%g%rofi, wavefunction, ps%ur(ii, is), ps%projectors_sphere_threshold)
1358 end do
1359
1360 if (.not. ps_has_density(ps) .or. ps%ispin == 2) then
1361 do is = 1, ps%ispin
1362 do ip = 1, ps_xml%grid_size
1363 dens(ip, is) = dens(ip, is) + ps%conf%occ(ii, is)*wavefunction(ip)**2/(m_four*m_pi)
1364 end do
1365 end do
1366 end if
1367
1368 end do
1369
1370 !If we assigned all the valence electrons, we can compute the (spin-resolved) atomic density
1371 if ((.not. ps_has_density(ps) .or. ps%ispin == 2) .and. ps_xml%nwavefunctions > 0) then
1372 do is = 1, ps%ispin
1373 call spline_fit(ps%g%nrval, ps%g%rofi, dens(:,is), ps%density(is), ps%projectors_sphere_threshold)
1374 end do
1375 safe_deallocate_a(dens)
1376 density_is_known = .true.
1377 ps%has_density = .true.
1378 end if
1379
1380 safe_deallocate_a(cmap)
1381
1382 else !Not ps_xml%kleinman_bylander
1383
1384 !Get the occupations from the valence charge of the atom
1385 ps%conf%occ = m_zero
1386 ps%conf%symbol = ps%label(1:2)
1387 call ps_guess_atomic_occupations(namespace, ps%z, ps%z_val, ps%ispin, ps%conf)
1388
1389 ! In order to work in the following, we need to sort the occupations by angular momentum
1390 safe_allocate(order(1:ps_xml%lmax+1))
1391 safe_allocate(occ_tmp(1:ps_xml%lmax+1, 1:2))
1392 occ_tmp(:,:) = ps%conf%occ(1:ps_xml%lmax+1,1:2)
1393 call sort(ps%conf%l(1:ps_xml%lmax+1), order)
1394 do ll = 0, ps_xml%lmax
1395 ps%conf%occ(ll+1, 1:2) = occ_tmp(order(ll+1), 1:2)
1396 end do
1397 safe_deallocate_a(order)
1398 safe_deallocate_a(occ_tmp)
1399
1400 !If we assigned all the valence electrons, we can compute the (spin-resolved) atomic density
1401 if (abs(sum(ps%conf%occ) - ps%z_val ) < m_epsilon) then
1402 safe_allocate(dens(1:ps%g%nrval, 1:ps%ispin))
1403 dens = m_zero
1404 end if
1405
1406 do ll = 0, ps_xml%lmax
1407 ! we need to build the KB projectors
1408 ! the procedure was copied from ps_in_grid.F90 (r12967)
1409 dnrm = m_zero
1410 avgv = m_zero
1411 do ip = 1, ps_xml%grid_size
1412 rr = ps_xml%grid(ip)
1413 volume_element = rr**2*ps_xml%weights(ip)
1414 kbprojector(ip) = (ps_xml%potential(ip, ll) - ps_xml%potential(ip, ps%llocal))*ps_xml%wavefunction(ip, ll)
1415 dnrm = dnrm + kbprojector(ip)**2*volume_element
1416 avgv = avgv + kbprojector(ip)*ps_xml%wavefunction(ip, ll)*volume_element
1417 end do
1418
1419 kbcos = dnrm/(safe_tol(avgv,1.0e-20_real64))
1420 kbnorm = m_one/(safe_tol(sqrt(dnrm),1.0e-20_real64))
1421
1422 if (ll /= ps%llocal) then
1423 ps%h(ll, 1, 1) = kbcos
1424 kbprojector = kbprojector*kbnorm
1425 else
1426 ps%h(ll, 1, 1) = m_zero
1427 end if
1428
1429 call spline_fit(ps%g%nrval, ps%g%rofi, kbprojector, ps%kb(ll, 1), ps%projectors_sphere_threshold)
1430
1431 ! wavefunctions, for the moment we pad them with zero
1432 call lalg_copy(ps_xml%grid_size, ps_xml%wavefunction(:, ll), wavefunction)
1433 wavefunction(ps_xml%grid_size+1:ps%g%nrval) = m_zero
1434
1435 do is = 1, ps%ispin
1436 call spline_fit(ps%g%nrval, ps%g%rofi, wavefunction, ps%ur(ll + 1, is), ps%projectors_sphere_threshold)
1437 end do
1438
1439 !If we assigned all the valence electrons, we can compute the (spin-resolved) atomic density
1440 if (abs(sum(ps%conf%occ) - ps%z_val) < m_epsilon) then
1441 do is = 1, ps%ispin
1442 do ip = 1, ps_xml%grid_size
1443 dens(ip, is) = dens(ip, is) + ps%conf%occ(ll+1, is)*wavefunction(ip)**2/(m_four*m_pi)
1444 end do
1445 end do
1446 end if
1447 end do
1448
1449 !If we assigned all the valence electrons, we can compute the (spin-resolved) atomic density
1450 if (abs(sum(ps%conf%occ) - ps%z_val) < m_epsilon) then
1451 do is = 1, ps%ispin
1452 call spline_fit(ps%g%nrval, ps%g%rofi, dens(:,is), ps%density(is), ps%projectors_sphere_threshold)
1453 end do
1454 safe_deallocate_a(dens)
1455 density_is_known = .true.
1456 ps%has_density = .true.
1457 end if
1458 end if ! ps_xml%kleinman_bylander
1459
1460
1461 if (ps_has_density(ps) .and. .not. density_is_known) then
1462
1463 safe_allocate(dens(1:ps%g%nrval, 1))
1464
1465 dens(1:ps_xml%grid_size, 1) = ps_xml%density(1:ps_xml%grid_size)/ps%ispin
1466 dens(ps_xml%grid_size + 1:ps%g%nrval, 1) = m_zero
1467
1468 do is = 1, ps%ispin
1469 call spline_fit(ps%g%nrval, ps%g%rofi, dens(:, 1), ps%density(is), ps%projectors_sphere_threshold)
1470 end do
1471
1472 safe_deallocate_a(dens)
1473 end if
1474
1475 ! Non-linear core-corrections
1476 ! We truncate the NLCC when below eps
1477 if (ps_xml%nlcc) then
1478
1479 safe_allocate(nlcc_density(1:ps%g%nrval))
1480
1481 call lalg_copy(ps_xml%grid_size, ps_xml%nlcc_density, nlcc_density)
1482
1483 ! find cutoff radius
1484 do ir = ps_xml%grid_size - 1, 1, -1
1485 if (nlcc_density(ir) > eps) then
1486 nrc = ir + 1
1487 exit
1488 end if
1489 end do
1490
1491 nlcc_density(nrc:ps%g%nrval) = m_zero
1492
1493 call spline_fit(ps%g%nrval, ps%g%rofi, nlcc_density, ps%core, ps%projectors_sphere_threshold)
1494
1495 safe_deallocate_a(nlcc_density)
1496 end if
1497
1498 call ps_getradius(ps)
1499
1500 !To be consistent with the other pseudopotentials, we are increasing the radius here
1501 ps%rc_max = ps%rc_max*1.05_real64
1502
1503 safe_deallocate_a(kbprojector)
1504 safe_deallocate_a(wavefunction)
1505
1506 pop_sub(ps_xml_load)
1507 end subroutine ps_xml_load
1508
1509
1510 ! ---------------------------------------------------------
1512 pure integer function ps_niwfs(ps)
1513 type(ps_t), intent(in) :: ps
1514
1515 integer :: i, l
1516
1517 ps_niwfs = 0
1518 do i = 1, ps%conf%p
1519 l = ps%conf%l(i)
1520 ps_niwfs = ps_niwfs + (2*l+1)
1521 end do
1522
1523 end function ps_niwfs
1524
1525 ! ---------------------------------------------------------
1527 pure integer function ps_bound_niwfs(ps)
1528 type(ps_t), intent(in) :: ps
1529
1530 integer :: i, l
1531
1532 ps_bound_niwfs = 0
1533 do i = 1, ps%conf%p
1534 l = ps%conf%l(i)
1535 if (any(.not. ps%bound(i,:))) cycle
1536 ps_bound_niwfs = ps_bound_niwfs + (2*l+1)
1537 end do
1538
1539 end function ps_bound_niwfs
1540
1541 !---------------------------------------
1542
1543 pure logical function ps_has_density(ps) result(has_density)
1544 type(ps_t), intent(in) :: ps
1545
1546 has_density = ps%has_density
1547
1548 end function ps_has_density
1549
1550 !---------------------------------------
1551
1552 pure logical function ps_has_nlcc(ps) result(has_nlcc)
1553 type(ps_t), intent(in) :: ps
1554
1555 has_nlcc = ps%nlcc
1556
1557 end function ps_has_nlcc
1558
1559 !---------------------------------------
1560 real(real64) function ps_density_volume(ps, namespace) result(volume)
1561 type(ps_t), intent(in) :: ps
1562 type(namespace_t), intent(in) :: namespace
1563
1564 integer :: ip, ispin
1565 real(real64) :: rr
1566 real(real64), allocatable ::vol(:)
1567 type(spline_t) :: volspl
1568
1569 push_sub(ps_density_volume)
1570
1571 if (.not. ps_has_density(ps)) then
1572 message(1) = "The pseudopotential does not contain an atomic density"
1573 call messages_fatal(1, namespace=namespace)
1574 end if
1575
1576 safe_allocate(vol(1:ps%g%nrval))
1577
1578 do ip = 1, ps%g%nrval
1579 rr = ps%g%rofi(ip)
1580 vol(ip) = m_zero
1581 do ispin = 1, ps%ispin
1582 vol(ip) = vol(ip) + spline_eval(ps%density(ispin), rr)*m_four*m_pi*rr**5
1583 end do
1584 end do
1585
1586 call spline_init(volspl)
1587 call spline_fit(ps%g%nrval, ps%g%rofi, vol, volspl, ps%projectors_sphere_threshold)
1588 volume = spline_integral(volspl)
1589 call spline_end(volspl)
1590
1591 safe_deallocate_a(vol)
1592
1593 pop_sub(ps_density_volume)
1594 end function ps_density_volume
1595
1599 subroutine ps_guess_atomic_occupations(namespace, zz, valcharge, ispin, conf)
1600 type(namespace_t), intent(in) :: namespace
1601 real(real64), intent(in) :: zz
1602 real(real64), intent(in) :: valcharge
1603 integer, intent(in) :: ispin
1604 type(valconf_t), intent(inout) :: conf
1606 real(real64) :: val
1607
1609
1610 val = valcharge
1611 conf%p = 0
1612 conf%l = 0
1613 conf%j = m_zero
1614 conf%n = 0
1615
1616 write(message(1), '(a,a,a)') 'Debug: Guessing the atomic occupations for ', trim(conf%symbol), "."
1617 call messages_info(1, namespace=namespace, debug_only=.true.)
1618
1619 assert(valcharge <= zz)
1621 ! Here we populate the core states
1622 ! 1s state for all atoms after He
1623 if(int(zz) > 2 .and. val > zz - 2) then
1624 call fill_s_orbs(val, 2, 1)
1625 end if
1626 ! 2s state for all atoms after Be
1627 if(int(zz) > 4 .and. val > zz - 4) then
1628 call fill_s_orbs(val, 2, 2)
1629 end if
1630 ! 2p state for all atoms after Ne
1631 ! For pseudopotentials Al-Ar, we fill the 2s but not the 2p
1632 if(int(zz) > 18 .and. val > zz - 10) then
1633 call fill_p_orbs(val, 6, 2)
1634 end if
1635 ! 3s state for all atoms after Mg
1636 if(int(zz) > 12 .and. val > zz - 12) then
1637 call fill_s_orbs(val, 2, 3)
1638 end if
1639 ! 3p state for all atoms after Ar
1640 if(int(zz) > 18 .and. val > zz - 18) then
1641 call fill_p_orbs(val, 6, 3)
1642 end if
1643 ! 3d states for all atoms after Ni
1644 if(int(zz) > 28 .and. val > zz - 28) then
1645 call fill_d_orbs(val, 10, 3)
1646 end if
1647 ! 4s states for all atoms after Zn
1648 if(int(zz) > 30 .and. val > zz - 30) then
1649 call fill_s_orbs(val, 2, 4)
1650 end if
1651 ! 4p states for all atoms after Kr
1652 if(int(zz) > 36 .and. val > zz - 36) then
1653 call fill_p_orbs(val, 6, 4)
1654 end if
1655 ! 4d states for all atoms after Pd
1656 if(int(zz) > 46 .and. val > zz - 46) then
1657 call fill_d_orbs(val, 10, 4)
1658 end if
1659 ! 5s states for all atoms after Cd
1660 ! For Z=71 ot Z=80, the 4f is filled before the 5s/5p
1661 if((int(zz) > 48 .and. val > zz - 48) .or. &
1662 (int(zz) > 70 .and. int(zz) <= 81 .and. val > zz - 62)) then
1663 call fill_s_orbs(val, 2, 5)
1664 end if
1665 ! 5p states for all atoms after Xe
1666 ! For Z=71 ot Z=80, the 4f is filled before the 5s/5p
1667 if((int(zz) > 54 .and. val > zz - 54) .or. &
1668 (int(zz) > 70 .and. int(zz) <= 81 .and. val > zz - 68) ) then
1669 call fill_p_orbs(val, 6, 5)
1670 end if
1671 ! 4f states for all atoms after Yb
1672 ! Only after Z=80 for pseudopotentials
1673 if(int(zz) > 80 .and. val > zz - 68) then
1674 call fill_f_orbs(val, 14, 4)
1675 end if
1676
1677
1678 ! We leave here the valence states
1679 select case (int(zz))
1680 case (1)
1681 call fill_s_orbs(val, 1, 1) ! H 1s^1
1682 case (2)
1683 call fill_s_orbs(val, 2, 1) ! He 1s^2
1684 case (3)
1685 call fill_s_orbs(val, 1, 2) ! Li 1s^2 2s^1
1686 case (4)
1687 call fill_s_orbs(val, 2, 2) ! Be 1s^2 ; 2s^2
1688 case (5)
1689 call fill_p_orbs(val, 1, 2) ! B 1s^2 ; 2s^2 2p^1
1690 case (6)
1691 call fill_p_orbs(val, 2, 2) ! C 1s^2 ; 2s^2 2p^2
1692 case (7)
1693 call fill_p_orbs(val, 3, 2) ! N 1s^2 ; 2s^2 2
1694 case (8)
1695 call fill_p_orbs(val, 4, 2) ! O 1s^2 ; 2s^2 2p^4
1696 case (9)
1697 call fill_p_orbs(val, 5, 2) ! F 1s^2 ; 2s^2 2p^5
1698 case (10)
1699 call fill_p_orbs(val, 6, 2) ! Ne 1s^2 ; 2s^2 2p^6
1700 case (11)
1701 if(val > 6) call fill_p_orbs(val, 6, 2)
1702 call fill_s_orbs(val, 1, 3) ! Na 1s^2 2s^2 2p^6 ; 3s^1
1703 case (12)
1704 if(val > 6) call fill_p_orbs(val, 6, 2)
1705 call fill_s_orbs(val, 2, 3) ! Mg 1s^2 2s^2 2p^6 ; 3s^2
1706 case (13)
1707 if(val > 6) call fill_p_orbs(val, 6, 2)
1708 call fill_p_orbs(val, 1, 3) ! Al 1s^2 2s^2 2p^6 ; 3s^2 3p^1
1709 case (14)
1710 if(val > 6) call fill_p_orbs(val, 6, 2)
1711 call fill_p_orbs(val, 2, 3) ! Si 1s^2 2s^2 2p^6 ; 3s^2 3p^2
1712 case (15)
1713 if(val > 6) call fill_p_orbs(val, 6, 2)
1714 call fill_p_orbs(val, 3, 3) ! P 1s^2 2s^2 2p^6 ; 3s^2 3p^3
1715 case (16)
1716 if(val > 6) call fill_p_orbs(val, 6, 2)
1717 call fill_p_orbs(val, 4, 3) ! S 1s^2 2s^2 2p^6 ; 3s^2 3p^4
1718 case (17)
1719 if(val > 6) call fill_p_orbs(val, 6, 2)
1720 call fill_p_orbs(val, 5, 3) ! Cl 1s^2 2s^2 2p^6 ; 3s^2 3p^5
1721 case (18)
1722 if(val > 6) call fill_p_orbs(val, 6, 2)
1723 call fill_p_orbs(val, 6, 3) ! Ar 1s^2 2s^2 2p^6 ; 3s^2 3p^6
1724 case (19)
1725 call fill_s_orbs(val, 1, 4) ! K 1s^2 2s^2 2p^6 3s^2 3p^6 ; 4s^1
1726 case (20)
1727 call fill_s_orbs(val, 2, 4) ! Ca 1s^2 2s^2 2p^6 3s^2 3p^6 ; 4s^2
1728 case (21)
1729 if (val > 1) call fill_s_orbs(val, 2, 4)
1730 call fill_d_orbs(val, 1, 3) ! Sc 1s^2 2s^2 2p^6 3s^2 3p^6 ; 4s^2 3d^1
1731 case (22)
1732 if (val > 2) call fill_s_orbs(val, 2, 4)
1733 call fill_d_orbs(val, 2, 3) ! Ti 1s^2 2s^2 2p^6 3s^2 3p^6 ; 4s^2 3d^2
1734 case (23)
1735 if (val > 3) call fill_s_orbs(val, 2, 4)
1736 call fill_d_orbs(val, 3, 3) ! V 1s^2 2s^2 2p^6 3s^2 3p^6 ; 4s^2 3d^3
1737 case (24)
1738 if (val > 4) call fill_s_orbs(val, 2, 4)
1739 call fill_d_orbs(val, 4, 3) ! Cr 1s^2 2s^2 2p^6 3s^2 3p^6 ; 4s^2 3d^4
1740 case (25)
1741 if (val > 5) call fill_s_orbs(val, 2, 4)
1742 call fill_d_orbs(val, 5, 3) ! Mn 1s^2 2s^2 2p^6 3s^2 3p^6 ; 4s^2 3d^5
1743 case (26)
1744 if (val > 6) call fill_s_orbs(val, 2, 4)
1745 call fill_d_orbs(val, 6, 3) ! Fe 1s^2 2s^2 2p^6 3s^2 3p^6 ; 4s^2 3d^6
1746 case (27)
1747 if (val > 7) call fill_s_orbs(val, 2, 4)
1748 call fill_d_orbs(val, 7, 3) ! Co 1s^2 2s^2 2p^6 3s^2 3p^6 ; 4s^2 3d^7
1749 case (28)
1750 if (val > 8 ) call fill_s_orbs(val, 2, 4)
1751 call fill_d_orbs(val, 8, 3) ! Ni 1s^2 2s^2 2p^6 3s^2 3p^6 ; 4s^2 3d^8
1752 case (29)
1753 call fill_s_orbs(val, 1, 4) ! Cu 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 ; 4s^1
1754 case (30)
1755 call fill_s_orbs(val, 2, 4) ! Zn 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 ; 4s^2
1756 case (31)
1757 call fill_p_orbs(val, 1, 4) ! Ga 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 ; 4s^2 4p^1
1758 case (32)
1759 call fill_p_orbs(val, 2, 4) ! Ge 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 ; 4s^2 4p^2
1760 case (33)
1761 call fill_p_orbs(val, 3, 4) ! As 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 ; 4s^2 4p^3
1762 case (34)
1763 call fill_p_orbs(val, 4, 4) ! Se 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 ; 4s^2 4p^4
1764 case (35)
1765 call fill_p_orbs(val, 5, 4) ! Br 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 ; 4s^2 4p^5
1766 case (36)
1767 call fill_p_orbs(val, 6, 4) ! Kr 1s^2 2s^2 2p^6 3s^2 3p^6 3d^10 ; 4s^2 4p^6
1768 case (37)
1769 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
1770 case (38)
1771 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
1772 case (39)
1773 if (val > 2) call fill_d_orbs(val, 1, 4)
1774 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
1775 case (40)
1776 if (val > 2) call fill_d_orbs(val, 2, 4)
1777 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
1778 case (41)
1779 if (val > 1) call fill_d_orbs(val, 4, 4)
1780 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
1781 case (42)
1782 if (val > 1) call fill_d_orbs(val, 5, 4)
1783 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
1784 case (43)
1785 if (val > 2) call fill_d_orbs(val, 5, 4)
1786 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
1787 case (44)
1788 if (val > 1) call fill_d_orbs(val, 7, 4)
1789 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
1790 case (45)
1791 if (val > 1) call fill_d_orbs(val, 8, 4)
1792 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
1793 case (46)
1794 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
1795 case (47)
1796 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
1797 case (48)
1798 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
1799 case (49)
1800 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
1801 case (50)
1802 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
1803 case (51)
1804 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
1805 case (52)
1806 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
1807 case (53)
1808 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
1809 case (54)
1810 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
1811 case (55)
1812 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
1813 case (56)
1814 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
1815 case (57)
1816 if (val > 2) call fill_d_orbs(val, 1, 5)
1817 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
1818 case (58)
1819 if (val > 3) call fill_f_orbs(val, 1, 4)
1820 if (val > 2) call fill_d_orbs(val, 1, 5)
1821 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
1822 case (59)
1823 if (val > 2) call fill_f_orbs(val, 3, 4)
1824 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
1825 case (60)
1826 if (val > 2) call fill_f_orbs(val, 4, 4)
1827 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
1828 case (61)
1829 if (val > 2) call fill_f_orbs(val, 5, 4)
1830 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
1831 case (62)
1832 if (val > 2) call fill_f_orbs(val, 6, 4)
1833 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
1834 case (63)
1835 if (val > 2) call fill_f_orbs(val, 7, 4)
1836 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
1837 case (64)
1838 if (val > 3) call fill_f_orbs(val, 7, 4)
1839 if (val > 2) call fill_d_orbs(val, 1, 5)
1840 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
1841 case (65)
1842 if (val > 2) call fill_f_orbs(val, 9, 4)
1843 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
1844 case (66)
1845 if (val > 2) call fill_f_orbs(val, 10, 4)
1846 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
1847 case (67)
1848 if (val > 2) call fill_f_orbs(val, 11, 4)
1849 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
1850 case (68)
1851 if (val > 2) call fill_f_orbs(val, 12, 4)
1852 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
1853 case (69)
1854 if (val > 2) call fill_f_orbs(val, 13, 4)
1855 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
1856 case (70)
1857 if (val > 2) call fill_f_orbs(val, 14, 4)
1858 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
1859 case (71)
1860 if (val > 3) call fill_f_orbs(val, 14, 4)
1861 if (val > 2) call fill_d_orbs(val, 1, 5)
1862 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
1863 case (72)
1864 if (val > 4) call fill_f_orbs(val, 14, 4)
1865 if (val > 2) call fill_d_orbs(val, 2, 5)
1866 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
1867 case (73)
1868 if (val > 5) call fill_f_orbs(val, 14, 4)
1869 if (val > 2) call fill_d_orbs(val, 3, 5)
1870 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
1871 case (74)
1872 if (val > 6) call fill_f_orbs(val, 14, 4)
1873 if (val > 2) call fill_d_orbs(val, 4, 5)
1874 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
1875 case (75)
1876 if (val > 7) call fill_f_orbs(val, 14, 4)
1877 if (val > 2) call fill_d_orbs(val, 5, 5)
1878 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
1879 case (76)
1880 if (val > 8) call fill_f_orbs(val, 14, 4)
1881 if (val > 2) call fill_d_orbs(val, 6, 5)
1882 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
1883 case (77)
1884 if (val > 9) call fill_f_orbs(val, 14, 4)
1885 if (val > 2) call fill_d_orbs(val, 7, 5)
1886 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
1887 case (78)
1888 if (val > 10) call fill_f_orbs(val, 14, 4)
1889 if (val > 1) call fill_d_orbs(val, 9, 5)
1890 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
1891 case (79)
1892 if (val > 21) call fill_f_orbs(val, 14, 4)
1893 if (val > 1) call fill_d_orbs(val, 10, 5)
1894 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
1895 case (80)
1896 if (val > 12) call fill_f_orbs(val, 14, 4)
1897 if (val > 2) call fill_d_orbs(val, 10, 5)
1898 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
1899 case (81)
1900 if (val > 10) call fill_d_orbs(val, 10, 5)
1901 if (val > 1) call fill_s_orbs(val, 2, 6)
1902 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
1903 case (82)
1904 if (val > 10) call fill_d_orbs(val, 10, 5)
1905 if (val > 2) call fill_s_orbs(val, 2, 6)
1906 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
1907 case (83)
1908 if (val > 10) call fill_d_orbs(val, 10, 5)
1909 if (val > 3) call fill_s_orbs(val, 2, 6)
1910 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
1911 case (84)
1912 if (val > 10) call fill_d_orbs(val, 10, 5)
1913 if (val > 4) call fill_s_orbs(val, 2, 6)
1914 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
1915 case (85)
1916 if (val > 10) call fill_d_orbs(val, 10, 5)
1917 if (val > 5) call fill_s_orbs(val, 2, 6)
1918 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
1919 case (86)
1920 if (val > 10) call fill_d_orbs(val, 10, 5)
1921 if (val > 6) call fill_s_orbs(val, 2, 6)
1922 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
1923 case (87)
1924 if (val > 1) call fill_p_orbs(val, 6, 6)
1925 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
1926 case (88)
1927 if (val > 2) call fill_p_orbs(val, 6, 6)
1928 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
1929 case (89)
1930 if (val > 3) call fill_p_orbs(val, 6, 6)
1931 if (val > 2) call fill_d_orbs(val, 1, 6)
1932 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
1933 case (90)
1934 if (val > 4) call fill_p_orbs(val, 6, 6)
1935 if (val > 2) call fill_d_orbs(val, 2, 6)
1936 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
1937 case (91)
1938 if (val > 5) call fill_p_orbs(val, 6, 6)
1939 if (val > 3) call fill_f_orbs(val, 2, 5)
1940 if (val > 2) call fill_d_orbs(val, 1, 6)
1941 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
1942 case (92)
1943 if (val > 6) call fill_p_orbs(val, 6, 6)
1944 if (val > 3) call fill_f_orbs(val, 3, 5)
1945 if (val > 2) call fill_d_orbs(val, 1, 6)
1946 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
1947 case (93)
1948 if (val > 7) call fill_p_orbs(val, 6, 6)
1949 if (val > 3) call fill_f_orbs(val, 4, 5)
1950 if (val > 2) call fill_d_orbs(val, 1, 6)
1951 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
1952 case (94)
1953 if (val > 8) call fill_p_orbs(val, 6, 6)
1954 if (val > 2) call fill_f_orbs(val, 6, 5)
1955 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
1956 case (95)
1957 if (val > 9) call fill_p_orbs(val, 6, 6)
1958 if (val > 2) call fill_f_orbs(val, 7, 5)
1959 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
1960 case (96)
1961 if (val > 10) call fill_p_orbs(val, 6, 6)
1962 if (val > 3) call fill_f_orbs(val, 7, 5)
1963 if (val > 2) call fill_d_orbs(val, 1, 6)
1964 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
1965 case (97)
1966 if (val > 11) call fill_p_orbs(val, 6, 6)
1967 if (val > 2) call fill_f_orbs(val, 9, 5)
1968 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
1969 case (98)
1970 if (val > 12) call fill_p_orbs(val, 6, 6)
1971 if (val > 2) call fill_f_orbs(val, 10, 5)
1972 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
1973 case (99)
1974 if (val > 13) call fill_p_orbs(val, 6, 6)
1975 if (val > 2) call fill_f_orbs(val, 11, 5)
1976 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
1977 case (100)
1978 if (val > 14) call fill_p_orbs(val, 6, 6)
1979 if (val > 2) call fill_f_orbs(val, 12, 5)
1980 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
1981 case (101)
1982 if (val > 15) call fill_p_orbs(val, 6, 6)
1983 if (val > 2) call fill_f_orbs(val, 13, 5)
1984 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
1985 case (102)
1986 if (val > 16) call fill_p_orbs(val, 6, 6)
1987 if (val > 2) call fill_f_orbs(val, 14, 5)
1988 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
1989 case (103)
1990 if (val > 3) call fill_f_orbs(val, 14, 5)
1991 if (val > 1) call fill_s_orbs(val, 2, 7)
1992 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
1993 end select
1994
1995 !If we attributed all the electrons, everything went fine
1996 if (val < m_epsilon) then
1997 !In case of spin-polarized calculations, we properly distribute the electrons
1998 if (ispin == 2) then
1999 call valconf_unpolarized_to_polarized(conf)
2000 end if
2001 else
2002 conf%occ = m_zero
2003 message(1) = "Error in attributing atomic occupations"
2004 call messages_warning(1, namespace=namespace)
2005 end if
2006
2008
2009 contains
2010 subroutine fill_s_orbs(val, max_occ, nn)
2011 real(real64), intent(inout) :: val
2012 integer, intent(in) :: max_occ, nn
2013
2014 conf%p = conf%p + 1
2015 conf%occ(conf%p,1) = min(val, real(max_occ, real64) )
2016 val = val - conf%occ(conf%p,1)
2017 conf%l(conf%p) = 0
2018 conf%n(conf%p) = nn
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:
constant times a vector plus a vector
Definition: lalg_basic.F90:171
Copies a vector x, to a vector y.
Definition: lalg_basic.F90:186
This is the common interface to a sorting routine. It performs the shell algorithm,...
Definition: sort.F90:149
Some operations may be done for one spline-function, or for an array of them.
Definition: splines.F90:179
The integral may be done with or without integration limits, but we want the interface to be common.
Definition: splines.F90:173
double exp(double __x) __attribute__((__nothrow__
double sqrt(double __x) __attribute__((__nothrow__
subroutine, public valconf_unpolarized_to_polarized(conf)
Definition: atomic.F90:236
subroutine, public valconf_copy(cout, cin)
Definition: atomic.F90:166
real(real64), parameter, public m_two
Definition: global.F90:190
real(real64), parameter, public m_zero
Definition: global.F90:188
real(real64), parameter, public m_four
Definition: global.F90:192
real(real64), parameter, public m_pi
some mathematical constants
Definition: global.F90:186
real(real64), parameter, public m_min_exp_arg
Definition: global.F90:207
real(real64), parameter, public m_epsilon
Definition: global.F90:204
real(real64), parameter, public m_half
Definition: global.F90:194
real(real64), parameter, public m_one
Definition: global.F90:189
real(real64), parameter, public m_three
Definition: global.F90:191
Definition: io.F90:114
subroutine, public io_close(iunit, grp)
Definition: io.F90:418
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
Definition: io.F90:352
subroutine, public logrid_copy(grid_in, grid_out)
Definition: logrid.F90:228
integer function, public logrid_index(grid, rofi)
Definition: logrid.F90:256
subroutine, public logrid_end(grid)
Definition: logrid.F90:213
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:115
logical pure function, public is_diagonal(dim, matrix, tol)
Returns true is the matrix of size dim x dim is diagonal, given a relative tolerance.
Definition: math.F90:1475
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:537
subroutine, public messages_obsolete_variable(namespace, name, rep)
Definition: messages.F90:1045
subroutine, public messages_new_line()
Definition: messages.F90:1134
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:160
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:414
subroutine, public messages_input_error(namespace, var, details, row, column)
Definition: messages.F90:713
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:616
subroutine, public profiling_out(label)
Increment out counter and sum up difference between entry and exit time.
Definition: profiling.F90:623
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
Definition: profiling.F90:552
subroutine, public ps_cpi_end(ps_cpi)
Definition: ps_cpi.F90:183
subroutine, public ps_cpi_init(ps_cpi, filename, namespace)
Definition: ps_cpi.F90:148
subroutine, public ps_cpi_process(ps_cpi, lloc, namespace)
Definition: ps_cpi.F90:221
subroutine, public ps_fhi_end(ps_fhi)
Definition: ps_fhi.F90:183
subroutine, public ps_fhi_init(ps_fhi, filename, namespace)
Definition: ps_fhi.F90:149
subroutine, public ps_fhi_process(ps_fhi, lmax, lloc, namespace)
Definition: ps_fhi.F90:202
subroutine, public hgh_end(psp)
Definition: ps_hgh.F90:228
subroutine, public hgh_get_eigen(psp, eigen)
Definition: ps_hgh.F90:289
subroutine, public hgh_process(psp, namespace)
Definition: ps_hgh.F90:247
subroutine, public hgh_init(psp, filename, namespace)
Definition: ps_hgh.F90:183
real(real64) function, public first_point_extrapolate(x, y, high_order)
Definition: ps_in_grid.F90:437
Definition: ps.F90:114
pure logical function, public ps_has_density(ps)
Definition: ps.F90:1637
integer, parameter, public ps_filter_ts
Definition: ps.F90:162
subroutine, public ps_end(ps)
Definition: ps.F90:1057
subroutine hgh_load(ps, ps_hgh)
Definition: ps.F90:1097
subroutine, public ps_init(ps, namespace, label, z, user_lmax, user_llocal, ispin, filename)
Definition: ps.F90:260
subroutine ps_info(ps, filename, namespace)
Definition: ps.F90:579
subroutine, public ps_separate(ps)
Separate the local potential into (soft) long-ranged and (hard) short-ranged parts.
Definition: ps.F90:690
subroutine ps_check_bound(ps, eigen)
Definition: ps.F90:881
subroutine ps_grid_load(ps, ps_grid)
Definition: ps.F90:1168
subroutine, public ps_debug(ps, dir, namespace, gmax)
Definition: ps.F90:926
integer, parameter, public proj_hgh
Definition: ps.F90:167
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:1621
subroutine, public ps_getradius(ps)
Definition: ps.F90:758
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:1693
pure logical function, public ps_has_nlcc(ps)
Definition: ps.F90:1646
pure integer function, public ps_niwfs(ps)
Returns the number of atomic orbitals taking into account then m quantum number multiplicity.
Definition: ps.F90:1606
subroutine, public ps_pspio_init(ps, namespace, label, z, lmax, ispin, filename)
Definition: ps.F90:2170
integer, parameter, public proj_j_average
Fully-relativistic pseudopotentials with separate j-average and SOC terms.
Definition: ps.F90:173
integer, parameter, public proj_rkb
Definition: ps.F90:167
integer, parameter, public proj_j_dependent
Fully-relativistic j-dependent pseudopotentials.
Definition: ps.F90:173
subroutine, public ps_derivatives(ps)
Definition: ps.F90:778
subroutine, public ps_filter(ps, filter, gmax)
Definition: ps.F90:795
real(real64) function, public ps_density_volume(ps, namespace)
Definition: ps.F90:1654
integer, parameter, public ps_filter_bsb
Definition: ps.F90:162
integer, parameter, public proj_kb
Definition: ps.F90:167
subroutine ps_xml_load(ps, ps_xml, namespace)
Loads XML files for QSO, UPF1, UPF2, PSML, and PSP8 formats.
Definition: ps.F90:1270
subroutine, public ps_psf_init(pstm, ispin, filename, namespace)
Definition: ps_psf.F90:151
subroutine, public ps_psf_process(ps_psf, namespace, lmax, lloc)
Definition: ps_psf.F90:301
subroutine, public ps_psf_get_eigen(ps_psf, eigen)
Definition: ps_psf.F90:335
subroutine, public ps_psf_end(ps_psf)
Definition: ps_psf.F90:202
subroutine, public ps_xml_end(this)
Definition: ps_xml.F90:337
subroutine, public ps_xml_init(this, namespace, filename, fmt, ierr)
Definition: ps_xml.F90:163
integer, parameter, public pseudo_format_file_not_found
Definition: pseudo.F90:173
integer, parameter, public pseudo_type_kleinman_bylander
Definition: pseudo.F90:158
integer, parameter, public pseudo_format_qso
Definition: pseudo.F90:173
integer, parameter, public pseudo_format_hgh
Definition: pseudo.F90:173
integer, parameter, public pseudo_format_upf2
Definition: pseudo.F90:173
integer, parameter, public pseudo_type_paw
Definition: pseudo.F90:158
integer, parameter, public pseudo_format_cpi
Definition: pseudo.F90:173
integer, parameter, public pseudo_exchange_unknown
Definition: pseudo.F90:188
integer, parameter, public pseudo_type_semilocal
Definition: pseudo.F90:158
integer, parameter, public pseudo_type_ultrasoft
Definition: pseudo.F90:158
logical pure function, public pseudo_has_total_angular_momentum(pseudo)
Definition: pseudo.F90:548
integer, parameter, public pseudo_correlation_unknown
Definition: pseudo.F90:192
integer, parameter, public pseudo_format_psf
Definition: pseudo.F90:173
integer, parameter, public pseudo_format_fhi
Definition: pseudo.F90:173
integer, parameter, public pseudo_format_unknown
Definition: pseudo.F90:173
integer, parameter, public pseudo_format_psml
Definition: pseudo.F90:173
This module is intended to contain "only mathematical" functions and procedures.
Definition: sort.F90:117
subroutine, public spline_filter_mask(spl, l, rmax, qmax, alpha, gamma, threshold)
subroutine, public spline_filter_bessel(spl, l, qmax, alpha, beta_fs, rcut, beta_rs, threshold)
subroutine, public spline_der(spl, dspl, threshold)
Definition: splines.F90:921
subroutine, public spline_force_pos(spl, threshold)
Definition: splines.F90:841
subroutine, public spline_fit(nrc, rofi, ffit, spl, threshold)
Definition: splines.F90:413
subroutine, public spline_3dft(spl, splw, threshold, gmax)
Definition: splines.F90:599
subroutine, public spline_besselft(spl, splw, l, threshold, gmax)
Definition: splines.F90:677
real(real64) function, public spline_eval(spl, x)
Definition: splines.F90:441
subroutine fill_s_orbs(val, max_occ, nn)
Definition: ps.F90:2104
subroutine fill_p_orbs(val, max_occ, nn)
Definition: ps.F90:2115
subroutine fill_d_orbs(val, max_occ, nn)
Definition: ps.F90:2126
subroutine get_splines()
Definition: ps.F90:1121
subroutine fill_f_orbs(val, max_occ, nn)
Definition: ps.F90:2137
remember that the FHI format is basically the CPI format with a header
Definition: ps_fhi.F90:137
The following data type contains: (a) the pseudopotential parameters, as read from a *....
Definition: ps_hgh.F90:144
A type storing the information and data about a pseudopotential.
Definition: ps.F90:184
the basic spline datatype
Definition: splines.F90:156
int true(void)