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
26 use math_oct_m
28 use mpi_oct_m
31 use ps_oct_m
32 use pseudo_oct_m
36 use unit_oct_m
38
39 implicit none
40
41 private
42 public :: &
49
50 integer, public, parameter :: &
51 SPECIES_PSEUDO = 7, & !< pseudopotential
52 species_pspio = 110
53
54 type, extends(species_t) :: pseudopotential_t
55 private
56
57 integer :: type
58 logical :: nlcc
59
60 type(ps_t), allocatable, public :: ps
61
62 integer :: user_lmax
63 integer :: user_llocal
64
65 integer, public :: pseudopotential_set_id
66 logical, public :: pseudopotential_set_initialized
67 type(pseudo_set_t), public :: pseudopotential_set
68
69 contains
70 procedure :: has_nlcc => pseudopotential_has_nlcc
71 procedure :: x_functional => pseudopotential_x_functional
72 procedure :: c_functional => pseudopotential_c_functional
73 procedure :: get_radius => pseudopotential_get_radius
74 procedure :: iwf_fix_qn => pseudopotential_iwf_fix_qn
75 procedure :: get_iwf_radius => pseudopotential_get_iwf_radius
76 procedure :: is_local => pseudopotential_is_local
77 procedure :: debug => pseudopotential_debug
78 procedure :: build => pseudopotential_build
79 procedure :: init_potential => pseudopotential_init_potential
80 procedure :: get_user_lmax => pseudopotential_get_user_lmax
81 procedure :: get_user_lloc => pseudopotential_get_user_lloc
82 procedure :: set_user_lmax => pseudopotential_set_user_lmax
83 procedure :: set_user_lloc => pseudopotential_set_user_lloc
84 procedure :: is_ps => pseudopotential_is_ps
85 procedure :: is_ps_with_nlcc => pseudopotential_is_ps_with_nlcc
86 procedure :: represents_real_atom => pseudopotential_represents_real_atom
88 end type pseudopotential_t
89
90 interface pseudopotential_t
91 procedure pseudopotential_constructor
92 end interface pseudopotential_t
93
94contains
95
96
97 ! ---------------------------------------------------------
102 function pseudopotential_constructor(label, index) result(spec)
103 class(pseudopotential_t), pointer :: spec
104 character(len=*), intent(in) :: label
105 integer, intent(in) :: index
106
108
109 safe_allocate(spec)
110
111 call species_init(spec, label, index)
112
113 spec%nlcc = .false.
115 spec%user_lmax = invalid_l
116 spec%user_llocal = invalid_l
117
118 spec%type = 0
119
120 spec%pseudopotential_set_id = option__pseudopotentialset__none
121 spec%pseudopotential_set_initialized = .false.
122 call pseudo_set_nullify(spec%pseudopotential_set)
123
125 end function pseudopotential_constructor
126
127
128 ! ---------------------------------------------------------
129 subroutine pseudopotential_finalize(spec)
130 type(pseudopotential_t), intent(inout) :: spec
131
133
134 if (spec%pseudopotential_set_initialized) then
135 call pseudo_set_end(spec%pseudopotential_set)
136 spec%pseudopotential_set_initialized = .false.
137 end if
138
139 if (allocated(spec%ps)) call ps_end(spec%ps)
140 safe_deallocate_a(spec%ps)
141
142 call species_end(spec)
145 end subroutine pseudopotential_finalize
146
147 ! ---------------------------------------------------------
148 logical pure function pseudopotential_has_nlcc(spec)
149 class(pseudopotential_t), intent(in) :: spec
150 pseudopotential_has_nlcc = spec%nlcc
152
154 ! ---------------------------------------------------------
155 integer pure function pseudopotential_x_functional(spec)
156 class(pseudopotential_t), intent(in) :: spec
157
158 pseudopotential_x_functional = spec%ps%exchange_functional
160 ! if we do not know, try the pseudpotential set
161 if (pseudopotential_x_functional == pseudo_exchange_unknown) then
162 select case (spec%pseudopotential_set_id)
163 case ( &
164 option__pseudopotentialset__standard, &
165 option__pseudopotentialset__hgh_lda, &
166 option__pseudopotentialset__hgh_lda_sc, &
167 option__pseudopotentialset__hscv_lda)
169 pseudopotential_x_functional = option__xcfunctional__lda_x
171 case (option__pseudopotentialset__hscv_pbe)
172 pseudopotential_x_functional = option__xcfunctional__gga_x_pbe
173 end select
174 end if
178 ! ---------------------------------------------------------
180 integer pure function pseudopotential_c_functional(spec)
181 class(pseudopotential_t), intent(in) :: spec
182
183 pseudopotential_c_functional = spec%ps%correlation_functional
184
185 ! if we do not know, try the pseudpotential set
186 if (pseudopotential_c_functional == pseudo_correlation_unknown) then
187 select case (spec%pseudopotential_set_id)
188 case ( &
189 option__pseudopotentialset__standard, &
190 option__pseudopotentialset__hgh_lda, &
191 option__pseudopotentialset__hgh_lda_sc, &
192 option__pseudopotentialset__hscv_lda)
193
194 pseudopotential_c_functional = option__xcfunctional__lda_c_pz_mod/1000
196 case (option__pseudopotentialset__hscv_pbe)
197 pseudopotential_c_functional = option__xcfunctional__gga_c_pbe/1000
198 end select
199 end if
200
202
203 ! ---------------------------------------------------------
205 real(real64) pure function pseudopotential_get_radius(spec) result(radius)
206 class(pseudopotential_t), intent(in) :: spec
207
208 radius = spec%ps%rc_max
209 end function pseudopotential_get_radius
210
211 ! ---------------------------------------------------------
212 integer pure function pseudopotential_get_user_lloc(spec) result(lloc)
213 class(pseudopotential_t), intent(in) :: spec
214 lloc = spec%user_llocal
216
217 ! ---------------------------------------------------------
218 integer pure function pseudopotential_get_user_lmax(spec) result(lmax)
219 class(pseudopotential_t), intent(in) :: spec
220 lmax = spec%user_lmax
223 ! ---------------------------------------------------------
224 pure subroutine pseudopotential_set_user_lmax(spec, ll)
225 class(pseudopotential_t), intent(inout) :: spec
226 integer, intent(in) :: ll
227 spec%user_lmax = ll
228 end subroutine pseudopotential_set_user_lmax
229
230 ! ---------------------------------------------------------
231 pure subroutine pseudopotential_set_user_lloc(spec, ll)
232 class(pseudopotential_t), intent(inout) :: spec
233 integer, intent(in) :: ll
234 spec%user_llocal = ll
235 end subroutine pseudopotential_set_user_lloc
236
237 ! ---------------------------------------------------------
238
239 character(len=MAX_PATH_LEN) function get_set_directory(set_id) result(filename)
240 integer, intent(in) :: set_id
242 push_sub(get_set_directory)
243
244 select case (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)
266 filename = ''
267 case default
268 assert(.false.)
269 end select
270
271 pop_sub(get_set_directory)
272 end function get_set_directory
274 ! ---------------------------------------------------------
276 subroutine read_from_set(spec, set_id, set, read_data)
277 type(pseudopotential_t), intent(inout) :: spec
278 integer, intent(in) :: set_id
279 type(pseudo_set_t), intent(in) :: set
280 integer, intent(out) :: read_data
281
282 type(element_t) :: el
283
284 push_sub(read_from_set)
285
286 spec%pseudopotential_set_id = set_id
287 spec%pseudopotential_set = set
288
289 call element_init(el, get_symbol(spec%get_label()))
290
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))
293
294 ! these might have been set before
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))
300 read_data = 8
301 else
302 read_data = 0
303 end if
304
305 call element_end(el)
306
307 pop_sub(read_from_set)
308 end subroutine read_from_set
309
310
311 ! ---------------------------------------------------------
312 subroutine read_from_default_file(iunit, read_data, spec)
313 integer, intent(in) :: iunit
314 integer, intent(inout) :: read_data
315 class(pseudopotential_t), intent(inout) :: spec
316
317 character(len=LABEL_LEN) :: label
318 character(len=MAX_PATH_LEN) :: filename
319 type(element_t) :: element
320 integer :: lmax, llocal
321 real(real64) :: zz
322
323 push_sub(read_from_default_file)
325 backspace(iunit)
326
327 read(iunit,*) label, filename, zz, lmax, llocal
328
329 call spec%set_filename(trim(conf%share)//'/pseudopotentials/'//trim(filename))
330
331 assert(trim(label) == trim(spec%get_label()))
333 read_data = 8
334
335 ! get the mass, vdw radius and atomic number for this element
336 call element_init(element, get_symbol(label))
337
338 assert(element_valid(element))
339
340 ! these might have been set before
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))
347
348 call element_end(element)
349
351 end subroutine read_from_default_file
352
353 ! ---------------------------------------------------------
356 subroutine pseudopotential_real_nl_projector(spec, np, x, r, l, lm, i, uV)
357 class(pseudopotential_t), intent(in) :: spec
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(:)
363
364 integer :: ip
365 real(real64) :: ylm
366
368
369 if (np > 0) then
370 uv(1:np) = r(1:np)
371 call spline_eval_vec(spec%ps%kb(l, i), np, uv)
372
373 do ip = 1, np
374 call ylmr_real(x(1:3, ip), l, lm, ylm)
375 uv(ip) = uv(ip) * ylm
376 end do
377 end if
378
381
382 ! ---------------------------------------------------------
385 subroutine pseudopotential_nl_projector(spec, np, x, r, l, lm, i, uV)
386 class(pseudopotential_t), intent(in) :: spec
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(:)
392
393 integer :: ip
394 complex(real64) :: ylm
395
397
398 if (np > 0) then
399 uv(1:np) = r(1:np)
400 call spline_eval_vec(spec%ps%kb(l, i), np, uv)
401
402 do ip = 1, np
403 call ylmr_cmplx(x(1:3, ip), l, lm, ylm)
404 uv(ip) = uv(ip) * ylm
405 end do
406 end if
407
409 end subroutine pseudopotential_nl_projector
410
411 ! ---------------------------------------------------------
413 subroutine pseudopotential_iwf_fix_qn(spec, namespace, nspin, dim)
414 class(pseudopotential_t), intent(inout) :: spec
415 type(namespace_t), intent(in) :: namespace
416 integer, intent(in) :: nspin
417 integer, intent(in) :: dim
418
419 integer :: is, n, i, l, m
420
422
423 do is = 1, nspin
424 n = 1
425 do i = 1, spec%ps%conf%p
426 if (n > spec%niwfs) exit
427 l = spec%ps%conf%l(i)
428
429 if (.not. spec%ps%bound(i,is)) cycle
430
431 do m = -l, l
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)
437 n = n + 1
438 end do
439
440 end do
441 end do
442
444 end subroutine pseudopotential_iwf_fix_qn
445
446 ! ---------------------------------------------------------
448 real(real64) function pseudopotential_get_iwf_radius(spec, ii, is, threshold) result(radius)
449 class(pseudopotential_t), intent(in) :: spec
450 integer, intent(in) :: ii
451 integer, intent(in) :: is
452 real(real64), optional, intent(in) :: threshold
453
454 real(real64) threshold_
455
457
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_)
461
464
465 ! ---------------------------------------------------------
466 logical function pseudopotential_is_local(spec) result(is_local)
467 class(pseudopotential_t), intent(in) :: spec
468
470
471 is_local = .false.
472 if (spec%ps%lmax == 0 .and. spec%ps%llocal == 0) is_local = .true.
473
475 end function pseudopotential_is_local
476
477 ! ---------------------------------------------------------
481 ! ---------------------------------------------------------
482 subroutine pseudopotential_init_potential(this, namespace, grid_cutoff, filter)
483 class(pseudopotential_t), intent(inout) :: this
484 type(namespace_t), intent(in) :: namespace
485 real(real64), intent(in) :: grid_cutoff
486 integer, intent(in) :: filter
487
488 character(len=256) :: dirname
489 integer :: iorb
490 real(real64) :: local_radius, orbital_radius
491
493
494 call ps_separate(this%ps)
495
496 call ps_getradius(this%ps)
497
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) ! radius may have changed
501 end if
502
503 call ps_derivatives(this%ps)
504
505 local_radius = this%ps%vl%x_threshold
507 orbital_radius = m_zero
508 ! FIXME: should take max over spins too here.
509 do iorb = 1, this%get_niwfs()
510 orbital_radius = max(orbital_radius, this%get_iwf_radius(this%iwf_i(iorb, 1), is = 1))
511 end do
512
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)
522
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)
527 end if
528
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)
533 end if
534
535 if (debug%info) then
536 write(dirname, '(a)') 'debug/geometry'
537 call io_mkdir(dirname, namespace)
538 call this%debug(trim(dirname), namespace, grid_cutoff)
539 end if
540
542 end subroutine pseudopotential_init_potential
543
544 ! ---------------------------------------------------------
545 subroutine pseudopotential_debug(spec, dir, namespace, gmax)
546 class(pseudopotential_t), intent(in) :: spec
547 character(len=*), intent(in) :: dir
548 type(namespace_t), intent(in) :: namespace
549 real(real64), intent(in) :: gmax
550
551
552 character(len=256) :: dirname
553 integer :: iunit
554
555 if (.not. mpi_grp_is_root(mpi_world)) then
556 return
557 end if
558
560
561 dirname = trim(dir)//'/'//trim(spec%get_label())
562
563 call io_mkdir(dirname, namespace)
564
565 iunit = io_open(trim(dirname)//'/info', namespace, action='write')
566
567 write(iunit, '(a,i3)') 'Index = ', spec%get_index()
568 write(iunit, '(2a)') 'Label = ', trim(spec%get_label())
569 !write(iunit, '(a,i3)') 'Type = ', spec%type
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()
580
581 if (debug%info) call ps_debug(spec%ps, trim(dirname), namespace, gmax)
582
583 call io_close(iunit)
584 pop_sub(pseudopotential_debug)
585 end subroutine pseudopotential_debug
586
587 ! ---------------------------------------------------------
588 subroutine pseudopotential_build(spec, namespace, ispin, dim, print_info)
589 class(pseudopotential_t), intent(inout) :: spec
590 type(namespace_t), intent(in) :: namespace
591 integer, intent(in) :: ispin
592 integer, intent(in) :: dim
593 logical, optional, intent(in) :: print_info
594
595 logical :: print_info_
596
597 push_sub(pseudopotential_build)
598
599 print_info_ = optional_default(print_info, .true.)
600
601 ! masses are always in amu, so convert them to a.u.
602 call spec%set_mass(units_to_atomic(unit_amu, spec%get_mass()))
603
604 spec%has_density = .false.
605
606 ! allocate structure
607 safe_allocate(spec%ps)
608 if (spec%type == species_pspio) then
609 call ps_pspio_init(spec%ps, namespace, spec%get_label(), spec%get_z(), spec%user_lmax, &
610 ispin, spec%get_filename())
611 else
612 call ps_init(spec%ps, namespace, spec%get_label(), spec%get_z(), spec%user_lmax, &
613 spec%user_llocal, ispin, spec%get_filename())
614 end if
615 call spec%set_zval(spec%ps%z_val)
616 spec%nlcc = spec%ps%nlcc
617 spec%niwfs = ps_bound_niwfs(spec%ps)
618
619 ! invalidate these variables as they should not be used after
620 spec%user_lmax = invalid_l
621 spec%user_llocal = invalid_l
622
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))
628
629 call spec%iwf_fix_qn(namespace, ispin, dim)
630
631 pop_sub(pseudopotential_build)
632 end subroutine pseudopotential_build
633
634 ! ---------------------------------------------------------
636 logical pure function pseudopotential_is_ps(this)
637 class(pseudopotential_t), intent(in) :: this
640 end function pseudopotential_is_ps
641
642 ! ---------------------------------------------------------
644 logical pure function pseudopotential_is_ps_with_nlcc(this)
645 class(pseudopotential_t), intent(in) :: this
646
647 pseudopotential_is_ps_with_nlcc = this%has_nlcc()
649
650 ! ---------------------------------------------------------
652 logical pure function pseudopotential_represents_real_atom(spec)
653 class(pseudopotential_t), intent(in) :: spec
654
657
658
659
660end module pseudopotential_oct_m
661
662!! Local Variables:
663!! mode: f90
664!! coding: utf-8
665!! 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)