Octopus
jellium.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
21module jellium_oct_m
22 use debug_oct_m
23 use global_oct_m
24 use io_oct_m
26 use mpi_oct_m
28 use parser_oct_m
31 use unit_oct_m
33
34 implicit none
35
36 private
37 public :: &
38 jellium_t, &
45
46 integer, public, parameter :: &
47 SPECIES_JELLIUM = 3, & !< jellium sphere.
50 species_usdef = 123, &
53
54 type, abstract, extends(species_t) :: jellium_t
55 private
56
57 real(real64) :: omega
58
59 contains
60 procedure :: iwf_fix_qn => jellium_iwf_fix_qn
61 procedure :: get_iwf_radius => jellium_get_iwf_radius
62 procedure :: is_local => jellium_is_local
63 procedure :: init_potential => jellium_init_potential
64 procedure :: debug => jellium_debug
65 procedure :: build => jellium_build
66 procedure :: get_omega => jellium_get_omega
67 procedure :: is_user_defined => jellium_user_defined
68 end type jellium_t
69
70 type, extends(jellium_t) :: jellium_sphere_t
71 private
72
73 real(real64) :: jradius
74
75 contains
76 procedure :: radius => jellium_radius
77 procedure :: set_radius => jellium_set_radius
79 end type jellium_sphere_t
80
81 type, extends(jellium_t) :: jellium_slab_t
82 private
83
84 real(real64) :: jthick
85
86 contains
87 procedure :: thickness => jellium_thick
88 procedure :: set_thickness => jellium_set_thickness
90 end type jellium_slab_t
91
92 type, extends(jellium_t) :: jellium_charge_t
93 private
94
95 character(len=200), public :: density_formula
96
97 contains
98 procedure :: rho_string => jellium_rho_string
100 end type jellium_charge_t
101
102 type, extends (jellium_t) :: species_from_file_t
103 private
104
105 contains
107 end type species_from_file_t
108
109 type, extends (jellium_t) :: species_user_defined_t
110 private
111
112 character(len=1024), public :: potential_formula
113 contains
114 procedure :: user_pot => jellium_userdef_pot
117
118 type, extends (jellium_t) :: species_charge_density_t
119 private
120
121 character(len=200), public :: density_formula
122 contains
123 procedure :: rho_string => species_rho_string
126
127 interface jellium_sphere_t
128 procedure jellium_sphere_constructor
129 end interface jellium_sphere_t
130
131 interface jellium_slab_t
132 procedure jellium_slab_constructor
133 end interface jellium_slab_t
134
135 interface jellium_charge_t
136 procedure jellium_charge_constructor
137 end interface jellium_charge_t
138
141 end interface species_from_file_t
142
143 interface species_user_defined_t
145 end interface species_user_defined_t
146
148 procedure species_charge_density_constructor
149 end interface species_charge_density_t
151
152contains
155 ! ---------------------------------------------------------
156 function jellium_slab_constructor(label, index) result(spec)
157 class(jellium_slab_t), pointer :: spec
158 character(len=*), intent(in) :: label
159 integer, intent(in) :: index
162
163 safe_allocate(spec)
164
165 call species_init(spec, label, index)
167 spec%omega = m_zero
168 spec%jthick = -m_one
172
173 ! ---------------------------------------------------------
174 subroutine jellium_slab_finalize(this)
175 type(jellium_slab_t), intent(inout) :: this
176
178
179 call species_end(this)
183
184 ! ---------------------------------------------------------
185 function jellium_sphere_constructor(label, index) result(spec)
186 class(jellium_sphere_t), pointer :: spec
187 character(len=*), intent(in) :: label
188 integer, intent(in) :: index
189
192 safe_allocate(spec)
193
194 call species_init(spec, label, index)
196 spec%omega = m_zero
197 spec%jradius = m_half
198
200 end function jellium_sphere_constructor
201
203 ! ---------------------------------------------------------
204 subroutine jellium_sphere_finalize(this)
205 type(jellium_sphere_t), intent(inout) :: this
206
209 call species_end(this)
210
212 end subroutine jellium_sphere_finalize
213
214 ! ---------------------------------------------------------
215 function jellium_charge_constructor(label, index) result(spec)
216 class(jellium_charge_t), pointer :: spec
217 character(len=*), intent(in) :: label
218 integer, intent(in) :: index
219
221
222 safe_allocate(spec)
223
224 call species_init(spec, label, index)
225
226 spec%omega = m_zero
227 spec%density_formula = ""
228
230 end function jellium_charge_constructor
231
232
233 ! ---------------------------------------------------------
234 subroutine jellium_charge_finalize(this)
235 type(jellium_charge_t), intent(inout) :: this
236
238
239 call species_end(this)
240
242 end subroutine jellium_charge_finalize
243
244
245 ! ---------------------------------------------------------
246 function species_from_file_constructor(label, index) result(spec)
247 class(species_from_file_t), pointer :: spec
248 character(len=*), intent(in) :: label
249 integer, intent(in) :: index
250
252
253 safe_allocate(spec)
254
255 call species_init(spec, label, index)
256
257 spec%omega = m_zero
258
261
262 ! ---------------------------------------------------------
263 subroutine species_from_file_finalize(this)
264 type(species_from_file_t), intent(inout) :: this
265
268 call species_end(this)
269
271 end subroutine species_from_file_finalize
272
273
274 ! ---------------------------------------------------------
275 function species_user_defined_constructor(label, index) result(spec)
276 class(species_user_defined_t), pointer :: spec
277 character(len=*), intent(in) :: label
278 integer, intent(in) :: index
279
281
282 safe_allocate(spec)
283
284 call species_init(spec, label, index)
285
286 spec%potential_formula = ""
287 spec%omega = m_zero
288
291
292
293 ! ---------------------------------------------------------
294 subroutine species_user_defined_finalize(this)
295 type(species_user_defined_t), intent(inout) :: this
296
298
299 call species_end(this)
300
302 end subroutine species_user_defined_finalize
303
304
305 ! ---------------------------------------------------------
306 function species_charge_density_constructor(label, index) result(spec)
307 class(species_charge_density_t), pointer :: spec
308 character(len=*), intent(in) :: label
309 integer, intent(in) :: index
310
312
313 safe_allocate(spec)
314
315 call species_init(spec, label, index)
316
317 spec%omega = m_zero
318 spec%density_formula = ""
319
322
323
324 ! ---------------------------------------------------------
325 subroutine species_charge_density_finalize(this)
326 type(species_charge_density_t), intent(inout) :: this
329
330 call species_end(this)
331
334
335 ! ---------------------------------------------------------
336 real(real64) pure function jellium_get_omega(spec)
337 class(jellium_t), intent(in) :: spec
338 jellium_get_omega = spec%omega
339 end function jellium_get_omega
340
341
342 ! ---------------------------------------------------------
343 real(real64) pure function jellium_radius(spec)
344 class(jellium_sphere_t), intent(in) :: spec
345 jellium_radius = spec%jradius
346 end function jellium_radius
347
348 ! ---------------------------------------------------------
349 pure subroutine jellium_set_radius(spec, radius)
350 class(jellium_sphere_t), intent(inout) :: spec
351 real(real64), intent(in) :: radius
352 spec%jradius = radius
353 end subroutine jellium_set_radius
354
355 ! ---------------------------------------------------------
356 real(real64) pure function jellium_thick(spec)
357 class(jellium_slab_t), intent(in) :: spec
358 jellium_thick = spec%jthick
359 end function jellium_thick
360
361 ! ---------------------------------------------------------
362 pure subroutine jellium_set_thickness(spec, thick)
363 class(jellium_slab_t), intent(inout) :: spec
364 real(real64), intent(in) :: thick
365 spec%jthick = thick
366 end subroutine jellium_set_thickness
367
368 ! ---------------------------------------------------------
369 character(len=200) pure function jellium_rho_string(spec)
370 class(jellium_charge_t), intent(in) :: spec
371 jellium_rho_string = trim(spec%density_formula)
372 end function jellium_rho_string
373
374 ! ---------------------------------------------------------
375 character(len=200) pure function species_rho_string(spec)
376 class(species_charge_density_t), intent(in) :: spec
377 species_rho_string = trim(spec%density_formula)
378 end function species_rho_string
379
380
381 ! ---------------------------------------------------------
382 complex(real64) function jellium_userdef_pot(spec, dim, xx, r)
383 class(species_user_defined_t), intent(in) :: spec
384 integer, intent(in) :: dim
385 real(real64), intent(in) :: xx(:)
386 real(real64), intent(in) :: r
388 real(real64) :: pot_re, pot_im
389
390 push_sub(jellium_userdef_pot)
391
392 call parse_expression(pot_re, pot_im, dim, xx, r, m_zero, spec%potential_formula)
393 jellium_userdef_pot = pot_re + m_zi * pot_im
394
395 pop_sub(jellium_userdef_pot)
396 end function jellium_userdef_pot
397
398 ! ---------------------------------------------------------
400 subroutine jellium_iwf_fix_qn(spec, namespace, nspin, dim)
401 class(jellium_t), intent(inout) :: spec
402 type(namespace_t), intent(in) :: namespace
403 integer, intent(in) :: nspin
404 integer, intent(in) :: dim
405
406 integer :: is, i, n1, n2, n3
407
408 push_sub(jellium_iwf_fix_qn)
409
410 select case (dim)
411 case (1)
412 do is = 1, nspin
413 do i = 1, spec%niwfs
414 spec%iwf_i(i, is) = i
415 spec%iwf_n(i, is) = 0
416 spec%iwf_l(i, is) = 0
417 spec%iwf_m(i, is) = 0
418 spec%iwf_j(i) = m_zero
419 end do
420 end do
421
422 case (2)
423 do is = 1, nspin
424 i = 1
425 n1 = 1
426 n2 = 1
427 do
428 spec%iwf_i(i, is) = n1
429 spec%iwf_n(i, is) = 1
430 spec%iwf_l(i, is) = n2
431 spec%iwf_m(i, is) = 0
432 spec%iwf_j(i) = m_zero
433 i = i + 1
434 if (i>spec%niwfs) exit
435
436 spec%iwf_i(i, is) = n1+1
437 spec%iwf_n(i, is) = 1
438 spec%iwf_l(i, is) = n2
439 spec%iwf_m(i, is) = 0
440 spec%iwf_j(i) = m_zero
441 i = i + 1
442 if (i>spec%niwfs) exit
443
444 spec%iwf_i(i, is) = n1
445 spec%iwf_n(i, is) = 1
446 spec%iwf_l(i, is) = n2+1
447 spec%iwf_m(i, is) = 0
448 spec%iwf_j(i) = m_zero
449 i = i + 1
450 if (i>spec%niwfs) exit
451
452 n1 = n1 + 1; n2 = n2 + 1
453 end do
454 end do
456 case (3)
457 do is = 1, nspin
458 i = 1
459 n1 = 1
460 n2 = 1
461 n3 = 1
462 do
463 spec%iwf_i(i, is) = n1
464 spec%iwf_n(i, is) = 1
465 spec%iwf_l(i, is) = n2
466 spec%iwf_m(i, is) = n3
467 spec%iwf_j(i) = m_zero
468 i = i + 1
469 if (i>spec%niwfs) exit
470
471 spec%iwf_i(i, is) = n1+1
472 spec%iwf_n(i, is) = 1
473 spec%iwf_l(i, is) = n2
474 spec%iwf_m(i, is) = n3
475 spec%iwf_j(i) = m_zero
476 i = i + 1
477 if (i>spec%niwfs) exit
478
479 spec%iwf_i(i, is) = n1
480 spec%iwf_n(i, is) = 1
481 spec%iwf_l(i, is) = n2+1
482 spec%iwf_m(i, is) = 0
483 spec%iwf_j(i) = m_zero
484 i = i + 1
485 if (i>spec%niwfs) exit
486
487 spec%iwf_i(i, is) = n1
488 spec%iwf_n(i, is) = 1
489 spec%iwf_l(i, is) = n2
490 spec%iwf_m(i, is) = n3+1
491 spec%iwf_j(i) = m_zero
492 i = i + 1
493 if (i>spec%niwfs) exit
494
495 spec%iwf_i(i, is) = n1+1
496 spec%iwf_n(i, is) = 1
497 spec%iwf_l(i, is) = n2+1
498 spec%iwf_m(i, is) = n3
499 spec%iwf_j(i) = m_zero
500 i = i + 1
501 if (i>spec%niwfs) exit
502
503 spec%iwf_i(i, is) = n1+1
504 spec%iwf_n(i, is) = 1
505 spec%iwf_l(i, is) = n2
506 spec%iwf_m(i, is) = n3+1
507 spec%iwf_j(i) = m_zero
508 i = i + 1
509 if (i>spec%niwfs) exit
510
511 spec%iwf_i(i, is) = n1
512 spec%iwf_n(i, is) = 1
513 spec%iwf_l(i, is) = n2+1
514 spec%iwf_m(i, is) = n3+1
515 spec%iwf_j(i) = m_zero
516 i = i + 1
517 if (i>spec%niwfs) exit
518
519 n1 = n1 + 1
520 n2 = n2 + 1
521 n3 = n3 + 1
522 end do
523 end do
524 case default
525 ! Not doing anything to allow for N-D simulations
526 end select
527
528 pop_sub(jellium_iwf_fix_qn)
529 end subroutine jellium_iwf_fix_qn
530
531 ! ---------------------------------------------------------
533 real(real64) pure function jellium_get_iwf_radius(spec, ii, is, threshold) result(radius)
534 class(jellium_t), intent(in) :: spec
535 integer, intent(in) :: ii
536 integer, intent(in) :: is
537 real(real64), optional, intent(in) :: threshold
538
539 real(real64) threshold_
540
541 threshold_ = optional_default(threshold, 0.001_real64)
542
543 radius = sqrt(-m_two*log(threshold_)/spec%omega)
544
545 ! The values for hydrogenic and harmonic-oscillator wavefunctions
546 ! come from taking the exponential part (i.e. the one that controls
547 ! the asymptotic behavior at large r), and setting it equal to
548 ! the threshold.
549 end function jellium_get_iwf_radius
550
551 ! ---------------------------------------------------------
552 logical pure function jellium_is_local(spec) result(is_local)
553 class(jellium_t), intent(in) :: spec
554
555 is_local = .true.
556 end function jellium_is_local
557
558 ! ---------------------------------------------------------
560 ! ---------------------------------------------------------
561 subroutine jellium_init_potential(this, namespace, grid_cutoff, filter)
562 class(jellium_t), intent(inout) :: this
563 type(namespace_t), intent(in) :: namespace
564 real(real64), intent(in) :: grid_cutoff
565 integer, intent(in) :: filter
566
567 push_sub(jellium_init_potential)
568
569
571 end subroutine jellium_init_potential
572
573 ! ---------------------------------------------------------
574 subroutine jellium_debug(spec, dir, namespace, gmax)
575 class(jellium_t), intent(in) :: spec
576 character(len=*), intent(in) :: dir
577 type(namespace_t), intent(in) :: namespace
578 real(real64), intent(in) :: gmax
579
580 character(len=256) :: dirname
581 integer :: iunit
582
583 if (.not. mpi_grp_is_root(mpi_world)) then
584 return
585 end if
586
587 push_sub(jellium_debug)
588
589 dirname = trim(dir)//'/'//trim(spec%get_label())
590
591 call io_mkdir(dirname, namespace)
592
593 iunit = io_open(trim(dirname)//'/info', namespace, action='write')
594
595 write(iunit, '(a,i3)') 'Index = ', spec%get_index()
596 write(iunit, '(2a)') 'Label = ', trim(spec%get_label())
597 write(iunit, '(a,f15.2)') 'z_val = ', spec%get_zval()
598 write(iunit, '(a,f15.2)') 'mass = ', spec%get_mass()
599 write(iunit, '(a,f15.2)') 'vdw_radius = ', spec%get_vdw_radius()
600 write(iunit, '(a,l1)') 'local = ', spec%is_local()
601
602 select type(spec)
603 type is(species_from_file_t)
604 write(iunit, '(a,f15.2)') 'z = ', spec%get_z()
605 write(iunit,'(a)') 'Species read from file "'//trim(spec%get_filename())//'".'
606 type is(jellium_sphere_t)
607 write(iunit, '(a,f15.2)') 'z = ', spec%get_z()
608 write(iunit, '(a,f15.2)') 'jradius= ', spec%radius()
609 type is(jellium_slab_t)
610 write(iunit, '(a,f15.2)') 'z = ', spec%get_z()
611 write(iunit, '(a,f15.2)') 'jthick= ', spec%thickness()
613 write(iunit, '(2a)') 'usdef = ', trim(spec%potential_formula)
614 end select
615
616 write(iunit, '(a,i3)') 'hubbard_l = ', spec%get_hubbard_l()
617 write(iunit, '(a,f15.2)') 'hubbard_U = ', spec%get_hubbard_U()
618 write(iunit, '(a,f15.2)') 'hubbard_j = ', spec%get_hubbard_j()
619 write(iunit, '(a,f15.2)') 'hubbard_alpha = ', spec%get_hubbard_alpha()
620
621 call io_close(iunit)
622 pop_sub(jellium_debug)
623 end subroutine jellium_debug
624
625 ! ---------------------------------------------------------
626 subroutine jellium_build(spec, namespace, ispin, dim, print_info)
627 class(jellium_t), intent(inout) :: spec
628 type(namespace_t), intent(in) :: namespace
629 integer, intent(in) :: ispin
630 integer, intent(in) :: dim
631 logical, optional, intent(in) :: print_info
632
633 logical :: print_info_
634 integer :: i
635 real(real64) :: pot_re, pot_im, xx(dim), rr
636
637 push_sub(jellium_build)
638
639 print_info_ = optional_default(print_info, .true.)
640
641 ! masses are always in amu, so convert them to a.u.
642 call spec%set_mass(units_to_atomic(unit_amu, spec%get_mass()))
643
644 spec%has_density = .false.
646 select type (spec)
648 if (print_info_) then
649 write(message(1),'(a,a,a)') 'Species "',trim(spec%get_label()),'" is a user-defined potential.'
650 i = min(237, len_trim(spec%potential_formula)-1) ! I subtract 1 to avoid the non-printable C "end-of-string" character.
651 write(message(2),'(a,a)') ' Potential = ', trim(spec%potential_formula(1:i))
652 if (len(trim(spec%potential_formula)) > 237) then
653 message(2) = trim(message(2))//'...'
654 end if
655 call messages_info(2, namespace=namespace)
656 end if
657 spec%niwfs = int(max(2*spec%get_zval(), m_one))
658
659 xx = m_zero
660 xx(1) = 0.01_real64
661 rr = norm2(xx)
662 call parse_expression(pot_re, pot_im, dim, xx, rr, m_zero, spec%potential_formula)
663 spec%omega = sqrt(abs(m_two / 1.0e-4_real64 * pot_re)) ! why...?
664 ! To avoid problems with constant potentials.
665 if (spec%omega <= m_zero) spec%omega = 0.1_real64
666
668 if (print_info_) then
669 write(message(1),'(a)') 'Species read from file "'//trim(spec%get_filename())//'".'
670 call messages_info(1, namespace=namespace)
671 end if
672 spec%niwfs = 2*nint(spec%get_zval()+m_half)
673 spec%omega = 0.1_real64
674
675 type is(jellium_sphere_t)
676 if (print_info_) then
677 write(message(1),'(a,a,a)') 'Species "', trim(spec%get_label()), &
678 '" is a jellium sphere / approximated point particle.'
679 write(message(2),'(a,f11.6)') ' Valence charge = ', spec%get_zval()
680 write(message(3),'(a,f11.6)') ' Radius [a.u] = ', spec%jradius
681 write(message(4),'(a,f11.6)') ' Rs [a.u] = ', spec%jradius * spec%get_zval() ** (-m_one/m_three)
682 call messages_info(4, namespace=namespace)
683 end if
684 spec%niwfs = species_closed_shell_size(2*nint(spec%get_zval()+m_half))
685 spec%omega = 0.1_real64
686
687 type is(jellium_slab_t)
688 if (print_info_) then
689 write(message(1),'(a,a,a)') 'Species "',trim(spec%get_label()),'" is a jellium slab.'
690 write(message(2),'(a,f11.6)') ' Valence charge = ', spec%get_zval()
691 write(message(3),'(a,f11.6)') ' Thickness [a.u] = ', spec%jthick
692 !write(message(4),'(a,f11.6)') ' Rs [a.u] = ', ( M_THREE /( M_FOUR *M_PI ) &
693 !& *spec%get_zval() /( *sb%lsize(1) *sb%lsize(2) ) )**(1.0/3.0)
694 call messages_info(3, namespace=namespace)
695 end if
696 spec%niwfs = 2*nint(spec%get_zval()+m_half)
697 spec%omega = 0.1_real64
698
699 type is(jellium_charge_t)
700 spec%niwfs = int(max(2*spec%get_zval(), m_one))
701 spec%omega = spec%get_zval()
702 spec%has_density = .true.
703 if (print_info_) then
704 write(message(1),'(a,a,a)') 'Species "', trim(spec%get_label()), '" is a distribution of charge:'
705 write(message(2),'(a,a,a)') ' rho is enclosed in volume defined by the "', &
706 trim(spec%density_formula), '" block'
707 write(message(3),'(a,f11.6)') ' Z = ', spec%get_zval()
708 call messages_info(3, namespace=namespace)
709 end if
710
712 spec%niwfs = int(max(2*spec%get_zval(), m_one))
713 spec%omega = spec%get_zval()
714 spec%has_density = .true.
715 if (print_info_) then
716 write(message(1),'(a,a,a)') 'Species "', trim(spec%get_label()), '" is a distribution of charge:'
717 write(message(2),'(a,a)') ' rho = ', trim(spec%density_formula)
718 write(message(3),'(a,f11.6)') ' Z = ', spec%get_zval()
719 call messages_info(3, namespace=namespace)
720 end if
721 class default
722 call messages_input_error(namespace, 'Species', 'Unknown species type')
723 end select
724
725 ! since there is no real cap, make sure there are at least a few available
726 spec%niwfs = max(5, spec%niwfs)
727
728 safe_allocate(spec%iwf_n(1:spec%niwfs, 1:ispin))
729 safe_allocate(spec%iwf_l(1:spec%niwfs, 1:ispin))
730 safe_allocate(spec%iwf_m(1:spec%niwfs, 1:ispin))
731 safe_allocate(spec%iwf_i(1:spec%niwfs, 1:ispin))
732 safe_allocate(spec%iwf_j(1:spec%niwfs))
733
734 call spec%iwf_fix_qn(namespace, ispin, dim)
735
736 write(message(1),'(a,i6,a,i6)') 'Number of orbitals: ', spec%niwfs
737 if (print_info_) call messages_info(1, namespace=namespace)
738
739 pop_sub(jellium_build)
740 end subroutine jellium_build
741
742 ! ---------------------------------------------------------
744 logical pure function jellium_user_defined(spec)
745 class(jellium_t), intent(in) :: spec
746
747 select type(spec)
748 class is(species_user_defined_t)
752 class is(species_from_file_t)
754 class default
755 jellium_user_defined = .false.
756 end select
757
758 end function jellium_user_defined
759
760
761end module jellium_oct_m
762
763!! Local Variables:
764!! mode: f90
765!! coding: utf-8
766!! End:
double log(double __x) __attribute__((__nothrow__
double sqrt(double __x) __attribute__((__nothrow__
real(real64), parameter, public m_zero
Definition: global.F90:187
real(real64), parameter, public m_half
Definition: global.F90:193
real(real64), parameter, public m_one
Definition: global.F90:188
Definition: io.F90:114
subroutine species_user_defined_finalize(this)
Definition: jellium.F90:388
real(real64) pure function jellium_get_omega(spec)
Definition: jellium.F90:430
real(real64) pure function jellium_radius(spec)
Definition: jellium.F90:437
logical pure function jellium_user_defined(spec)
Is the species user-defined or not.
Definition: jellium.F90:838
logical pure function jellium_is_local(spec)
Definition: jellium.F90:646
pure subroutine jellium_set_thickness(spec, thick)
Definition: jellium.F90:456
subroutine jellium_iwf_fix_qn(spec, namespace, nspin, dim)
set up quantum numbers of orbitals
Definition: jellium.F90:494
integer, parameter, public species_charge_density
user-defined function for charge density
Definition: jellium.F90:139
integer, parameter, public species_jellium_charge_density
jellium volume read from file
Definition: jellium.F90:139
subroutine jellium_charge_finalize(this)
Definition: jellium.F90:328
class(jellium_slab_t) function, pointer jellium_slab_constructor(label, index)
Definition: jellium.F90:250
class(species_user_defined_t) function, pointer species_user_defined_constructor(label, index)
Definition: jellium.F90:369
subroutine jellium_slab_finalize(this)
Definition: jellium.F90:268
subroutine jellium_build(spec, namespace, ispin, dim, print_info)
Definition: jellium.F90:720
real(real64) pure function jellium_get_iwf_radius(spec, ii, is, threshold)
Return radius outside which orbital is less than threshold value 0.001.
Definition: jellium.F90:627
class(jellium_charge_t) function, pointer jellium_charge_constructor(label, index)
Definition: jellium.F90:309
class(jellium_sphere_t) function, pointer jellium_sphere_constructor(label, index)
Definition: jellium.F90:279
pure subroutine jellium_set_radius(spec, radius)
Definition: jellium.F90:443
subroutine species_charge_density_finalize(this)
Definition: jellium.F90:419
character(len=200) pure function jellium_rho_string(spec)
Definition: jellium.F90:463
subroutine species_from_file_finalize(this)
Definition: jellium.F90:357
real(real64) pure function jellium_thick(spec)
Definition: jellium.F90:450
complex(real64) function jellium_userdef_pot(spec, dim, xx, r)
Definition: jellium.F90:476
integer, parameter, public species_from_file
Definition: jellium.F90:139
class(species_from_file_t) function, pointer species_from_file_constructor(label, index)
Definition: jellium.F90:340
integer, parameter, public species_usdef
user-defined function for local potential
Definition: jellium.F90:139
subroutine jellium_sphere_finalize(this)
Definition: jellium.F90:298
class(species_charge_density_t) function, pointer species_charge_density_constructor(label, index)
Definition: jellium.F90:400
integer, parameter, public species_jellium_slab
jellium slab.
Definition: jellium.F90:139
subroutine jellium_init_potential(this, namespace, grid_cutoff, filter)
Some operations like filtering of the potentials.
Definition: jellium.F90:655
character(len=200) pure function species_rho_string(spec)
Definition: jellium.F90:469
subroutine jellium_debug(spec, dir, namespace, gmax)
Definition: jellium.F90:668
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)