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