51 integer,
public,
parameter :: &
52 SPECIES_PSEUDO = 7, & !< pseudopotential
61 type(ps_t),
allocatable,
public :: ps
64 integer :: user_llocal
66 integer,
public :: pseudopotential_set_id
67 logical,
public :: pseudopotential_set_initialized
68 type(pseudo_set_t),
public :: pseudopotential_set
92 procedure pseudopotential_constructor
104 class(pseudopotential_t),
pointer :: spec
105 character(len=*),
intent(in) :: label
106 integer,
intent(in) :: index
121 spec%pseudopotential_set_id = option__pseudopotentialset__none
122 spec%pseudopotential_set_initialized = .false.
131 type(pseudopotential_t),
intent(inout) :: spec
135 if (spec%pseudopotential_set_initialized)
then
137 spec%pseudopotential_set_initialized = .false.
140 if (
allocated(spec%ps))
call ps_end(spec%ps)
141 safe_deallocate_a(spec%ps)
149 logical pure function pseudopotential_has_nlcc(spec)
150 class(pseudopotential_t),
intent(in) :: spec
151 pseudopotential_has_nlcc = spec%nlcc
163 select case (spec%pseudopotential_set_id)
165 option__pseudopotentialset__standard, &
166 option__pseudopotentialset__hgh_lda, &
167 option__pseudopotentialset__hgh_lda_sc, &
168 option__pseudopotentialset__hscv_lda)
172 case (option__pseudopotentialset__hscv_pbe)
188 select case (spec%pseudopotential_set_id)
190 option__pseudopotentialset__standard, &
191 option__pseudopotentialset__hgh_lda, &
192 option__pseudopotentialset__hgh_lda_sc, &
193 option__pseudopotentialset__hscv_lda)
197 case (option__pseudopotentialset__hscv_pbe)
209 radius = spec%ps%rc_max
215 lloc = spec%user_llocal
221 lmax = spec%user_lmax
227 integer,
intent(in) :: ll
234 integer,
intent(in) :: ll
235 spec%user_llocal = ll
241 integer,
intent(in) :: set_id
246 case (option__pseudopotentialset__standard)
247 filename = trim(conf%share)//
'/pseudopotentials/PSF'
248 case (option__pseudopotentialset__sg15)
249 filename = trim(conf%share)//
'/pseudopotentials/quantum-simulation.org/sg15/'
250 case (option__pseudopotentialset__hgh_lda)
251 filename = trim(conf%share)//
'/pseudopotentials/HGH/lda/'
252 case (option__pseudopotentialset__hgh_lda_sc)
253 filename = trim(conf%share)//
'/pseudopotentials/HGH/lda_sc/'
254 case (option__pseudopotentialset__hscv_lda)
255 filename = trim(conf%share)//
'/pseudopotentials/quantum-simulation.org/hscv/lda/'
256 case (option__pseudopotentialset__hscv_pbe)
257 filename = trim(conf%share)//
'/pseudopotentials/quantum-simulation.org/hscv/pbe/'
258 case (option__pseudopotentialset__pseudodojo_lda)
259 filename = trim(conf%share)//
'/pseudopotentials/pseudo-dojo.org/nc-sr-04_pw_standard/'
260 case (option__pseudopotentialset__pseudodojo_pbe)
261 filename = trim(conf%share)//
'/pseudopotentials/pseudo-dojo.org/nc-sr-05_pbe_standard/'
262 case (option__pseudopotentialset__pseudodojo_pbesol)
263 filename = trim(conf%share)//
'/pseudopotentials/pseudo-dojo.org/nc-sr-04_pbesol_standard/'
264 case (option__pseudopotentialset__pseudodojo_pbe_fr)
265 filename = trim(conf%share)//
'/pseudopotentials/pseudo-dojo.org/nc-fr-04_pbe_standard/'
266 case (option__pseudopotentialset__none)
279 integer,
intent(in) :: set_id
280 type(pseudo_set_t),
intent(in) :: set
281 integer,
intent(out) :: read_data
283 type(element_t) :: el
287 spec%pseudopotential_set_id = set_id
288 spec%pseudopotential_set = set
290 call element_init(el, get_symbol(spec%get_label()))
292 if (spec%pseudopotential_set_id /= option__pseudopotentialset__none .and. pseudo_set_has(spec%pseudopotential_set, el))
then
293 call spec%set_filename(pseudo_set_file_path(spec%pseudopotential_set, el))
296 if (spec%get_z() < 0)
call spec%set_z(real(element_atomic_number(el), real64))
297 if (spec%user_lmax == invalid_l) spec%user_lmax = pseudo_set_lmax(spec%pseudopotential_set, el)
298 if (spec%user_llocal == invalid_l) spec%user_llocal = pseudo_set_llocal(spec%pseudopotential_set, el)
299 if (spec%get_mass() < 0)
call spec%set_mass(element_mass(el))
300 if (spec%get_vdw_radius() < 0)
call spec%set_vdw_radius(element_vdw_radius(el))
314 integer,
intent(in) :: iunit
315 integer,
intent(inout) :: read_data
318 character(len=LABEL_LEN) :: label
319 character(len=MAX_PATH_LEN) :: filename
320 type(element_t) :: element
321 integer :: lmax, llocal
328 read(iunit,*) label, filename, zz, lmax, llocal
330 call spec%set_filename(trim(conf%share)//
'/pseudopotentials/'//trim(filename))
332 assert(trim(label) == trim(spec%get_label()))
337 call element_init(element, get_symbol(label))
339 assert(element_valid(element))
342 if (spec%get_z() < 0)
call spec%set_z(zz)
343 if (spec%get_z() < 0)
call spec%set_z(real(element_atomic_number(element), real64))
344 if (spec%user_lmax == invalid_l) spec%user_lmax = lmax
345 if (spec%user_llocal == invalid_l) spec%user_llocal = llocal
346 if (spec%get_mass() < 0)
call spec%set_mass(element_mass(element))
347 if (spec%get_vdw_radius() < 0)
call spec%set_vdw_radius(element_vdw_radius(element))
349 call element_end(element)
359 integer,
intent(in) :: np
360 real(real64),
contiguous,
intent(in) :: x(:,:)
361 real(real64),
contiguous,
intent(in) :: r(:)
362 integer,
intent(in) :: l, lm, i
363 real(real64),
contiguous,
intent(out) :: uv(:)
371 call lalg_copy(np, r, uv)
372 call spline_eval_vec(spec%ps%kb(l, i), np, uv)
375 call ylmr_real(x(1:3, ip), l, lm, ylm)
376 uv(ip) = uv(ip) * ylm
388 integer,
intent(in) :: np
389 real(real64),
intent(in) :: x(:,:)
390 real(real64),
intent(in) :: r(:)
391 integer,
intent(in) :: l, lm, i
392 complex(real64),
intent(out) :: uv(:)
395 complex(real64) :: ylm
401 call spline_eval_vec(spec%ps%kb(l, i), np, uv)
404 call ylmr_cmplx(x(1:3, ip), l, lm, ylm)
405 uv(ip) = uv(ip) * ylm
416 type(namespace_t),
intent(in) :: namespace
417 integer,
intent(in) :: nspin
418 integer,
intent(in) :: dim
420 integer :: is, n, i, l, m
426 do i = 1, spec%ps%conf%p
427 if (n > spec%niwfs)
exit
428 l = spec%ps%conf%l(i)
430 if (.not. spec%ps%bound(i,is)) cycle
433 spec%iwf_i(n, is) = i
434 spec%iwf_n(n, is) = spec%ps%conf%n(i)
435 spec%iwf_l(n, is) = l
436 spec%iwf_m(n, is) = m
437 spec%iwf_j(n) = spec%ps%conf%j(i)
451 integer,
intent(in) :: ii
452 integer,
intent(in) :: is
453 real(real64),
optional,
intent(in) :: threshold
455 real(real64) threshold_
459 threshold_ = optional_default(threshold, spec%ps%projectors_sphere_threshold)
460 assert(ii <= spec%ps%conf%p)
461 radius = spline_x_threshold(spec%ps%ur(ii, is), threshold_)
473 if (spec%ps%lmax == 0 .and. spec%ps%llocal == 0) is_local = .
true.
485 type(namespace_t),
intent(in) :: namespace
486 real(real64),
intent(in) :: grid_cutoff
487 integer,
intent(in) :: filter
489 character(len=256) :: dirname
491 real(real64) :: local_radius, orbital_radius
495 call ps_separate(this%ps)
497 call ps_getradius(this%ps)
499 if (filter /= ps_filter_none .and. this%ps%projector_type /= proj_hgh)
then
500 call ps_filter(this%ps, filter, grid_cutoff)
501 call ps_getradius(this%ps)
504 call ps_derivatives(this%ps)
506 local_radius = this%ps%vl%x_threshold
508 orbital_radius = m_zero
510 do iorb = 1, this%get_niwfs()
511 orbital_radius = max(orbital_radius, this%get_iwf_radius(this%iwf_i(iorb, 1), is = 1))
514 call messages_write(
'Info: Pseudopotential for '//trim(this%get_label()), new_line = .
true.)
515 call messages_write(
' Radii for localized parts:', new_line = .
true.)
516 call messages_write(
' local part = ')
517 call messages_write(local_radius, fmt =
'f5.1', units = units_out%length, new_line = .
true.)
518 call messages_write(
' non-local part = ')
519 call messages_write(this%ps%rc_max, fmt =
'f5.1', units = units_out%length, new_line = .
true.)
520 call messages_write(
' orbitals = ')
521 call messages_write(orbital_radius, fmt =
'f5.1', units = units_out%length, new_line = .
true.)
522 call messages_info(namespace=namespace)
524 if (max(local_radius, this%ps%rc_max) > 6.0_real64)
then
525 call messages_write(
"One of the radii of your pseudopotential's localized parts seems", new_line = .
true.)
526 call messages_write(
"unusually large; check that your pseudopotential is correct.")
527 call messages_warning(namespace=namespace)
530 if (orbital_radius > 20.0_real64)
then
531 call messages_write(
"The radius of the atomic orbitals given by your pseudopotential seems", new_line = .
true.)
532 call messages_write(
"unusually large; check that your pseudopotential is correct.")
533 call messages_warning(namespace=namespace)
537 write(dirname,
'(a)')
'debug/geometry'
538 call io_mkdir(dirname, namespace)
539 call this%debug(trim(dirname), namespace, grid_cutoff)
548 character(len=*),
intent(in) :: dir
549 type(namespace_t),
intent(in) :: namespace
550 real(real64),
intent(in) :: gmax
553 character(len=256) :: dirname
556 if (.not. mpi_grp_is_root(mpi_world))
then
562 dirname = trim(dir)//
'/'//trim(spec%get_label())
564 call io_mkdir(dirname, namespace)
566 iunit = io_open(trim(dirname)//
'/info', namespace, action=
'write')
568 write(iunit,
'(a,i3)')
'Index = ', spec%get_index()
569 write(iunit,
'(2a)')
'Label = ', trim(spec%get_label())
571 write(iunit,
'(a,f15.2)')
'z = ', spec%get_z()
572 write(iunit,
'(a,f15.2)')
'z_val = ', spec%get_zval()
573 write(iunit,
'(a,f15.2)')
'mass = ', spec%get_mass()
574 write(iunit,
'(a,f15.2)')
'vdw_radius = ', spec%get_vdw_radius()
575 write(iunit,
'(a,l1)')
'local = ', spec%is_local()
576 write(iunit,
'(a,l1)')
'nlcc = ', spec%nlcc
577 write(iunit,
'(a,i3)')
'hubbard_l = ', spec%get_hubbard_l()
578 write(iunit,
'(a,f15.2)')
'hubbard_U = ', spec%get_hubbard_U()
579 write(iunit,
'(a,f15.2)')
'hubbard_j = ', spec%get_hubbard_j()
580 write(iunit,
'(a,f15.2)')
'hubbard_alpha = ', spec%get_hubbard_alpha()
582 if (debug%info)
call ps_debug(spec%ps, trim(dirname), namespace, gmax)
591 type(namespace_t),
intent(in) :: namespace
592 integer,
intent(in) :: ispin
593 integer,
intent(in) :: dim
594 logical,
optional,
intent(in) :: print_info
596 logical :: print_info_
600 print_info_ = optional_default(print_info, .
true.)
603 call spec%set_mass(units_to_atomic(unit_amu, spec%get_mass()))
605 spec%has_density = .false.
608 safe_allocate(spec%ps)
610 call ps_pspio_init(spec%ps, namespace, spec%get_label(), spec%get_z(), spec%user_lmax, &
611 ispin, spec%get_filename())
613 call ps_init(spec%ps, namespace, spec%get_label(), spec%get_z(), spec%user_lmax, &
614 spec%user_llocal, ispin, spec%get_filename())
616 call spec%set_zval(spec%ps%z_val)
617 spec%nlcc = spec%ps%nlcc
618 spec%niwfs = ps_bound_niwfs(spec%ps)
621 spec%user_lmax = invalid_l
622 spec%user_llocal = invalid_l
624 safe_allocate(spec%iwf_n(1:spec%niwfs, 1:ispin))
625 safe_allocate(spec%iwf_l(1:spec%niwfs, 1:ispin))
626 safe_allocate(spec%iwf_m(1:spec%niwfs, 1:ispin))
627 safe_allocate(spec%iwf_i(1:spec%niwfs, 1:ispin))
628 safe_allocate(spec%iwf_j(1:spec%niwfs))
630 call spec%iwf_fix_qn(namespace, ispin, dim)
637 logical pure function pseudopotential_is_ps(this)
645 logical pure function pseudopotential_is_ps_with_nlcc(this)
653 logical pure function pseudopotential_represents_real_atom(spec)
This module is intended to contain "only mathematical" functions and procedures.
subroutine, public ps_end(ps)
integer, parameter, public invalid_l
real(real64) pure function pseudopotential_get_radius(spec)
Return radius of the pseudopotential if this is a pseudo, zero otherwise.
subroutine pseudopotential_build(spec, namespace, ispin, dim, print_info)
subroutine, public pseudopotential_nl_projector(spec, np, x, r, l, lm, i, uV)
This routine returns the non-local projector, built using spherical harmonics.
logical pure function pseudopotential_has_nlcc(spec)
integer pure function pseudopotential_get_user_lmax(spec)
logical pure function pseudopotential_is_ps_with_nlcc(this)
Is the species a pseudopotential derived class or not with nlcc.
subroutine pseudopotential_iwf_fix_qn(spec, namespace, nspin, dim)
set up quantum numbers of orbitals
subroutine pseudopotential_finalize(spec)
real(real64) function pseudopotential_get_iwf_radius(spec, ii, is, threshold)
Return radius outside which orbital is less than threshold value 0.001.
character(len=max_path_len) function, public get_set_directory(set_id)
subroutine, public read_from_default_file(iunit, read_data, spec)
subroutine, public read_from_set(spec, set_id, set, read_data)
Creates a pseudopotential type from a set.
logical pure function pseudopotential_represents_real_atom(spec)
Is the species representing an atomic species or not.
logical function pseudopotential_is_local(spec)
integer pure function pseudopotential_x_functional(spec)
subroutine pseudopotential_init_potential(this, namespace, grid_cutoff, filter)
This routine performs some operations on the pseudopotential functions (filtering,...
integer pure function pseudopotential_c_functional(spec)
logical pure function pseudopotential_is_ps(this)
Is the species a pseudopotential derived class or not.
integer pure function pseudopotential_get_user_lloc(spec)
pure subroutine pseudopotential_set_user_lloc(spec, ll)
pure subroutine pseudopotential_set_user_lmax(spec, ll)
subroutine, public pseudopotential_real_nl_projector(spec, np, x, r, l, lm, i, uV)
This routine returns the non-local projector and its derivative, built using real spherical harmonics...
subroutine pseudopotential_debug(spec, dir, namespace, gmax)
class(pseudopotential_t) function, pointer pseudopotential_constructor(label, index)
The factory routine (or constructor) allocates a pointer of the corresponding type and then calls the...
integer, parameter, public species_pspio
pseudopotential parsed by pspio library
subroutine, public species_end(species)
subroutine, public species_init(this, label, index)
Initializes a species object. This should be the first routine to be called (before species_read and ...
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
This module defines the unit system, used for input and output.
An abstract class for species. Derived classes include jellium, all electron, and pseudopotential spe...