24 use,
intrinsic :: iso_fortran_env
42 interface operator(==)
44 end interface operator(==)
47 integer,
public,
parameter :: LABEL_LEN=15
54 character(len=LABEL_LEN) :: label
59 real(real64) :: vdw_radius
61 logical,
public :: has_density
63 character(len=MAX_PATH_LEN) :: filename
66 integer,
public :: niwfs
67 integer,
allocatable,
public :: iwf_l(:, :), iwf_m(:, :), iwf_i(:, :), iwf_n(:, :)
68 real(real64),
allocatable,
public :: iwf_j(:)
71 real(real64) :: hubbard_U
72 real(real64) :: hubbard_j
73 real(real64) :: hubbard_alpha
76 procedure(species_is_local),
deferred :: is_local
77 procedure(species_get_iwf_radius),
deferred :: get_iwf_radius
78 procedure(species_iwf_fix_qn),
deferred :: iwf_fix_qn
79 procedure(species_init_potential),
deferred :: init_potential
80 procedure(species_build),
deferred :: build
83 procedure(species_debug),
deferred :: debug
118 class(species_t),
pointer :: s
125 class(species_t),
intent(in) :: spec
129 real(real64) function species_get_iwf_radius(spec, ii, is, threshold) result(radius)
132 class(species_t),
intent(in) :: spec
133 integer,
intent(in) :: ii
134 integer,
intent(in) :: is
135 real(real64),
optional,
intent(in) :: threshold
143 type(namespace_t),
intent(in) :: namespace
144 integer,
intent(in) :: nspin
145 integer,
intent(in) :: dim
155 real(real64),
intent(in) :: grid_cutoff
156 integer,
intent(in) :: filter
165 integer,
intent(in) :: ispin
166 integer,
intent(in) :: dim
167 logical,
optional,
intent(in) :: print_info
176 character(len=*),
intent(in) :: dir
178 real(real64),
intent(in) :: gmax
195 character(len=*),
intent(in) :: label
196 integer,
intent(in) :: index
200 this%label = trim(label)
206 this%vdw_radius = -m_one
207 this%has_density = .false.
210 this%hubbard_U = m_zero
211 this%hubbard_j = m_zero
212 this%hubbard_alpha = m_zero
221 integer,
intent(in) :: min_niwfs
227 if (
size >= min_niwfs)
exit
234 character(len=LABEL_LEN) pure function species_label(species)
240 character(len=2) pure function species_label_short(species)
246 integer pure function species_index(species)
252 real(real64)
pure function species_zval(species)
259 class(
species_t),
intent(inout) :: species
260 real(real64),
intent(in) :: zval
265 real(real64)
pure function species_z(species)
272 class(
species_t),
intent(inout) :: species
273 real(real64),
intent(in) :: z
285 class(
species_t),
intent(inout) :: species
286 real(real64),
intent(in) :: mass
298 class(
species_t),
intent(inout) :: species
299 real(real64),
intent(in) :: radius
300 species%vdw_radius = radius
335 class(
species_t),
intent(inout) :: species
336 integer,
intent(in) :: hubbard_l
337 species%hubbard_l = hubbard_l
342 class(
species_t),
intent(inout) :: species
343 real(real64),
intent(in) :: hubbard_u
344 species%hubbard_u = hubbard_u
349 class(
species_t),
intent(inout) :: species
350 real(real64),
intent(in) :: hubbard_j
351 species%hubbard_j = hubbard_j
356 class(
species_t),
intent(inout) :: species
357 real(real64),
intent(in) :: hubbard_alpha
358 species%hubbard_alpha = hubbard_alpha
369 class(
species_t),
intent(inout) :: species
370 character(len=*),
intent(in) :: filename
371 species%filename = trim(filename)
377 integer,
intent(in) :: j, is
378 integer,
intent(out) :: i, l, m
380 i = species%iwf_i(j, is)
381 l = species%iwf_l(j, is)
382 m = species%iwf_m(j, is)
388 integer,
intent(in) :: j, is
389 integer,
intent(out) :: n
391 n = species%iwf_n(j, is)
397 integer,
intent(in) :: iorb
398 real(real64),
intent(out) :: j
400 j = species%iwf_j(iorb)
405 class(
species_t),
intent(inout) :: species
409 safe_deallocate_a(species%iwf_n)
410 safe_deallocate_a(species%iwf_l)
411 safe_deallocate_a(species%iwf_m)
412 safe_deallocate_a(species%iwf_i)
413 safe_deallocate_a(species%iwf_j)
424 character(len=LABEL_LEN) function get_symbol(label)
result(symbol)
425 character(len=*),
intent(in) :: label
430 do ilend = 1, len(label)
431 if (iachar(label(ilend:ilend)) >= iachar(
'a') .and. iachar(label(ilend:ilend)) <= iachar(
'z')) cycle
432 if (iachar(label(ilend:ilend)) >= iachar(
'A') .and. iachar(label(ilend:ilend)) <= iachar(
'Z')) cycle
437 symbol = label(1:ilend)
443 class(
species_t),
intent(in) :: spec1, spec2
447 same = same_type_as(spec1, spec2)
448 if(abs(spec1%z - spec2%z) > m_epsilon) same = .false.
449 if(abs(spec1%z_val - spec2%z_val) > m_epsilon) same = .false.
450 if(abs(spec1%mass - spec2%mass) > m_epsilon) same = .false.
real(real64) pure function species_zval(species)
integer, parameter, private libxc_c_index
pure subroutine species_set_hubbard_l(species, hubbard_l)
pure subroutine species_set_z(species, z)
logical pure function species_user_defined(spec)
Is the species user-defined or not.
real(real64) pure function species_hubbard_u(species)
character(len=label_len) function, public get_symbol(label)
character(len=2) pure function, public species_label_short(species)
pure subroutine species_set_hubbard_j(species, hubbard_j)
integer pure function species_niwfs(species)
pure subroutine species_set_filename(species, filename)
pure subroutine species_set_zval(species, zval)
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 ...
logical pure function species_represents_real_atom(spec)
Is the species representing an atomic species or not.
logical pure function species_is_ps_with_nlcc(this)
Is the species a pseudopotential derived class or not with nlcc.
pure subroutine species_set_hubbard_alpha(species, hubbard_alpha)
character(len=200) pure function species_filename(species)
pure subroutine species_iwf_ilm(species, j, is, i, l, m)
logical function species_is_same_species(spec1, spec2)
integer pure function species_index(species)
logical pure function species_is_ps(this)
Is the species a pseudopotential derived class or not.
real(real64) pure function species_vdw_radius(species)
pure subroutine species_iwf_n(species, j, is, n)
pure subroutine species_set_hubbard_u(species, hubbard_u)
real(real64) pure function species_hubbard_alpha(species)
real(real64) pure function species_z(species)
pure subroutine species_iwf_j(species, iorb, j)
character(len=label_len) pure function species_label(species)
logical pure function species_is_full(this)
Is the species an all-electron derived class or not.
pure subroutine species_set_vdw_radius(species, radius)
integer pure function species_hubbard_l(species)
pure subroutine species_set_mass(species, mass)
real(real64) pure function species_hubbard_j(species)
integer pure function, public species_closed_shell_size(min_niwfs)
find size of closed shell for hydrogenic atom with size at least min_niwfs
real(real64) pure function species_mass(species)
An abstract class for species. Derived classes include jellium, all electron, and pseudopotential spe...
Needed for having an array of pointers See for instance https: