24 use,
intrinsic :: iso_fortran_env
40 logical :: initialized
41 logical :: kleinman_bylander
43 logical :: has_density
44 integer :: atomic_number
46 real(real64) :: valence_charge
51 integer :: nwavefunctions
52 real(real64),
allocatable :: grid(:)
53 real(real64),
allocatable :: weights(:)
54 real(real64),
allocatable :: potential(:, :)
55 real(real64),
allocatable :: wavefunction(:, :)
56 real(real64),
allocatable :: projector(:, :, :)
57 real(real64),
allocatable :: dij(:, :, :)
58 real(real64),
allocatable :: nlcc_density(:)
59 real(real64),
allocatable :: density(:)
60 integer,
allocatable :: wf_n(:)
61 integer,
allocatable :: wf_l(:)
62 real(real64),
allocatable :: wf_occ(:)
63 type(pseudo_t) :: pseudo
69 subroutine ps_xml_init(this, namespace, filename, fmt, ierr)
70 type(ps_xml_t),
intent(inout) :: this
71 type(namespace_t),
intent(in) :: namespace
72 character(len=*),
intent(in) :: filename
73 integer,
intent(in) :: fmt
74 integer,
intent(out) :: ierr
76 integer :: ll, ii, ic, jc
77 type(pseudo_t) :: pseudo
81 this%initialized = .false.
86 call messages_write(
"Pseudopotential file '" // trim(filename) //
"' not found")
91 call messages_write(
"Cannot determine the format for pseudopotential file '" // trim(filename) //
"'")
96 call messages_write(
"Ultrasoft pseudopotential file '" // trim(filename) //
"' not supported")
101 call messages_write(
"PAW pseudopotential file '" // trim(filename) //
"' not supported")
106 call messages_write(
"Pseudopotential file '" // trim(filename) //
"' not supported")
111 call messages_write(
"Pseudopotential file '" // trim(filename) //
"' not supported")
115 this%initialized = .
true.
126 safe_allocate(this%grid(1:this%grid_size))
127 safe_allocate(this%weights(1:this%grid_size))
132 if (.not. this%kleinman_bylander)
then
134 safe_allocate(this%potential(1:this%grid_size, 0:this%lmax))
135 safe_allocate(this%wavefunction(1:this%grid_size, 0:this%lmax))
139 call messages_write(
"The pseudopotential file '"//trim(filename)//
"' does not contain")
151 safe_allocate(this%potential(1:this%grid_size, this%llocal:this%llocal))
155 safe_allocate(this%projector(1:this%grid_size, 0:this%lmax, 1:this%nchannels))
157 this%projector = 0.0_real64
160 if (this%llocal == ll) cycle
166 safe_allocate(this%dij(0:this%lmax, 1:this%nchannels, 1:this%nchannels))
168 this%dij = 0.0_real64
170 do ic = 1, this%nchannels
171 do jc = 1, this%nchannels
172 this%dij(ll, ic, jc) =
pseudo_dij(pseudo, ll, ic, jc)
179 safe_allocate(this%wavefunction(1:this%grid_size, 1:this%nwavefunctions))
180 safe_allocate(this%wf_n(1:this%nwavefunctions))
181 safe_allocate(this%wf_l(1:this%nwavefunctions))
182 safe_allocate(this%wf_occ(1:this%nwavefunctions))
184 do ii = 1, this%nwavefunctions
185 call pseudo_wavefunction(pseudo, ii, this%wf_n(ii), this%wf_l(ii), this%wf_occ(ii), this%wavefunction(1, ii))
192 if (this%has_density)
then
193 safe_allocate(this%density(1:this%grid_size))
199 safe_allocate(this%nlcc_density(1:this%grid_size))
217 real(real64) :: nrm, rr
224 do ip = 1, this%grid_size
226 nrm = nrm + this%wavefunction(ip, ll)**2*this%weights(ip)*rr**2
230 nrm = abs(nrm -
m_one)
231 if (nrm > 1.0e-5_real64)
then
232 write(
message(1),
'(a,i2,a)')
"Eigenstate for l = ", ll,
' is not normalized'
233 write(
message(2),
'(a, f12.6,a)')
'(abs(1 - norm) = ', nrm,
')'
244 type(
ps_xml_t),
intent(inout) :: this
250 safe_deallocate_a(this%grid)
251 safe_deallocate_a(this%weights)
252 safe_deallocate_a(this%potential)
253 safe_deallocate_a(this%wavefunction)
254 safe_deallocate_a(this%projector)
255 safe_deallocate_a(this%dij)
256 safe_deallocate_a(this%nlcc_density)
257 safe_deallocate_a(this%density)
258 safe_deallocate_a(this%wf_n)
259 safe_deallocate_a(this%wf_l)
260 safe_deallocate_a(this%wf_occ)
double sqrt(double __x) __attribute__((__nothrow__
real(real64), parameter, public m_zero
real(real64), parameter, public m_one
subroutine, public messages_warning(no_lines, all_nodes, namespace)
subroutine, public messages_new_line()
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
subroutine, public ps_xml_end(this)
subroutine ps_xml_check_normalization(this, namespace)
checks normalization of the pseudo wavefunctions
subroutine, public ps_xml_init(this, namespace, filename, fmt, ierr)
integer, parameter, public pseudo_type_kleinman_bylander
logical function, public pseudo_has_radial_function(pseudo, l)
logical function, public pseudo_has_density(pseudo)
logical function, public pseudo_has_nlcc(pseudo)
integer, parameter, public pseudo_status_unsupported_type
integer, parameter, public pseudo_status_file_not_found
integer, parameter, public pseudo_status_unsupported_type_ultrasoft
integer, parameter, public pseudo_status_unknown_format
integer, parameter, public pseudo_status_unsupported_type_paw
integer, parameter, public pseudo_status_format_not_supported