50 integer,
public,
parameter :: &
51 SPECIES_PSEUDO = 7, & !< pseudopotential
60 type(ps_t),
allocatable,
public :: ps
63 integer :: user_llocal
65 integer,
public :: pseudopotential_set_id
66 logical,
public :: pseudopotential_set_initialized
67 type(pseudo_set_t),
public :: pseudopotential_set
91 procedure pseudopotential_constructor
103 class(pseudopotential_t),
pointer :: spec
104 character(len=*),
intent(in) :: label
105 integer,
intent(in) :: index
120 spec%pseudopotential_set_id = option__pseudopotentialset__none
121 spec%pseudopotential_set_initialized = .false.
130 type(pseudopotential_t),
intent(inout) :: spec
134 if (spec%pseudopotential_set_initialized)
then
136 spec%pseudopotential_set_initialized = .false.
139 if (
allocated(spec%ps))
call ps_end(spec%ps)
140 safe_deallocate_a(spec%ps)
148 logical pure function pseudopotential_has_nlcc(spec)
149 class(pseudopotential_t),
intent(in) :: spec
150 pseudopotential_has_nlcc = spec%nlcc
162 select case (spec%pseudopotential_set_id)
164 option__pseudopotentialset__standard, &
165 option__pseudopotentialset__hgh_lda, &
166 option__pseudopotentialset__hgh_lda_sc, &
167 option__pseudopotentialset__hscv_lda)
171 case (option__pseudopotentialset__hscv_pbe)
187 select case (spec%pseudopotential_set_id)
189 option__pseudopotentialset__standard, &
190 option__pseudopotentialset__hgh_lda, &
191 option__pseudopotentialset__hgh_lda_sc, &
192 option__pseudopotentialset__hscv_lda)
196 case (option__pseudopotentialset__hscv_pbe)
208 radius = spec%ps%rc_max
214 lloc = spec%user_llocal
220 lmax = spec%user_lmax
226 integer,
intent(in) :: ll
233 integer,
intent(in) :: ll
234 spec%user_llocal = ll
240 integer,
intent(in) :: set_id
245 case (option__pseudopotentialset__standard)
246 filename = trim(conf%share)//
'/pseudopotentials/PSF'
247 case (option__pseudopotentialset__sg15)
248 filename = trim(conf%share)//
'/pseudopotentials/quantum-simulation.org/sg15/'
249 case (option__pseudopotentialset__hgh_lda)
250 filename = trim(conf%share)//
'/pseudopotentials/HGH/lda/'
251 case (option__pseudopotentialset__hgh_lda_sc)
252 filename = trim(conf%share)//
'/pseudopotentials/HGH/lda_sc/'
253 case (option__pseudopotentialset__hscv_lda)
254 filename = trim(conf%share)//
'/pseudopotentials/quantum-simulation.org/hscv/lda/'
255 case (option__pseudopotentialset__hscv_pbe)
256 filename = trim(conf%share)//
'/pseudopotentials/quantum-simulation.org/hscv/pbe/'
257 case (option__pseudopotentialset__pseudodojo_lda)
258 filename = trim(conf%share)//
'/pseudopotentials/pseudo-dojo.org/nc-sr-04_pw_standard/'
259 case (option__pseudopotentialset__pseudodojo_pbe)
260 filename = trim(conf%share)//
'/pseudopotentials/pseudo-dojo.org/nc-sr-04_pbe_standard/'
261 case (option__pseudopotentialset__pseudodojo_pbesol)
262 filename = trim(conf%share)//
'/pseudopotentials/pseudo-dojo.org/nc-sr-04_pbesol_standard/'
263 case (option__pseudopotentialset__pseudodojo_pbe_fr)
264 filename = trim(conf%share)//
'/pseudopotentials/pseudo-dojo.org/nc-fr-04_pbe_standard/'
265 case (option__pseudopotentialset__none)
278 integer,
intent(in) :: set_id
279 type(pseudo_set_t),
intent(in) :: set
280 integer,
intent(out) :: read_data
282 type(element_t) :: el
286 spec%pseudopotential_set_id = set_id
287 spec%pseudopotential_set = set
289 call element_init(el, get_symbol(spec%get_label()))
291 if (spec%pseudopotential_set_id /= option__pseudopotentialset__none .and. pseudo_set_has(spec%pseudopotential_set, el))
then
292 call spec%set_filename(pseudo_set_file_path(spec%pseudopotential_set, el))
295 if (spec%get_z() < 0)
call spec%set_z(real(element_atomic_number(el), real64))
296 if (spec%user_lmax == invalid_l) spec%user_lmax = pseudo_set_lmax(spec%pseudopotential_set, el)
297 if (spec%user_llocal == invalid_l) spec%user_llocal = pseudo_set_llocal(spec%pseudopotential_set, el)
298 if (spec%get_mass() < 0)
call spec%set_mass(element_mass(el))
299 if (spec%get_vdw_radius() < 0)
call spec%set_vdw_radius(element_vdw_radius(el))
313 integer,
intent(in) :: iunit
314 integer,
intent(inout) :: read_data
317 character(len=LABEL_LEN) :: label
318 character(len=MAX_PATH_LEN) :: filename
319 type(element_t) :: element
320 integer :: lmax, llocal
327 read(iunit,*) label, filename, zz, lmax, llocal
329 call spec%set_filename(trim(conf%share)//
'/pseudopotentials/'//trim(filename))
331 assert(trim(label) == trim(spec%get_label()))
336 call element_init(element, get_symbol(label))
338 assert(element_valid(element))
341 if (spec%get_z() < 0)
call spec%set_z(zz)
342 if (spec%get_z() < 0)
call spec%set_z(real(element_atomic_number(element), real64))
343 if (spec%user_lmax == invalid_l) spec%user_lmax = lmax
344 if (spec%user_llocal == invalid_l) spec%user_llocal = llocal
345 if (spec%get_mass() < 0)
call spec%set_mass(element_mass(element))
346 if (spec%get_vdw_radius() < 0)
call spec%set_vdw_radius(element_vdw_radius(element))
348 call element_end(element)
358 integer,
intent(in) :: np
359 real(real64),
intent(in) :: x(:,:)
360 real(real64),
intent(in) :: r(:)
361 integer,
intent(in) :: l, lm, i
362 real(real64),
intent(out) :: uv(:)
371 call spline_eval_vec(spec%ps%kb(l, i), np, uv)
374 call ylmr_real(x(1:3, ip), l, lm, ylm)
375 uv(ip) = uv(ip) * ylm
387 integer,
intent(in) :: np
388 real(real64),
intent(in) :: x(:,:)
389 real(real64),
intent(in) :: r(:)
390 integer,
intent(in) :: l, lm, i
391 complex(real64),
intent(out) :: uv(:)
394 complex(real64) :: ylm
400 call spline_eval_vec(spec%ps%kb(l, i), np, uv)
403 call ylmr_cmplx(x(1:3, ip), l, lm, ylm)
404 uv(ip) = uv(ip) * ylm
415 type(namespace_t),
intent(in) :: namespace
416 integer,
intent(in) :: nspin
417 integer,
intent(in) :: dim
419 integer :: is, n, i, l, m
425 do i = 1, spec%ps%conf%p
426 if (n > spec%niwfs)
exit
427 l = spec%ps%conf%l(i)
429 if (.not. spec%ps%bound(i,is)) cycle
432 spec%iwf_i(n, is) = i
433 spec%iwf_n(n, is) = spec%ps%conf%n(i)
434 spec%iwf_l(n, is) = l
435 spec%iwf_m(n, is) = m
436 spec%iwf_j(n) = spec%ps%conf%j(i)
450 integer,
intent(in) :: ii
451 integer,
intent(in) :: is
452 real(real64),
optional,
intent(in) :: threshold
454 real(real64) threshold_
458 threshold_ = optional_default(threshold, spec%ps%projectors_sphere_threshold)
459 assert(ii <= spec%ps%conf%p)
460 radius = spline_x_threshold(spec%ps%ur(ii, is), threshold_)
472 if (spec%ps%lmax == 0 .and. spec%ps%llocal == 0) is_local = .
true.
484 type(namespace_t),
intent(in) :: namespace
485 real(real64),
intent(in) :: grid_cutoff
486 integer,
intent(in) :: filter
488 character(len=256) :: dirname
490 real(real64) :: local_radius, orbital_radius
494 call ps_separate(this%ps)
496 call ps_getradius(this%ps)
498 if (filter /= ps_filter_none .and. this%ps%projector_type /= proj_hgh)
then
499 call ps_filter(this%ps, filter, grid_cutoff)
500 call ps_getradius(this%ps)
503 call ps_derivatives(this%ps)
505 local_radius = this%ps%vl%x_threshold
507 orbital_radius = m_zero
509 do iorb = 1, this%get_niwfs()
510 orbital_radius = max(orbital_radius, this%get_iwf_radius(this%iwf_i(iorb, 1), is = 1))
513 call messages_write(
'Info: Pseudopotential for '//trim(this%get_label()), new_line = .
true.)
514 call messages_write(
' Radii for localized parts:', new_line = .
true.)
515 call messages_write(
' local part = ')
516 call messages_write(local_radius, fmt =
'f5.1', units = units_out%length, new_line = .
true.)
517 call messages_write(
' non-local part = ')
518 call messages_write(this%ps%rc_max, fmt =
'f5.1', units = units_out%length, new_line = .
true.)
519 call messages_write(
' orbitals = ')
520 call messages_write(orbital_radius, fmt =
'f5.1', units = units_out%length, new_line = .
true.)
521 call messages_info(namespace=namespace)
523 if (max(local_radius, this%ps%rc_max) > 6.0_real64)
then
524 call messages_write(
"One of the radii of your pseudopotential's localized parts seems", new_line = .
true.)
525 call messages_write(
"unusually large; check that your pseudopotential is correct.")
526 call messages_warning(namespace=namespace)
529 if (orbital_radius > 20.0_real64)
then
530 call messages_write(
"The radius of the atomic orbitals given by your pseudopotential seems", new_line = .
true.)
531 call messages_write(
"unusually large; check that your pseudopotential is correct.")
532 call messages_warning(namespace=namespace)
536 write(dirname,
'(a)')
'debug/geometry'
537 call io_mkdir(dirname, namespace)
538 call this%debug(trim(dirname), namespace, grid_cutoff)
547 character(len=*),
intent(in) :: dir
548 type(namespace_t),
intent(in) :: namespace
549 real(real64),
intent(in) :: gmax
552 character(len=256) :: dirname
555 if (.not. mpi_grp_is_root(mpi_world))
then
561 dirname = trim(dir)//
'/'//trim(spec%get_label())
563 call io_mkdir(dirname, namespace)
565 iunit = io_open(trim(dirname)//
'/info', namespace, action=
'write')
567 write(iunit,
'(a,i3)')
'Index = ', spec%get_index()
568 write(iunit,
'(2a)')
'Label = ', trim(spec%get_label())
570 write(iunit,
'(a,f15.2)')
'z = ', spec%get_z()
571 write(iunit,
'(a,f15.2)')
'z_val = ', spec%get_zval()
572 write(iunit,
'(a,f15.2)')
'mass = ', spec%get_mass()
573 write(iunit,
'(a,f15.2)')
'vdw_radius = ', spec%get_vdw_radius()
574 write(iunit,
'(a,l1)')
'local = ', spec%is_local()
575 write(iunit,
'(a,l1)')
'nlcc = ', spec%nlcc
576 write(iunit,
'(a,i3)')
'hubbard_l = ', spec%get_hubbard_l()
577 write(iunit,
'(a,f15.2)')
'hubbard_U = ', spec%get_hubbard_U()
578 write(iunit,
'(a,f15.2)')
'hubbard_j = ', spec%get_hubbard_j()
579 write(iunit,
'(a,f15.2)')
'hubbard_alpha = ', spec%get_hubbard_alpha()
581 if (debug%info)
call ps_debug(spec%ps, trim(dirname), namespace, gmax)
590 type(namespace_t),
intent(in) :: namespace
591 integer,
intent(in) :: ispin
592 integer,
intent(in) :: dim
593 logical,
optional,
intent(in) :: print_info
595 logical :: print_info_
599 print_info_ = optional_default(print_info, .
true.)
602 call spec%set_mass(units_to_atomic(unit_amu, spec%get_mass()))
604 spec%has_density = .false.
607 safe_allocate(spec%ps)
609 call ps_pspio_init(spec%ps, namespace, spec%get_label(), spec%get_z(), spec%user_lmax, &
610 ispin, spec%get_filename())
612 call ps_init(spec%ps, namespace, spec%get_label(), spec%get_z(), spec%user_lmax, &
613 spec%user_llocal, ispin, spec%get_filename())
615 call spec%set_zval(spec%ps%z_val)
616 spec%nlcc = spec%ps%nlcc
617 spec%niwfs = ps_bound_niwfs(spec%ps)
620 spec%user_lmax = invalid_l
621 spec%user_llocal = invalid_l
623 safe_allocate(spec%iwf_n(1:spec%niwfs, 1:ispin))
624 safe_allocate(spec%iwf_l(1:spec%niwfs, 1:ispin))
625 safe_allocate(spec%iwf_m(1:spec%niwfs, 1:ispin))
626 safe_allocate(spec%iwf_i(1:spec%niwfs, 1:ispin))
627 safe_allocate(spec%iwf_j(1:spec%niwfs))
629 call spec%iwf_fix_qn(namespace, ispin, dim)
636 logical pure function pseudopotential_is_ps(this)
644 logical pure function pseudopotential_is_ps_with_nlcc(this)
652 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...