Octopus
pseudopotential.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-2006 M. Marques, A. Castro, A. Rubio, G. Bertsch
2!! Copyright (C) 2023-2024 N. Tancogne-Dejean
3!!
4!! This program is free software; you can redistribute it and/or modify
5!! it under the terms of the GNU General Public License as published by
6!! the Free Software Foundation; either version 2, or (at your option)
7!! any later version.
8!!
9!! This program is distributed in the hope that it will be useful,
10!! but WITHOUT ANY WARRANTY; without even the implied warranty of
11!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12!! GNU General Public License for more details.
13!!
14!! You should have received a copy of the GNU General Public License
15!! along with this program; if not, write to the Free Software
16!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
17!! 02110-1301, USA.
18!!
19#include "global.h"
20
22 use debug_oct_m
24 use global_oct_m
25 use io_oct_m
27 use math_oct_m
29 use mpi_oct_m
32 use ps_oct_m
33 use pseudo_oct_m
37 use unit_oct_m
39
40 implicit none
41
42 private
43 public :: &
50
51 integer, public, parameter :: &
52 SPECIES_PSEUDO = 7, & !< pseudopotential
53 species_pspio = 110
54
55 type, extends(species_t) :: pseudopotential_t
56 private
57
58 integer :: type
59 logical :: nlcc
60
61 type(ps_t), allocatable, public :: ps
62
63 integer :: user_lmax
64 integer :: user_llocal
65
66 integer, public :: pseudopotential_set_id
67 logical, public :: pseudopotential_set_initialized
68 type(pseudo_set_t), public :: pseudopotential_set
69
70 contains
71 procedure :: has_nlcc => pseudopotential_has_nlcc
72 procedure :: x_functional => pseudopotential_x_functional
73 procedure :: c_functional => pseudopotential_c_functional
74 procedure :: get_radius => pseudopotential_get_radius
75 procedure :: iwf_fix_qn => pseudopotential_iwf_fix_qn
76 procedure :: get_iwf_radius => pseudopotential_get_iwf_radius
77 procedure :: is_local => pseudopotential_is_local
78 procedure :: debug => pseudopotential_debug
79 procedure :: build => pseudopotential_build
80 procedure :: init_potential => pseudopotential_init_potential
81 procedure :: get_user_lmax => pseudopotential_get_user_lmax
82 procedure :: get_user_lloc => pseudopotential_get_user_lloc
83 procedure :: set_user_lmax => pseudopotential_set_user_lmax
84 procedure :: set_user_lloc => pseudopotential_set_user_lloc
85 procedure :: is_ps => pseudopotential_is_ps
86 procedure :: is_ps_with_nlcc => pseudopotential_is_ps_with_nlcc
87 procedure :: represents_real_atom => pseudopotential_represents_real_atom
89 end type pseudopotential_t
90
91 interface pseudopotential_t
92 procedure pseudopotential_constructor
93 end interface pseudopotential_t
94
95contains
96
97
98 ! ---------------------------------------------------------
103 function pseudopotential_constructor(label, index) result(spec)
104 class(pseudopotential_t), pointer :: spec
105 character(len=*), intent(in) :: label
106 integer, intent(in) :: index
107
109
110 safe_allocate(spec)
111
112 call species_init(spec, label, index)
113
114 spec%nlcc = .false.
115
116 spec%user_lmax = invalid_l
117 spec%user_llocal = invalid_l
118
119 spec%type = 0
120
121 spec%pseudopotential_set_id = option__pseudopotentialset__none
122 spec%pseudopotential_set_initialized = .false.
123 call pseudo_set_nullify(spec%pseudopotential_set)
124
126 end function pseudopotential_constructor
127
128
129 ! ---------------------------------------------------------
130 subroutine pseudopotential_finalize(spec)
131 type(pseudopotential_t), intent(inout) :: spec
132
134
135 if (spec%pseudopotential_set_initialized) then
136 call pseudo_set_end(spec%pseudopotential_set)
137 spec%pseudopotential_set_initialized = .false.
138 end if
139
140 if (allocated(spec%ps)) call ps_end(spec%ps)
141 safe_deallocate_a(spec%ps)
142
143 call species_end(spec)
146 end subroutine pseudopotential_finalize
147
148 ! ---------------------------------------------------------
149 logical pure function pseudopotential_has_nlcc(spec)
150 class(pseudopotential_t), intent(in) :: spec
151 pseudopotential_has_nlcc = spec%nlcc
153
155 ! ---------------------------------------------------------
156 integer pure function pseudopotential_x_functional(spec)
157 class(pseudopotential_t), intent(in) :: spec
158
159 pseudopotential_x_functional = spec%ps%exchange_functional
161 ! if we do not know, try the pseudpotential set
162 if (pseudopotential_x_functional == pseudo_exchange_unknown) then
163 select case (spec%pseudopotential_set_id)
164 case ( &
165 option__pseudopotentialset__standard, &
166 option__pseudopotentialset__hgh_lda, &
167 option__pseudopotentialset__hgh_lda_sc, &
168 option__pseudopotentialset__hscv_lda)
170 pseudopotential_x_functional = option__xcfunctional__lda_x
172 case (option__pseudopotentialset__hscv_pbe)
173 pseudopotential_x_functional = option__xcfunctional__gga_x_pbe
174 end select
175 end if
179 ! ---------------------------------------------------------
181 integer pure function pseudopotential_c_functional(spec)
182 class(pseudopotential_t), intent(in) :: spec
183
184 pseudopotential_c_functional = spec%ps%correlation_functional
185
186 ! if we do not know, try the pseudpotential set
187 if (pseudopotential_c_functional == pseudo_correlation_unknown) then
188 select case (spec%pseudopotential_set_id)
189 case ( &
190 option__pseudopotentialset__standard, &
191 option__pseudopotentialset__hgh_lda, &
192 option__pseudopotentialset__hgh_lda_sc, &
193 option__pseudopotentialset__hscv_lda)
194
195 pseudopotential_c_functional = option__xcfunctional__lda_c_pz_mod/1000
197 case (option__pseudopotentialset__hscv_pbe)
198 pseudopotential_c_functional = option__xcfunctional__gga_c_pbe/1000
199 end select
200 end if
201
203
204 ! ---------------------------------------------------------
206 real(real64) pure function pseudopotential_get_radius(spec) result(radius)
207 class(pseudopotential_t), intent(in) :: spec
208
209 radius = spec%ps%rc_max
210 end function pseudopotential_get_radius
211
212 ! ---------------------------------------------------------
213 integer pure function pseudopotential_get_user_lloc(spec) result(lloc)
214 class(pseudopotential_t), intent(in) :: spec
215 lloc = spec%user_llocal
217
218 ! ---------------------------------------------------------
219 integer pure function pseudopotential_get_user_lmax(spec) result(lmax)
220 class(pseudopotential_t), intent(in) :: spec
221 lmax = spec%user_lmax
224 ! ---------------------------------------------------------
225 pure subroutine pseudopotential_set_user_lmax(spec, ll)
226 class(pseudopotential_t), intent(inout) :: spec
227 integer, intent(in) :: ll
228 spec%user_lmax = ll
229 end subroutine pseudopotential_set_user_lmax
230
231 ! ---------------------------------------------------------
232 pure subroutine pseudopotential_set_user_lloc(spec, ll)
233 class(pseudopotential_t), intent(inout) :: spec
234 integer, intent(in) :: ll
235 spec%user_llocal = ll
236 end subroutine pseudopotential_set_user_lloc
237
238 ! ---------------------------------------------------------
239
240 character(len=MAX_PATH_LEN) function get_set_directory(set_id) result(filename)
241 integer, intent(in) :: set_id
243 push_sub(get_set_directory)
244
245 select case (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)
267 filename = ''
268 case default
269 assert(.false.)
270 end select
271
272 pop_sub(get_set_directory)
273 end function get_set_directory
275 ! ---------------------------------------------------------
277 subroutine read_from_set(spec, set_id, set, read_data)
278 type(pseudopotential_t), intent(inout) :: spec
279 integer, intent(in) :: set_id
280 type(pseudo_set_t), intent(in) :: set
281 integer, intent(out) :: read_data
282
283 type(element_t) :: el
284
285 push_sub(read_from_set)
286
287 spec%pseudopotential_set_id = set_id
288 spec%pseudopotential_set = set
289
290 call element_init(el, get_symbol(spec%get_label()))
291
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))
294
295 ! these might have been set before
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))
301 read_data = 8
302 else
303 read_data = 0
304 end if
305
306 call element_end(el)
307
308 pop_sub(read_from_set)
309 end subroutine read_from_set
310
311
312 ! ---------------------------------------------------------
313 subroutine read_from_default_file(iunit, read_data, spec)
314 integer, intent(in) :: iunit
315 integer, intent(inout) :: read_data
316 class(pseudopotential_t), intent(inout) :: spec
317
318 character(len=LABEL_LEN) :: label
319 character(len=MAX_PATH_LEN) :: filename
320 type(element_t) :: element
321 integer :: lmax, llocal
322 real(real64) :: zz
323
324 push_sub(read_from_default_file)
326 backspace(iunit)
327
328 read(iunit,*) label, filename, zz, lmax, llocal
329
330 call spec%set_filename(trim(conf%share)//'/pseudopotentials/'//trim(filename))
331
332 assert(trim(label) == trim(spec%get_label()))
334 read_data = 8
335
336 ! get the mass, vdw radius and atomic number for this element
337 call element_init(element, get_symbol(label))
338
339 assert(element_valid(element))
340
341 ! these might have been set before
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))
348
349 call element_end(element)
350
352 end subroutine read_from_default_file
353
354 ! ---------------------------------------------------------
357 subroutine pseudopotential_real_nl_projector(spec, np, x, r, l, lm, i, uV)
358 class(pseudopotential_t), intent(in) :: spec
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(:)
364
365 integer :: ip
366 real(real64) :: ylm
367
369
370 if (np > 0) then
371 call lalg_copy(np, r, uv)
372 call spline_eval_vec(spec%ps%kb(l, i), np, uv)
373
374 do ip = 1, np
375 call ylmr_real(x(1:3, ip), l, lm, ylm)
376 uv(ip) = uv(ip) * ylm
377 end do
378 end if
379
382
383 ! ---------------------------------------------------------
386 subroutine pseudopotential_nl_projector(spec, np, x, r, l, lm, i, uV)
387 class(pseudopotential_t), intent(in) :: spec
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(:)
393
394 integer :: ip
395 complex(real64) :: ylm
396
398
399 if (np > 0) then
400 uv(1:np) = r(1:np)
401 call spline_eval_vec(spec%ps%kb(l, i), np, uv)
402
403 do ip = 1, np
404 call ylmr_cmplx(x(1:3, ip), l, lm, ylm)
405 uv(ip) = uv(ip) * ylm
406 end do
407 end if
408
410 end subroutine pseudopotential_nl_projector
411
412 ! ---------------------------------------------------------
414 subroutine pseudopotential_iwf_fix_qn(spec, namespace, nspin, dim)
415 class(pseudopotential_t), intent(inout) :: spec
416 type(namespace_t), intent(in) :: namespace
417 integer, intent(in) :: nspin
418 integer, intent(in) :: dim
419
420 integer :: is, n, i, l, m
421
423
424 do is = 1, nspin
425 n = 1
426 do i = 1, spec%ps%conf%p
427 if (n > spec%niwfs) exit
428 l = spec%ps%conf%l(i)
429
430 if (.not. spec%ps%bound(i,is)) cycle
431
432 do m = -l, l
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)
438 n = n + 1
439 end do
440
441 end do
442 end do
443
445 end subroutine pseudopotential_iwf_fix_qn
446
447 ! ---------------------------------------------------------
449 real(real64) function pseudopotential_get_iwf_radius(spec, ii, is, threshold) result(radius)
450 class(pseudopotential_t), intent(in) :: spec
451 integer, intent(in) :: ii
452 integer, intent(in) :: is
453 real(real64), optional, intent(in) :: threshold
454
455 real(real64) threshold_
456
458
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_)
462
465
466 ! ---------------------------------------------------------
467 logical function pseudopotential_is_local(spec) result(is_local)
468 class(pseudopotential_t), intent(in) :: spec
469
471
472 is_local = .false.
473 if (spec%ps%lmax == 0 .and. spec%ps%llocal == 0) is_local = .true.
474
476 end function pseudopotential_is_local
477
478 ! ---------------------------------------------------------
482 ! ---------------------------------------------------------
483 subroutine pseudopotential_init_potential(this, namespace, grid_cutoff, filter)
484 class(pseudopotential_t), intent(inout) :: this
485 type(namespace_t), intent(in) :: namespace
486 real(real64), intent(in) :: grid_cutoff
487 integer, intent(in) :: filter
488
489 character(len=256) :: dirname
490 integer :: iorb
491 real(real64) :: local_radius, orbital_radius
492
494
495 call ps_separate(this%ps)
496
497 call ps_getradius(this%ps)
498
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) ! radius may have changed
502 end if
503
504 call ps_derivatives(this%ps)
505
506 local_radius = this%ps%vl%x_threshold
508 orbital_radius = m_zero
509 ! FIXME: should take max over spins too here.
510 do iorb = 1, this%get_niwfs()
511 orbital_radius = max(orbital_radius, this%get_iwf_radius(this%iwf_i(iorb, 1), is = 1))
512 end do
513
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)
523
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)
528 end if
529
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)
534 end if
535
536 if (debug%info) then
537 write(dirname, '(a)') 'debug/geometry'
538 call io_mkdir(dirname, namespace)
539 call this%debug(trim(dirname), namespace, grid_cutoff)
540 end if
541
543 end subroutine pseudopotential_init_potential
544
545 ! ---------------------------------------------------------
546 subroutine pseudopotential_debug(spec, dir, namespace, gmax)
547 class(pseudopotential_t), intent(in) :: spec
548 character(len=*), intent(in) :: dir
549 type(namespace_t), intent(in) :: namespace
550 real(real64), intent(in) :: gmax
551
552
553 character(len=256) :: dirname
554 integer :: iunit
555
556 if (.not. mpi_grp_is_root(mpi_world)) then
557 return
558 end if
559
561
562 dirname = trim(dir)//'/'//trim(spec%get_label())
563
564 call io_mkdir(dirname, namespace)
565
566 iunit = io_open(trim(dirname)//'/info', namespace, action='write')
567
568 write(iunit, '(a,i3)') 'Index = ', spec%get_index()
569 write(iunit, '(2a)') 'Label = ', trim(spec%get_label())
570 !write(iunit, '(a,i3)') 'Type = ', spec%type
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()
581
582 if (debug%info) call ps_debug(spec%ps, trim(dirname), namespace, gmax)
583
584 call io_close(iunit)
585 pop_sub(pseudopotential_debug)
586 end subroutine pseudopotential_debug
587
588 ! ---------------------------------------------------------
589 subroutine pseudopotential_build(spec, namespace, ispin, dim, print_info)
590 class(pseudopotential_t), intent(inout) :: spec
591 type(namespace_t), intent(in) :: namespace
592 integer, intent(in) :: ispin
593 integer, intent(in) :: dim
594 logical, optional, intent(in) :: print_info
595
596 logical :: print_info_
597
598 push_sub(pseudopotential_build)
599
600 print_info_ = optional_default(print_info, .true.)
601
602 ! masses are always in amu, so convert them to a.u.
603 call spec%set_mass(units_to_atomic(unit_amu, spec%get_mass()))
604
605 spec%has_density = .false.
606
607 ! allocate structure
608 safe_allocate(spec%ps)
609 if (spec%type == species_pspio) then
610 call ps_pspio_init(spec%ps, namespace, spec%get_label(), spec%get_z(), spec%user_lmax, &
611 ispin, spec%get_filename())
612 else
613 call ps_init(spec%ps, namespace, spec%get_label(), spec%get_z(), spec%user_lmax, &
614 spec%user_llocal, ispin, spec%get_filename())
615 end if
616 call spec%set_zval(spec%ps%z_val)
617 spec%nlcc = spec%ps%nlcc
618 spec%niwfs = ps_bound_niwfs(spec%ps)
619
620 ! invalidate these variables as they should not be used after
621 spec%user_lmax = invalid_l
622 spec%user_llocal = invalid_l
623
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))
629
630 call spec%iwf_fix_qn(namespace, ispin, dim)
631
632 pop_sub(pseudopotential_build)
633 end subroutine pseudopotential_build
634
635 ! ---------------------------------------------------------
637 logical pure function pseudopotential_is_ps(this)
638 class(pseudopotential_t), intent(in) :: this
641 end function pseudopotential_is_ps
642
643 ! ---------------------------------------------------------
645 logical pure function pseudopotential_is_ps_with_nlcc(this)
646 class(pseudopotential_t), intent(in) :: this
647
648 pseudopotential_is_ps_with_nlcc = this%has_nlcc()
650
651 ! ---------------------------------------------------------
653 logical pure function pseudopotential_represents_real_atom(spec)
654 class(pseudopotential_t), intent(in) :: spec
655
658
659
660
661end module pseudopotential_oct_m
662
663!! Local Variables:
664!! mode: f90
665!! coding: utf-8
666!! End:
Delete the C++ object.
Definition: pseudo_set.F90:152
Nullify the C++ pointer.
Definition: pseudo_set.F90:142
Definition: io.F90:114
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:115
Definition: ps.F90:114
subroutine, public ps_end(ps)
Definition: ps.F90:1059
integer, parameter, public invalid_l
Definition: ps.F90:178
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)
Definition: species.F90:498
subroutine, public species_init(this, label, index)
Initializes a species object. This should be the first routine to be called (before species_read and ...
Definition: species.F90:287
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:132
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...
Definition: species.F90:143
int true(void)