37 use,
intrinsic :: iso_fortran_env
109 real(real64) :: total(3,3) =
m_zero
110 real(real64) :: kinetic(3,3) =
m_zero
111 real(real64) :: Hartree(3,3) =
m_zero
112 real(real64) :: xc(3,3) =
m_zero
114 real(real64) :: ps_nl(3,3) =
m_zero
115 real(real64) :: ion_ion(3,3) =
m_zero
116 real(real64) :: vdw(3,3) =
m_zero
117 real(real64) :: hubbard(3,3) =
m_zero
119 real(real64) :: kinetic_sumrule
120 real(real64) :: hartree_sumrule
135 type(states_elec_dim_t) :: d
138 logical :: only_userdef_istates
140 type(states_elec_group_t) :: group
141 integer :: block_size
145 logical :: pack_states
151 character(len=1024),
allocatable :: user_def_states(:,:,:)
157 real(real64),
allocatable :: rho(:,:)
158 real(real64),
allocatable :: rho_core(:)
160 real(real64),
allocatable :: current(:, :, :)
161 real(real64),
allocatable :: current_para(:, :, :)
162 real(real64),
allocatable :: current_dia(:, :, :)
163 real(real64),
allocatable :: current_mag(:, :, :)
164 real(real64),
allocatable :: current_kpt(:,:,:)
169 real(real64),
allocatable :: frozen_rho(:, :)
170 real(real64),
allocatable :: frozen_tau(:, :)
171 real(real64),
allocatable :: frozen_gdens(:,:,:)
172 real(real64),
allocatable :: frozen_ldens(:,:)
174 logical :: uniform_occ
176 real(real64),
allocatable :: eigenval(:,:)
178 logical :: restart_fixed_occ
179 logical :: restart_reorder_occs
180 real(real64),
allocatable :: occ(:,:)
181 real(real64),
allocatable :: kweights(:)
184 logical,
private :: fixed_spins
186 real(real64),
allocatable :: spin(:, :, :)
189 real(real64) :: val_charge
191 type(stress_t) :: stress_tensors
193 logical :: fromScratch
194 type(smear_t) :: smear
197 type(modelmb_particle_t) :: modelmbparticles
198 integer,
allocatable :: mmb_nspindown(:,:)
199 integer,
allocatable :: mmb_iyoung(:,:)
200 real(real64),
allocatable :: mmb_proj(:)
202 logical :: parallel_in_states = .false.
213 logical :: scalapack_compatible
217 integer :: st_start, st_end
218 integer,
allocatable :: node(:)
221 integer,
allocatable :: st_kpt_task(:,:)
224 logical :: symmetrize_density
225 integer :: randomization
226 integer :: orth_method = 0
228 real(real64) :: cl_states_mem
240 integer,
public,
parameter :: &
284 real(real64),
intent(in) :: valence_charge
287 real(real64) :: excess_charge, nempty_percent
288 integer :: nempty, ntot, default
289 integer :: nempty_conv
291 real(real64),
parameter :: tol = 1e-13_real64
300 st%d%ispin = space%ispin
305 if (st%d%ispin /=
unpolarized .and. kpoints%use_time_reversal)
then
306 message(1) =
"Time reversal symmetry is only implemented for unpolarized spins."
307 message(2) =
"Use KPointsUseTimeReversal = no."
339 write(
message(1),
'(a,i5,a)')
"Input: '", ntot,
"' is not a valid value for TotalStates."
361 write(
message(1),
'(a,i5,a)')
"Input: '", nempty,
"' is not a valid value for ExtraStates."
362 message(2) =
'(0 <= ExtraStates)'
366 if (ntot > 0 .and. nempty > 0)
then
367 message(1) =
'You cannot set TotalStates and ExtraStates at the same time.'
382 if (nempty_percent < 0)
then
383 write(
message(1),
'(a,f8.6,a)')
"Input: '", nempty_percent, &
384 "' should be a percentage value x (where x is parts in hundred) larger or equal 0"
388 if (nempty > 0 .and. nempty_percent > 0)
then
389 message(1) =
'You cannot set ExtraStates and ExtraStatesInPercent at the same time.'
405 call parse_variable(namespace,
'ExtraStatesToConverge', nempty, nempty_conv)
407 write(
message(1),
'(a,i5,a)')
"Input: '", nempty_conv,
"' is not a valid value for ExtraStatesToConverge."
408 message(2) =
'(0 <= ExtraStatesToConverge)'
412 if (nempty_conv > nempty)
then
413 message(1) =
'You cannot set ExtraStatesToConverge to a higher value than ExtraStates.'
420 st%val_charge = valence_charge
422 st%qtot = -(st%val_charge + excess_charge)
425 write(
message(1),
'(a,f12.6,a)')
'Total charge = ', st%qtot,
' < 0'
426 message(2) =
'Check Species and ExcessCharge.'
430 select case (st%d%ispin)
433 st%nst = nint(st%qtot/2)
434 if (st%nst*2 - st%qtot < -tol) st%nst = st%nst + 1
436 st%d%spin_channels = 1
439 st%nst = nint(st%qtot/2)
440 print *, st%nst, st%nst*2 - st%qtot
441 if (st%nst*2 - st%qtot < -tol) st%nst = st%nst + 1
443 st%d%spin_channels = 2
446 st%nst = nint(st%qtot)
447 if (st%nst - st%qtot < -tol) st%nst = st%nst + 1
449 st%d%spin_channels = 2
453 if (ntot < st%nst)
then
454 message(1) =
'TotalStates is smaller than the number of states required by the system.'
461 if (nempty_percent > 0)
then
462 nempty = ceiling(nempty_percent * st%nst / 100)
465 st%nst_conv = st%nst + nempty_conv
466 st%nst = st%nst + nempty
467 if (st%nst == 0)
then
468 message(1) =
"Cannot run with number of states = zero."
488 default = max(
accel%warp_size, 32)
497 call parse_variable(namespace,
'StatesBlockSize', default, st%block_size)
498 if (st%block_size < 1)
then
499 call messages_write(
"The variable 'StatesBlockSize' must be greater than 0.")
503 st%block_size = min(st%block_size, st%nst)
504 conf%target_states_block_size = st%block_size
506 safe_allocate(st%eigenval(1:st%nst, 1:st%nik))
507 st%eigenval = huge(st%eigenval)
511 if (.not. kpoints%gamma_only())
then
524 call parse_variable(namespace,
'OnlyUserDefinedInitialStates', .false., st%only_userdef_istates)
527 safe_allocate(st%occ (1:st%nst, 1:st%nik))
532 safe_allocate(st%user_def_states(1:st%d%dim, 1:st%nst, 1:st%nik))
534 st%user_def_states(1:st%d%dim, 1:st%nst, 1:st%nik) =
'undefined'
537 if (st%d%ispin ==
spinors)
then
538 safe_allocate(st%spin(1:3, 1:st%nst, 1:st%nik))
565 if (st%smear%photodop)
then
566 if (nempty == 0)
then
567 write(
message(1),
'(a,i5,a)')
"PhotoDoping requires to specify ExtraStates."
568 message(2) =
'(0 == ExtraStates)'
576 safe_allocate(st%node(1:st%nst))
577 st%node(1:st%nst) = 0
580 st%parallel_in_states = .false.
585 if (st%modelmbparticles%nparticle > 0)
then
587 safe_allocate(st%mmb_nspindown(1:st%modelmbparticles%ntype_of_particle, 1:st%nst))
588 st%mmb_nspindown(:,:) = -1
589 safe_allocate(st%mmb_iyoung(1:st%modelmbparticles%ntype_of_particle, 1:st%nst))
590 st%mmb_iyoung(:,:) = -1
591 safe_allocate(st%mmb_proj(1:st%nst))
603 call parse_variable(namespace,
'SymmetrizeDensity', kpoints%use_symmetries, st%symmetrize_density)
633 integer,
intent(out) :: nik
634 integer,
intent(out) :: dim
635 integer,
intent(out) :: nst
636 integer,
intent(out) :: ierr
638 character(len=256) :: lines(3)
639 character(len=20) :: char
649 read(lines(1), *) char, nst
650 read(lines(2), *) char, dim
651 read(lines(3), *) char, nik
671 real(real64),
intent(in) :: excess_charge
674 integer :: ik, ist, ispin, nspin, ncols, nrows, el_per_state, icol, start_pos, spin_n
676 real(real64) :: rr, charge
677 logical :: integral_occs, unoccupied_states
678 real(real64),
allocatable :: read_occs(:, :)
679 real(real64) :: charge_in_block
692 call parse_variable(namespace,
'RestartFixedOccupations', .
true., st%restart_fixed_occ)
763 integral_occs = .
true.
765 occ_fix:
if (
parse_block(namespace,
'Occupations', blk) == 0)
then
767 st%fixed_occ = .
true.
770 if (ncols > st%nst)
then
775 if (nrows /= st%nik)
then
776 call messages_input_error(namespace,
"Occupations",
"Wrong number of rows in block Occupations.")
779 do ik = 1, st%nik - 1
782 "All rows in block Occupations must have the same number of columns.")
793 safe_allocate(read_occs(1:ncols, 1:st%nik))
799 charge_in_block = charge_in_block + read_occs(icol, ik) * st%kweights(ik)
804 select case (st%d%ispin)
813 start_pos = nint((st%qtot - charge_in_block)/spin_n)
815 if (start_pos + ncols > st%nst)
then
816 message(1) =
"To balance charge, the first column in block Occupations is taken to refer to state"
817 write(
message(2),
'(a,i6,a)')
"number ", start_pos,
" but there are too many columns for the number of states."
818 write(
message(3),
'(a,i6,a)')
"Solution: set ExtraStates = ", start_pos + ncols - st%nst
823 do ist = 1, start_pos
824 st%occ(ist, ik) = el_per_state
829 do ist = start_pos + 1, start_pos + ncols
830 st%occ(ist, ik) = read_occs(ist - start_pos, ik)
831 integral_occs = integral_occs .and. &
832 abs((st%occ(ist, ik) - el_per_state) * st%occ(ist, ik)) <=
m_epsilon
837 do ist = start_pos + ncols + 1, st%nst
844 safe_deallocate_a(read_occs)
847 st%fixed_occ = .false.
848 integral_occs = .false.
855 st%qtot = -(st%val_charge + excess_charge)
858 if (st%d%nspin == 2) nspin = 2
860 do ik = 1, st%nik, nspin
863 do ispin = ik, ik + nspin - 1
864 st%occ(ist, ispin) = min(rr, -(st%val_charge + excess_charge) - charge)
865 charge = charge + st%occ(ist, ispin)
884 if (st%fixed_occ)
then
885 call parse_variable(namespace,
'RestartReorderOccs', .false., st%restart_reorder_occs)
887 st%restart_reorder_occs = .false.
890 call smear_init(st%smear, namespace, st%d%ispin, st%fixed_occ, integral_occs, kpoints)
892 unoccupied_states = (st%d%ispin /=
spinors .and. st%nst*2 > st%qtot) .or. (st%d%ispin ==
spinors .and. st%nst > st%qtot)
895 if (.not. unoccupied_states)
then
896 call messages_write(
'Smearing needs unoccupied states (via ExtraStates or TotalStates) to be useful.')
904 charge = charge + sum(st%occ(ist, 1:st%nik) * st%kweights(1:st%nik))
906 if (abs(charge - st%qtot) > 1e-6_real64)
then
907 message(1) =
"Initial occupations do not integrate to total charge."
908 write(
message(2),
'(6x,f12.6,a,f12.6)') charge,
' != ', st%qtot
935 st%fixed_spins = .false.
936 if (st%d%ispin /=
spinors)
then
975 spin_fix:
if (
parse_block(namespace,
'InitialSpins', blk) == 0)
then
983 st%fixed_spins = .
true.
985 st%spin(:, :, i) = st%spin(:, :, 1)
998 class(
mesh_t),
intent(in) :: mesh
999 type(
type_t),
optional,
intent(in) :: wfs_type
1000 logical,
optional,
intent(in) :: skip(:)
1001 logical,
optional,
intent(in) :: packed
1005 if (
present(wfs_type))
then
1007 st%wfs_type = wfs_type
1035 type(
mesh_t),
intent(in) :: mesh
1036 logical,
optional,
intent(in) :: verbose
1037 logical,
optional,
intent(in) :: skip(:)
1038 logical,
optional,
intent(in) :: packed
1040 integer :: ib, iqn, ist, istmin, istmax
1041 logical :: same_node, verbose_, packed_
1042 integer,
allocatable :: bstart(:), bend(:)
1046 safe_allocate(bstart(1:st%nst))
1047 safe_allocate(bend(1:st%nst))
1048 safe_allocate(st%group%iblock(1:st%nst))
1057 if (
present(skip))
then
1059 if (.not. skip(ist))
then
1067 if (
present(skip))
then
1068 do ist = st%nst, istmin, -1
1069 if (.not. skip(ist))
then
1076 if (
present(skip) .and. verbose_)
then
1086 st%group%nblocks = 0
1088 do ist = istmin, istmax
1091 st%group%iblock(ist) = st%group%nblocks + 1
1094 if (st%parallel_in_states .and. ist /= istmax)
then
1097 same_node = (st%node(ist + 1) == st%node(ist))
1100 if (ib == st%block_size .or. ist == istmax .or. .not. same_node)
then
1102 st%group%nblocks = st%group%nblocks + 1
1103 bend(st%group%nblocks) = ist
1104 if (ist /= istmax) bstart(st%group%nblocks + 1) = ist + 1
1108 safe_allocate(st%group%psib(1:st%group%nblocks, st%d%kpt%start:st%d%kpt%end))
1110 safe_allocate(st%group%block_is_local(1:st%group%nblocks, 1:st%nik))
1111 st%group%block_is_local = .false.
1112 st%group%block_start = -1
1113 st%group%block_end = -2
1115 do ib = 1, st%group%nblocks
1116 if (bstart(ib) >= st%st_start .and. bend(ib) <= st%st_end)
then
1117 if (st%group%block_start == -1) st%group%block_start = ib
1118 st%group%block_end = ib
1119 do iqn = st%d%kpt%start, st%d%kpt%end
1120 st%group%block_is_local(ib, iqn) = .
true.
1123 call dwfs_elec_init(st%group%psib(ib, iqn), st%d%dim, bstart(ib), bend(ib), mesh%np_part, iqn, &
1124 special=.
true., packed=packed_)
1126 call zwfs_elec_init(st%group%psib(ib, iqn), st%d%dim, bstart(ib), bend(ib), mesh%np_part, iqn, &
1127 special=.
true., packed=packed_)
1134 safe_allocate(st%group%block_range(1:st%group%nblocks, 1:2))
1135 safe_allocate(st%group%block_size(1:st%group%nblocks))
1137 st%group%block_range(1:st%group%nblocks, 1) = bstart(1:st%group%nblocks)
1138 st%group%block_range(1:st%group%nblocks, 2) = bend(1:st%group%nblocks)
1139 st%group%block_size(1:st%group%nblocks) = bend(1:st%group%nblocks) - bstart(1:st%group%nblocks) + 1
1141 st%group%block_initialized = .
true.
1143 safe_allocate(st%group%block_node(1:st%group%nblocks))
1145 assert(
allocated(st%node))
1146 assert(all(st%node >= 0) .and. all(st%node < st%mpi_grp%size))
1148 do ib = 1, st%group%nblocks
1149 st%group%block_node(ib) = st%node(st%group%block_range(ib, 1))
1150 assert(st%group%block_node(ib) == st%node(st%group%block_range(ib, 2)))
1156 do ib = 1, st%group%nblocks
1162 if (st%group%block_size(ib) > 0)
then
1192 safe_deallocate_a(bstart)
1193 safe_deallocate_a(bend)
1214 type(
grid_t),
intent(in) :: gr
1216 real(real64) :: fsize
1220 safe_allocate(st%rho(1:gr%np_part, 1:st%d%nspin))
1223 fsize = gr%np_part*8.0_real64*st%block_size
1235 class(
space_t),
intent(in) :: space
1236 class(
mesh_t),
intent(in) :: mesh
1240 if (.not.
allocated(st%current))
then
1241 safe_allocate(st%current(1:mesh%np_part, 1:space%dim, 1:st%d%nspin))
1245 if (.not.
allocated(st%current_para))
then
1246 safe_allocate(st%current_para(1:mesh%np_part, 1:space%dim, 1:st%d%nspin))
1250 if (.not.
allocated(st%current_dia))
then
1251 safe_allocate(st%current_dia(1:mesh%np_part, 1:space%dim, 1:st%d%nspin))
1255 if (.not.
allocated(st%current_mag))
then
1256 safe_allocate(st%current_mag(1:mesh%np_part, 1:space%dim, 1:st%d%nspin))
1260 if (.not.
allocated(st%current_kpt))
then
1261 safe_allocate(st%current_kpt(1:mesh%np,1:space%dim,st%d%kpt%start:st%d%kpt%end))
1339 default = option__statesorthogonalization__cholesky_serial
1340#ifdef HAVE_SCALAPACK
1342 default = option__statesorthogonalization__cholesky_parallel
1346 call parse_variable(namespace,
'StatesOrthogonalization', default, st%orth_method)
1367 call parse_variable(namespace,
'StatesCLDeviceMemory', -512.0_real64, st%cl_states_mem)
1376 subroutine states_elec_copy(stout, stin, exclude_wfns, exclude_eigenval, special)
1379 logical,
optional,
intent(in) :: exclude_wfns
1380 logical,
optional,
intent(in) :: exclude_eigenval
1381 logical,
optional,
intent(in) :: special
1383 logical :: exclude_wfns_
1392 safe_allocate_source_a(stout%kweights, stin%kweights)
1393 stout%nik = stin%nik
1396 if (stin%modelmbparticles%nparticle > 0)
then
1397 safe_allocate_source_a(stout%mmb_nspindown, stin%mmb_nspindown)
1398 safe_allocate_source_a(stout%mmb_iyoung, stin%mmb_iyoung)
1399 safe_allocate_source_a(stout%mmb_proj, stin%mmb_proj)
1402 stout%wfs_type = stin%wfs_type
1403 stout%nst = stin%nst
1404 stout%block_size = stin%block_size
1405 stout%orth_method = stin%orth_method
1407 stout%cl_states_mem = stin%cl_states_mem
1408 stout%pack_states = stin%pack_states
1411 stout%only_userdef_istates = stin%only_userdef_istates
1413 if (.not. exclude_wfns_)
then
1414 safe_allocate_source_a(stout%rho, stin%rho)
1417 stout%uniform_occ = stin%uniform_occ
1420 safe_allocate_source_a(stout%eigenval, stin%eigenval)
1421 safe_allocate_source_a(stout%occ, stin%occ)
1422 safe_allocate_source_a(stout%spin, stin%spin)
1427 stout%group%nblocks = stin%group%nblocks
1429 safe_allocate_source_a(stout%user_def_states, stin%user_def_states)
1431 safe_allocate_source_a(stout%current, stin%current)
1432 safe_allocate_source_a(stout%current_kpt, stin%current_kpt)
1433 safe_allocate_source_a(stout%rho_core, stin%rho_core)
1434 safe_allocate_source_a(stout%frozen_rho, stin%frozen_rho)
1435 safe_allocate_source_a(stout%frozen_tau, stin%frozen_tau)
1436 safe_allocate_source_a(stout%frozen_gdens, stin%frozen_gdens)
1437 safe_allocate_source_a(stout%frozen_ldens, stin%frozen_ldens)
1439 stout%fixed_occ = stin%fixed_occ
1440 stout%restart_fixed_occ = stin%restart_fixed_occ
1442 stout%fixed_spins = stin%fixed_spins
1444 stout%qtot = stin%qtot
1445 stout%val_charge = stin%val_charge
1449 stout%parallel_in_states = stin%parallel_in_states
1451 call mpi_grp_copy(stout%dom_st_kpt_mpi_grp, stin%dom_st_kpt_mpi_grp)
1452 call mpi_grp_copy(stout%st_kpt_mpi_grp, stin%st_kpt_mpi_grp)
1453 call mpi_grp_copy(stout%dom_st_mpi_grp, stin%dom_st_mpi_grp)
1454 safe_allocate_source_a(stout%node, stin%node)
1455 safe_allocate_source_a(stout%st_kpt_task, stin%st_kpt_task)
1457#ifdef HAVE_SCALAPACK
1463 stout%scalapack_compatible = stin%scalapack_compatible
1465 stout%lnst = stin%lnst
1466 stout%st_start = stin%st_start
1467 stout%st_end = stin%st_end
1471 stout%symmetrize_density = stin%symmetrize_density
1475 stout%packed = stin%packed
1477 stout%randomization = stin%randomization
1493 if (st%modelmbparticles%nparticle > 0)
then
1494 safe_deallocate_a(st%mmb_nspindown)
1495 safe_deallocate_a(st%mmb_iyoung)
1496 safe_deallocate_a(st%mmb_proj)
1503 safe_deallocate_a(st%user_def_states)
1505 safe_deallocate_a(st%rho)
1506 safe_deallocate_a(st%eigenval)
1508 safe_deallocate_a(st%current)
1509 safe_deallocate_a(st%current_para)
1510 safe_deallocate_a(st%current_dia)
1511 safe_deallocate_a(st%current_mag)
1512 safe_deallocate_a(st%current_kpt)
1513 safe_deallocate_a(st%rho_core)
1514 safe_deallocate_a(st%frozen_rho)
1515 safe_deallocate_a(st%frozen_tau)
1516 safe_deallocate_a(st%frozen_gdens)
1517 safe_deallocate_a(st%frozen_ldens)
1518 safe_deallocate_a(st%occ)
1519 safe_deallocate_a(st%spin)
1520 safe_deallocate_a(st%kweights)
1523#ifdef HAVE_SCALAPACK
1529 safe_deallocate_a(st%node)
1530 safe_deallocate_a(st%st_kpt_task)
1532 if (st%parallel_in_states)
then
1533 safe_deallocate_a(st%ap%schedule)
1545 class(
mesh_t),
intent(in) :: mesh
1547 integer,
optional,
intent(in) :: ist_start_
1548 integer,
optional,
intent(in) :: ist_end_
1549 integer,
optional,
intent(in) :: ikpt_start_
1550 integer,
optional,
intent(in) :: ikpt_end_
1551 logical,
optional,
intent(in) :: normalized
1554 integer :: ist, ik, id, ist_start, ist_end, jst, ikpt_start, ikpt_end
1555 complex(real64) :: alpha, beta
1556 real(real64),
allocatable :: dpsi(:, :)
1557 complex(real64),
allocatable :: zpsi(:, :), zpsi2(:)
1558 integer :: ikpoint, ip
1561 logical :: normalized_
1572 safe_allocate(dpsi(1:mesh%np, 1:st%d%dim))
1574 safe_allocate(zpsi(1:mesh%np, 1:st%d%dim))
1577 select case (st%d%ispin)
1580 do ik = ikpt_start, ikpt_end
1581 ikpoint = st%d%get_kpoint_index(ik)
1582 do ist = ist_start, ist_end
1586 pre_shift = mesh%pv%xlocal-1, &
1587 post_shift = mesh%pv%np_global - mesh%pv%xlocal - mesh%np + 1, &
1588 normalized = normalized)
1590 if(mesh%parallel_in_domains)
then
1596 call dmf_random(mesh, dpsi(:, 1), normalized = normalized)
1601 zpsi(ip,1) = cmplx(dpsi(ip,1),
m_zero, real64)
1610 pre_shift = mesh%pv%xlocal-1, &
1611 post_shift = mesh%pv%np_global - mesh%pv%xlocal - mesh%np + 1, &
1612 normalized = normalized)
1614 if(mesh%parallel_in_domains)
then
1620 call zmf_random(mesh, zpsi(:, 1), normalized = normalized)
1632 if (st%fixed_spins)
then
1634 do ik = ikpt_start, ikpt_end
1635 ikpoint = st%d%get_kpoint_index(ik)
1636 do ist = ist_start, ist_end
1640 pre_shift = mesh%pv%xlocal-1, &
1641 post_shift = mesh%pv%np_global - mesh%pv%xlocal - mesh%np + 1, &
1642 normalized = normalized)
1644 if(mesh%parallel_in_domains)
then
1650 call dmf_random(mesh, dpsi(:, 1), normalized = normalized)
1654 zpsi(ip,1) = cmplx(dpsi(ip,1),
m_zero, real64)
1660 pre_shift = mesh%pv%xlocal-1, &
1661 post_shift = mesh%pv%np_global - mesh%pv%xlocal - mesh%np + 1, &
1662 normalized = normalized)
1664 if(mesh%parallel_in_domains)
then
1670 call zmf_random(mesh, zpsi(:, 1), normalized = normalized)
1683 safe_allocate(zpsi2(1:mesh%np))
1684 do jst = ist_start, ist - 1
1686 zpsi(1:mesh%np, 1) = zpsi(1:mesh%np, 1) -
zmf_dotp(mesh, zpsi(:, 1), zpsi2)*zpsi2(1:mesh%np)
1688 safe_deallocate_a(zpsi2)
1691 zpsi(1:mesh%np, 2) = zpsi(1:mesh%np, 1)
1695 if (abs(alpha) >
m_zero)
then
1696 beta = cmplx(st%spin(1, ist, ik) / abs(alpha), st%spin(2, ist, ik) / abs(alpha), real64)
1698 zpsi(1:mesh%np, 1) = alpha*zpsi(1:mesh%np, 1)
1699 zpsi(1:mesh%np, 2) = beta*zpsi(1:mesh%np, 2)
1704 do ik = ikpt_start, ikpt_end
1705 do ist = ist_start, ist_end
1709 pre_shift = mesh%pv%xlocal-1, &
1710 post_shift = mesh%pv%np_global - mesh%pv%xlocal - mesh%np + 1, &
1711 normalized = .false.)
1713 if(mesh%parallel_in_domains)
then
1719 call zmf_random(mesh, zpsi(:, id), normalized = .false.)
1734 safe_deallocate_a(dpsi)
1735 safe_deallocate_a(zpsi)
1746 class(
mesh_t),
intent(in) :: mesh
1747 logical,
optional,
intent(in) :: compute_spin
1751 real(real64) :: charge
1752 complex(real64),
allocatable :: zpsi(:, :)
1757 st%nik, st%nst, st%kweights)
1764 charge = charge + sum(st%occ(ist, 1:st%nik) * st%kweights(1:st%nik))
1766 if (abs(charge-st%qtot) > 1e-6_real64)
then
1767 message(1) =
'Occupations do not integrate to total charge.'
1768 write(
message(2),
'(6x,f12.8,a,f12.8)') charge,
' != ', st%qtot
1771 message(1) =
"There don't seem to be any electrons at all!"
1781 safe_allocate(zpsi(1:mesh%np, st%d%dim))
1782 do ik = st%d%kpt%start, st%d%kpt%end
1783 do ist = st%st_start, st%st_end
1785 st%spin(1:3, ist, ik) =
state_spin(mesh, zpsi)
1788 safe_deallocate_a(zpsi)
1790 if (st%parallel_in_states .or. st%d%kpt%parallel)
then
1805 real(real64),
optional,
intent(in) :: alt_eig(st%st_start:, st%d%kpt%start:)
1814 do ik = st%d%kpt%start, st%d%kpt%end
1815 if (
present(alt_eig))
then
1816 tot = tot + st%kweights(ik) * sum(st%occ(st%st_start:st%st_end, ik) * &
1817 alt_eig(st%st_start:st%st_end, ik))
1819 tot = tot + st%kweights(ik) * sum(st%occ(st%st_start:st%st_end, ik) * &
1820 st%eigenval(st%st_start:st%st_end, ik))
1824 if (st%parallel_in_states .or. st%d%kpt%parallel)
call comm_allreduce(st%st_kpt_mpi_grp, tot)
1836 logical :: default_scalapack_compatible
1847 st%parallel_in_states = .false.
1850 call mpi_grp_init(st%dom_st_kpt_mpi_grp, mc%dom_st_kpt_comm)
1854 default_scalapack_compatible =
calc_mode_par%scalapack_compat() .and. .not. st%d%kpt%parallel
1868 call parse_variable(namespace,
'ScaLAPACKCompatible', default_scalapack_compatible, st%scalapack_compatible)
1870#ifdef HAVE_SCALAPACK
1871 if (default_scalapack_compatible .neqv. st%scalapack_compatible)
then
1872 call messages_experimental(
'Setting ScaLAPACKCompatible to other than default', namespace=namespace)
1875 if (st%scalapack_compatible)
then
1882 st%scalapack_compatible = .false.
1891 if (st%nst < st%mpi_grp%size)
then
1892 message(1) =
"Have more processors than necessary"
1893 write(
message(2),
'(i4,a,i4,a)') st%mpi_grp%size,
" processors and ", st%nst,
" states."
1897 call distributed_init(st%dist, st%nst, st%mpi_grp%comm,
"states", scalapack_compat = st%scalapack_compatible)
1899 st%parallel_in_states = st%dist%parallel
1902 st%st_start = st%dist%start
1903 st%st_end = st%dist%end
1904 st%lnst = st%dist%nlocal
1905 st%node(1:st%nst) = st%dist%node(1:st%nst)
1922 kinetic_energy_density, paramagnetic_current, density_gradient, density_laplacian, &
1923 gi_kinetic_energy_density, st_end)
1927 logical,
intent(in) :: nlcc
1928 real(real64),
optional,
target,
intent(out) :: kinetic_energy_density(:,:)
1929 real(real64),
optional,
target,
intent(out) :: paramagnetic_current(:,:,:)
1930 real(real64),
optional,
intent(out) :: density_gradient(:,:,:)
1931 real(real64),
contiguous,
optional,
intent(out) :: density_laplacian(:,:)
1932 real(real64),
optional,
intent(out) :: gi_kinetic_energy_density(:,:)
1933 integer,
optional,
intent(in) :: st_end
1935 real(real64),
pointer :: jp(:, :, :)
1936 real(real64),
pointer :: tau(:, :)
1937 complex(real64),
allocatable :: wf_psi(:,:), gwf_psi(:,:,:), wf_psi_conj(:,:), lwf_psi(:,:)
1938 real(real64),
allocatable :: abs_wf_psi(:,:), abs_gwf_psi(:,:)
1939 complex(real64),
allocatable :: psi_gpsi(:,:)
1940 complex(real64) :: c_tmp
1941 integer :: is, ik, ist, i_dim, st_dim, ii, st_end_
1942 real(real64) :: ww, kpoint(gr%der%dim)
1943 logical :: something_to_do
1951 something_to_do =
present(kinetic_energy_density) .or.
present(gi_kinetic_energy_density) .or. &
1952 present(paramagnetic_current) .or.
present(density_gradient) .or.
present(density_laplacian)
1953 assert(something_to_do)
1955 safe_allocate( wf_psi(1:gr%np_part, 1:st%d%dim))
1956 safe_allocate( wf_psi_conj(1:gr%np_part, 1:st%d%dim))
1957 safe_allocate(gwf_psi(1:gr%np, 1:gr%der%dim, 1:st%d%dim))
1958 safe_allocate(abs_wf_psi(1:gr%np, 1:st%d%dim))
1959 safe_allocate(abs_gwf_psi(1:gr%np, 1:st%d%dim))
1960 safe_allocate(psi_gpsi(1:gr%np, 1:st%d%dim))
1961 if (
present(density_laplacian))
then
1962 safe_allocate(lwf_psi(1:gr%np, 1:st%d%dim))
1966 if (
present(kinetic_energy_density)) tau => kinetic_energy_density
1969 if (
present(paramagnetic_current)) jp => paramagnetic_current
1973 if (
present(gi_kinetic_energy_density))
then
1975 safe_allocate(jp(1:gr%np, 1:gr%der%dim, 1:st%d%nspin))
1977 if (.not.
present(kinetic_energy_density))
then
1978 safe_allocate(tau(1:gr%np, 1:st%d%nspin))
1982 if (
associated(tau)) tau =
m_zero
1983 if (
associated(jp)) jp =
m_zero
1984 if (
present(density_gradient)) density_gradient(:,:,:) =
m_zero
1985 if (
present(density_laplacian)) density_laplacian(:,:) =
m_zero
1986 if (
present(gi_kinetic_energy_density)) gi_kinetic_energy_density =
m_zero
1988 do ik = st%d%kpt%start, st%d%kpt%end
1990 kpoint(1:gr%der%dim) = kpoints%get_point(st%d%get_kpoint_index(ik))
1991 is = st%d%get_spin_index(ik)
1993 do ist = st%st_start, st_end_
1994 ww = st%kweights(ik)*st%occ(ist, ik)
2000 do st_dim = 1, st%d%dim
2005 do st_dim = 1, st%d%dim
2006 call zderivatives_grad(gr%der, wf_psi(:,st_dim), gwf_psi(:,:,st_dim), set_bc = .false.)
2010 if (
present(density_laplacian))
then
2011 do st_dim = 1, st%d%dim
2012 call zderivatives_lapl(gr%der, wf_psi(:,st_dim), lwf_psi(:,st_dim), ghost_update = .false., set_bc = .false.)
2017 wf_psi_conj(1:gr%np, 1:st%d%dim) = conjg(wf_psi(1:gr%np,1:st%d%dim))
2018 abs_wf_psi(1:gr%np, 1:st%d%dim) = real(wf_psi_conj(1:gr%np, 1:st%d%dim) * wf_psi(1:gr%np, 1:st%d%dim), real64)
2020 if (
present(density_laplacian))
then
2021 density_laplacian(1:gr%np, is) = density_laplacian(1:gr%np, is) + &
2022 ww *
m_two*real(wf_psi_conj(1:gr%np, 1) * lwf_psi(1:gr%np, 1), real64)
2023 if (st%d%ispin ==
spinors)
then
2024 density_laplacian(1:gr%np, 2) = density_laplacian(1:gr%np, 2) + &
2025 ww *
m_two*real(wf_psi_conj(1:gr%np, 2) * lwf_psi(1:gr%np, 2), real64)
2028 c_tmp = ww*(lwf_psi(ii, 1) * wf_psi_conj(ii, 2) + wf_psi(ii, 1) * conjg(lwf_psi(ii, 2)))
2029 density_laplacian(ii, 3) = density_laplacian(ii, 3) + real(c_tmp, real64)
2030 density_laplacian(ii, 4) = density_laplacian(ii, 4) + aimag(c_tmp)
2035 if (
associated(tau))
then
2036 tau(1:gr%np, is) = tau(1:gr%np, is) &
2037 + ww * sum(kpoint(1:gr%der%dim)**2) * abs_wf_psi(1:gr%np, 1)
2038 if (st%d%ispin ==
spinors)
then
2039 tau(1:gr%np, 2) = tau(1:gr%np, 2) &
2040 + ww * sum(kpoint(1:gr%der%dim)**2) * abs_wf_psi(1:gr%np, 2)
2044 c_tmp = ww * sum(kpoint(1:gr%der%dim)**2) * wf_psi(ii, 1) * wf_psi_conj(ii, 2)
2045 tau(ii, 3) = tau(ii, 3) + real(c_tmp, real64)
2046 tau(ii, 4) = tau(ii, 4) + aimag(c_tmp)
2051 do i_dim = 1, gr%der%dim
2054 psi_gpsi(1:gr%np, 1:st%d%dim) = wf_psi_conj(1:gr%np, 1:st%d%dim) &
2055 * gwf_psi(1:gr%np,i_dim,1:st%d%dim)
2056 abs_gwf_psi(1:gr%np, 1:st%d%dim) = real(conjg(gwf_psi(1:gr%np, i_dim, 1:st%d%dim)) &
2057 * gwf_psi(1:gr%np, i_dim, 1:st%d%dim), real64)
2059 if (
present(density_gradient))
then
2060 density_gradient(1:gr%np, i_dim, is) = density_gradient(1:gr%np, i_dim, is) &
2061 + ww *
m_two * real(psi_gpsi(1:gr%np, 1), real64)
2062 if (st%d%ispin ==
spinors)
then
2063 density_gradient(1:gr%np, i_dim, 2) = density_gradient(1:gr%np, i_dim, 2) &
2064 + ww *
m_two*real(psi_gpsi(1:gr%np, 2), real64)
2067 c_tmp = ww * (gwf_psi(ii, i_dim, 1) * wf_psi_conj(ii, 2) + wf_psi(ii, 1) * conjg(gwf_psi(ii, i_dim, 2)))
2068 density_gradient(ii, i_dim, 3) = density_gradient(ii, i_dim, 3) + real(c_tmp, real64)
2069 density_gradient(ii, i_dim, 4) = density_gradient(ii, i_dim, 4) + aimag(c_tmp)
2074 if (
present(density_laplacian))
then
2075 call lalg_axpy(gr%np, ww*
m_two, abs_gwf_psi(:,1), density_laplacian(:,is))
2076 if (st%d%ispin ==
spinors)
then
2077 call lalg_axpy(gr%np, ww*
m_two, abs_gwf_psi(:,2), density_laplacian(:,2))
2080 c_tmp =
m_two * ww * gwf_psi(ii, i_dim, 1) * conjg(gwf_psi(ii, i_dim, 2))
2081 density_laplacian(ii, 3) = density_laplacian(ii, 3) + real(c_tmp, real64)
2082 density_laplacian(ii, 4) = density_laplacian(ii, 4) + aimag(c_tmp)
2090 if (
associated(jp))
then
2092 jp(1:gr%np, i_dim, is) = jp(1:gr%np, i_dim, is) + &
2093 ww*(aimag(psi_gpsi(1:gr%np, 1)) - abs_wf_psi(1:gr%np, 1)*kpoint(i_dim))
2094 if (st%d%ispin ==
spinors)
then
2095 jp(1:gr%np, i_dim, 2) = jp(1:gr%np, i_dim, 2) + &
2096 ww*( aimag(psi_gpsi(1:gr%np, 2)) - abs_wf_psi(1:gr%np, 2)*kpoint(i_dim))
2099 c_tmp = -ww*
m_half*
m_zi*(gwf_psi(ii, i_dim, 1)*wf_psi_conj(ii, 2) - wf_psi(ii, 1)*conjg(gwf_psi(ii, i_dim, 2)) &
2100 -
m_two *
m_zi*wf_psi(ii, 1)*wf_psi_conj(ii, 2)*kpoint(i_dim))
2101 jp(ii, i_dim, 3) = jp(ii, i_dim, 3) + real(c_tmp, real64)
2102 jp(ii, i_dim, 4) = jp(ii, i_dim, 4) + aimag(c_tmp)
2111 if (
associated(tau))
then
2112 tau(1:gr%np, is) = tau(1:gr%np, is) + ww*(abs_gwf_psi(1:gr%np,1) &
2113 -
m_two*aimag(psi_gpsi(1:gr%np, 1))*kpoint(i_dim))
2114 if (st%d%ispin ==
spinors)
then
2115 tau(1:gr%np, 2) = tau(1:gr%np, 2) + ww*(abs_gwf_psi(1:gr%np, 2) &
2116 -
m_two*aimag(psi_gpsi(1:gr%np, 2))*kpoint(i_dim))
2119 c_tmp = ww * ( gwf_psi(ii, i_dim, 1)*conjg(gwf_psi(ii, i_dim, 2)) &
2120 +
m_zi * (gwf_psi(ii,i_dim,1)*wf_psi_conj(ii, 2) &
2121 - wf_psi(ii, 1)*conjg(gwf_psi(ii,i_dim,2)))*kpoint(i_dim))
2122 tau(ii, 3) = tau(ii, 3) + real(c_tmp, real64)
2123 tau(ii, 4) = tau(ii, 4) + aimag(c_tmp)
2133 safe_deallocate_a(wf_psi_conj)
2134 safe_deallocate_a(abs_wf_psi)
2135 safe_deallocate_a(abs_gwf_psi)
2136 safe_deallocate_a(psi_gpsi)
2138 if (.not.
present(gi_kinetic_energy_density))
then
2139 if (.not.
present(paramagnetic_current))
then
2140 safe_deallocate_p(jp)
2142 if (.not.
present(kinetic_energy_density))
then
2143 safe_deallocate_p(tau)
2147 if (st%parallel_in_states .or. st%d%kpt%parallel)
call reduce_all(st%st_kpt_mpi_grp)
2152 if (st%symmetrize_density)
then
2153 do is = 1, st%d%nspin
2154 if (
associated(tau))
then
2158 if (
present(density_laplacian))
then
2162 if (
associated(jp))
then
2166 if (
present(density_gradient))
then
2173 if (
allocated(st%rho_core) .and. nlcc .and. (
present(density_laplacian) .or.
present(density_gradient)))
then
2175 wf_psi(ii, 1) = st%rho_core(ii)/st%d%spin_channels
2180 if (
present(density_gradient))
then
2183 do is = 1, st%d%spin_channels
2184 density_gradient(1:gr%np, 1:gr%der%dim, is) = density_gradient(1:gr%np, 1:gr%der%dim, is) + &
2185 real(gwf_psi(1:gr%np, 1:gr%der%dim,1))
2190 if (
present(density_laplacian))
then
2193 do is = 1, st%d%spin_channels
2194 density_laplacian(1:gr%np, is) = density_laplacian(1:gr%np, is) + real(lwf_psi(1:gr%np, 1))
2201 if (
allocated(st%frozen_tau) .and. .not.
present(st_end))
then
2204 if (
allocated(st%frozen_gdens) .and. .not.
present(st_end))
then
2205 do is = 1, st%d%nspin
2206 call lalg_axpy(gr%np, gr%der%dim,
m_one, st%frozen_gdens(:,:,is), density_gradient(:,:,is))
2209 if (
allocated(st%frozen_tau) .and. .not.
present(st_end))
then
2210 call lalg_axpy(gr%np, st%d%nspin,
m_one, st%frozen_ldens, density_laplacian)
2213 safe_deallocate_a(wf_psi)
2214 safe_deallocate_a(gwf_psi)
2215 safe_deallocate_a(lwf_psi)
2219 if (
present(gi_kinetic_energy_density))
then
2220 do is = 1, st%d%nspin
2221 assert(
associated(tau))
2222 gi_kinetic_energy_density(1:gr%np, is) = tau(1:gr%np, is)
2225 assert(
associated(jp))
2226 if (st%d%ispin /=
spinors)
then
2227 do is = 1, st%d%nspin
2230 if (st%rho(ii, is) < 1.0e-7_real64) cycle
2231 gi_kinetic_energy_density(ii, is) = &
2232 gi_kinetic_energy_density(ii, is) - sum(jp(ii,1:gr%der%dim, is)**2)/st%rho(ii, is)
2238 gi_kinetic_energy_density(ii, 1) = gi_kinetic_energy_density(ii, 1) &
2239 - sum(jp(ii,1:gr%der%dim, 1)**2)/(safe_tol(st%rho(ii, 1),
m_epsilon))
2240 gi_kinetic_energy_density(ii, 2) = gi_kinetic_energy_density(ii, 2) &
2241 - sum(jp(ii,1:gr%der%dim, 2)**2)/(safe_tol(st%rho(ii, 2),
m_epsilon))
2242 gi_kinetic_energy_density(ii, 3) = &
2243 gi_kinetic_energy_density(ii, 3) - sum(jp(ii,1:gr%der%dim, 3)**2 + jp(ii,1:gr%der%dim, 4)**2) &
2244 /(safe_tol((st%rho(ii, 3)**2 + st%rho(ii, 4)**2),
m_epsilon))*st%rho(ii, 3)
2245 gi_kinetic_energy_density(ii, 4) = &
2246 gi_kinetic_energy_density(ii, 4) + sum(jp(ii,1:gr%der%dim, 3)**2 + jp(ii,1:gr%der%dim, 4)**2) &
2247 /(safe_tol((st%rho(ii, 3)**2 + st%rho(ii, 4)**2),
m_epsilon))*st%rho(ii, 4)
2253 if (.not.
present(kinetic_energy_density))
then
2254 safe_deallocate_p(tau)
2256 if (.not.
present(paramagnetic_current))
then
2257 safe_deallocate_p(jp)
2272 if (
associated(tau))
call comm_allreduce(grp, tau, dim = (/gr%np, st%d%nspin/))
2274 if (
present(density_laplacian))
call comm_allreduce(grp, density_laplacian, dim = (/gr%np, st%d%nspin/))
2276 do is = 1, st%d%nspin
2277 if (
associated(jp))
call comm_allreduce(grp, jp(:, :, is), dim = (/gr%np, gr%der%dim/))
2279 if (
present(density_gradient))
then
2280 call comm_allreduce(grp, density_gradient(:, :, is), dim = (/gr%np, gr%der%dim/))
2294 type(
mesh_t),
intent(in) :: mesh
2295 complex(real64),
intent(in) :: f1(:, :)
2296 real(real64) :: spin(1:3)
2298 complex(real64) :: z
2302 z =
zmf_dotp(mesh, f1(:, 1) , f1(:, 2))
2304 spin(1) =
m_two*real(z, real64)
2305 spin(2) =
m_two*aimag(z)
2317 integer,
intent(in) :: ist
2331 integer,
intent(in) :: ist
2332 integer,
intent(in) :: ik
2337 ik >= st%d%kpt%start .and. ik <= st%d%kpt%end
2347 class(
mesh_t),
intent(in) :: mesh
2353 memory = memory + real64*real(mesh%np_part_global, real64) *st%d%dim*real(st%nst, real64) *st%d%kpt%nglobal
2363 logical,
optional,
intent(in) :: copy
2366 integer(int64) :: max_mem, mem
2378 if (accel_is_enabled())
then
2379 max_mem = accel_global_memory_size()
2381 if (st%cl_states_mem > m_one)
then
2382 max_mem = int(st%cl_states_mem, int64)*(1024_8)**2
2383 else if (st%cl_states_mem < 0.0_real64)
then
2384 max_mem = max_mem + int(st%cl_states_mem, int64)*(1024_8)**2
2386 max_mem = int(st%cl_states_mem*real(max_mem, real64) , int64)
2389 max_mem = huge(max_mem)
2393 qnloop:
do iqn = st%d%kpt%start, st%d%kpt%end
2394 do ib = st%group%block_start, st%group%block_end
2396 mem = mem + st%group%psib(ib, iqn)%pack_total_size()
2398 if (mem > max_mem)
then
2399 call messages_write(
'Not enough CL device memory to store all states simultaneously.', new_line = .
true.)
2400 call messages_write(
'Only ')
2401 call messages_write(ib - st%group%block_start)
2402 call messages_write(
' of ')
2403 call messages_write(st%group%block_end - st%group%block_start + 1)
2404 call messages_write(
' blocks will be stored in device memory.', new_line = .
true.)
2405 call messages_warning()
2409 call st%group%psib(ib, iqn)%do_pack(copy)
2421 logical,
optional,
intent(in) :: copy
2430 do iqn = st%d%kpt%start, st%d%kpt%end
2431 do ib = st%group%block_start, st%group%block_end
2432 if (st%group%psib(ib, iqn)%is_packed())
call st%group%psib(ib, iqn)%do_unpack(copy)
2445 type(namespace_t),
intent(in) :: namespace
2449 call messages_print_with_emphasis(msg=
"States", namespace=namespace)
2451 write(message(1),
'(a,f12.3)')
'Total electronic charge = ', st%qtot
2452 write(message(2),
'(a,i8)')
'Number of states = ', st%nst
2453 write(message(3),
'(a,i8)')
'States block-size = ', st%block_size
2454 call messages_info(3, namespace=namespace)
2456 call messages_print_with_emphasis(namespace=namespace)
2471 do iqn = st%d%kpt%start, st%d%kpt%end
2472 do ib = st%group%block_start, st%group%block_end
2473 call batch_set_zero(st%group%psib(ib, iqn))
2483 integer pure function states_elec_block_min(st, ib) result(range)
2485 integer,
intent(in) :: ib
2487 range = st%group%block_range(ib, 1)
2493 integer pure function states_elec_block_max(st, ib) result(range)
2495 integer,
intent(in) :: ib
2497 range = st%group%block_range(ib, 2)
2503 integer pure function states_elec_block_size(st, ib) result(size)
2505 integer,
intent(in) :: ib
2507 size = st%group%block_size(ib)
2515 type(namespace_t),
intent(in) :: namespace
2516 integer,
intent(out) :: n_pairs
2517 integer,
intent(out) :: n_occ(:)
2518 integer,
intent(out) :: n_unocc(:)
2519 logical,
allocatable,
intent(out) :: is_included(:,:,:)
2521 logical,
intent(out) :: is_frac_occ
2523 integer :: ik, ist, ast, n_filled, n_partially_filled, n_half_filled
2524 character(len=80) :: nst_string, default, wfn_list
2525 real(real64) :: energy_window
2529 is_frac_occ = .false.
2531 call occupied_states(st, namespace, ik, n_filled, n_partially_filled, n_half_filled)
2532 if (n_partially_filled > 0 .or. n_half_filled > 0) is_frac_occ = .
true.
2533 n_occ(ik) = n_filled + n_partially_filled + n_half_filled
2534 n_unocc(ik) = st%nst - n_filled
2547 call parse_variable(namespace,
'CasidaKSEnergyWindow', -m_one, energy_window, units_inp%energy)
2570 safe_allocate(is_included(maxval(n_occ), minval(n_occ) + 1:st%nst , st%nik))
2571 is_included(:,:,:) = .false.
2573 if (energy_window < m_zero)
then
2574 write(nst_string,
'(i6)') st%nst
2575 write(default,
'(a,a)')
"1-", trim(adjustl(nst_string))
2576 call parse_variable(namespace,
'CasidaKohnShamStates', default, wfn_list)
2578 write(message(1),
'(a,a)')
"Info: States that form the basis: ", trim(wfn_list)
2579 call messages_info(1, namespace=namespace)
2584 do ast = n_occ(ik) + 1, st%nst
2585 if (loct_isinstringlist(ast, wfn_list))
then
2586 do ist = 1, n_occ(ik)
2587 if (loct_isinstringlist(ist, wfn_list))
then
2588 n_pairs = n_pairs + 1
2589 is_included(ist, ast, ik) = .
true.
2598 write(message(1),
'(a,f12.6,a)')
"Info: including transitions with energy < ", &
2599 units_from_atomic(units_out%energy, energy_window), trim(units_abbrev(units_out%energy))
2600 call messages_info(1, namespace=namespace)
2605 do ast = n_occ(ik) + 1, st%nst
2606 do ist = 1, n_occ(ik)
2607 if (st%eigenval(ast, ik) - st%eigenval(ist, ik) < energy_window)
then
2608 n_pairs = n_pairs + 1
2609 is_included(ist, ast, ik) = .
true.
2637 subroutine occupied_states(st, namespace, ik, n_filled, n_partially_filled, n_half_filled, &
2638 filled, partially_filled, half_filled)
2640 type(namespace_t),
intent(in) :: namespace
2641 integer,
intent(in) :: ik
2642 integer,
intent(out) :: n_filled, n_partially_filled, n_half_filled
2643 integer,
optional,
intent(out) :: filled(:), partially_filled(:), half_filled(:)
2649 if (
present(filled)) filled(:) = 0
2650 if (
present(partially_filled)) partially_filled(:) = 0
2651 if (
present(half_filled)) half_filled(:) = 0
2653 n_partially_filled = 0
2656 select case (st%d%ispin)
2659 if (abs(st%occ(ist, ik) - m_two) < m_min_occ)
then
2660 n_filled = n_filled + 1
2661 if (
present(filled)) filled(n_filled) = ist
2662 else if (abs(st%occ(ist, ik) - m_one) < m_min_occ)
then
2663 n_half_filled = n_half_filled + 1
2664 if (
present(half_filled)) half_filled(n_half_filled) = ist
2665 else if (st%occ(ist, ik) > m_min_occ)
then
2666 n_partially_filled = n_partially_filled + 1
2667 if (
present(partially_filled)) partially_filled(n_partially_filled) = ist
2668 else if (abs(st%occ(ist, ik)) > m_min_occ)
then
2669 write(message(1),*)
'Internal error in occupied_states: Illegal occupation value ', st%occ(ist, ik)
2670 call messages_fatal(1, namespace=namespace)
2673 case (spin_polarized, spinors)
2675 if (abs(st%occ(ist, ik)-m_one) < m_min_occ)
then
2676 n_filled = n_filled + 1
2677 if (
present(filled)) filled(n_filled) = ist
2678 else if (st%occ(ist, ik) > m_min_occ)
then
2679 n_partially_filled = n_partially_filled + 1
2680 if (
present(partially_filled)) partially_filled(n_partially_filled) = ist
2681 else if (abs(st%occ(ist, ik)) > m_min_occ)
then
2682 write(message(1),*)
'Internal error in occupied_states: Illegal occupation value ', st%occ(ist, ik)
2683 call messages_fatal(1, namespace=namespace)
2695 type(multicomm_t),
intent(in) :: mc
2698 call distributed_init(this%d%kpt, this%nik, mc%group_comm(p_strategy_kpoints),
"k-points")
2713 if (.not.
allocated(st%st_kpt_task))
then
2714 safe_allocate(st%st_kpt_task(0:st%st_kpt_mpi_grp%size-1, 1:4))
2717 st%st_kpt_task(0:st%st_kpt_mpi_grp%size-1, :) = 0
2718 st%st_kpt_task(st%st_kpt_mpi_grp%rank, :) = [st%st_start, st%st_end, st%d%kpt%start, st%d%kpt%end]
2720 if (st%parallel_in_states .or. st%d%kpt%parallel)
then
2721 call comm_allreduce(st%st_kpt_mpi_grp, st%st_kpt_task)
2732 type(kpoints_t),
intent(in) :: kpoints
2733 type(namespace_t),
intent(in) :: namespace
2736 type(states_elec_dim_t),
pointer :: dd
2743 st%nik = kpoints_number(kpoints)
2745 if (dd%ispin == spin_polarized) st%nik = 2*st%nik
2747 safe_allocate(st%kweights(1:st%nik))
2750 ik = dd%get_kpoint_index(iq)
2751 st%kweights(iq) = kpoints%get_weight(ik)
2763 call io_mkdir(
'debug/', namespace)
2764 iunit = io_open(
'debug/kpoints', namespace, action =
'write')
2765 call kpoints%write_info(iunit=iunit, absolute_coordinates = .
true.)
2766 call io_close(iunit)
2780 class(mesh_t),
intent(in) :: gr
2781 real(real64) :: dipole(1:gr%box%dim)
2784 real(real64) :: e_dip(1:gr%box%dim, this%d%spin_channels)
2788 do ispin = 1, this%d%spin_channels
2789 call dmf_dipole(gr, this%rho(:, ispin), e_dip(:, ispin))
2792 dipole(:) = sum(e_dip(:, 1:this%d%spin_channels), 2)
2800#include "states_elec_inc.F90"
2803#include "complex.F90"
2804#include "states_elec_inc.F90"
initialize a batch with existing memory
constant times a vector plus a vector
Prints out to iunit a message in the form: ["InputVariable" = value] where "InputVariable" is given b...
pure logical function, public accel_is_enabled()
type(accel_t), public accel
This module implements batches of mesh functions.
This module implements common operations on batches of mesh functions.
This module provides the BLACS processor grid.
subroutine, public blacs_proc_grid_init(this, mpi_grp, procdim)
Initializes a blacs context from an MPI communicator with topological information.
subroutine, public blacs_proc_grid_end(this)
subroutine, public blacs_proc_grid_copy(cin, cout)
Module implementing boundary conditions in Octopus.
This module handles the calculation mode.
type(calc_mode_par_t), public calc_mode_par
Singleton instance of parallel calculation mode.
integer, parameter, public p_strategy_states
parallelization in states
This module calculates the derivatives (gradients, Laplacians, etc.) of a function.
subroutine, public zderivatives_grad(der, ff, op_ff, ghost_update, set_bc, to_cartesian)
apply the gradient to a mesh function
subroutine, public zderivatives_lapl(der, ff, op_ff, ghost_update, set_bc, factor)
apply the Laplacian to a mesh function
subroutine, public distributed_end(this)
subroutine, public distributed_nullify(this, total)
subroutine, public distributed_init(this, total, comm, tag, scalapack_compat)
Distribute N instances across M processes of communicator comm
subroutine, public distributed_copy(in, out)
@Brief Create a copy of a distributed instance
integer, parameter, public unpolarized
Parameters...
integer, parameter, public spinors
integer, parameter, public spin_polarized
real(real64), parameter, public m_two
real(real64), parameter, public m_zero
real(real64), parameter, public m_fourth
complex(real64), parameter, public m_zi
real(real64), parameter, public m_epsilon
type(conf_t), public conf
Global instance of Octopus configuration.
real(real64), parameter, public m_half
real(real64), parameter, public m_one
This module implements the underlying real-space grid.
subroutine, public dgrid_symmetrize_vector_field(gr, field, suppress_warning)
subroutine, public dgrid_symmetrize_scalar_field(gr, field, suppress_warning)
logical pure function, public kpoints_point_is_gamma(this, ik)
This module is intended to contain "only mathematical" functions and procedures.
integer pure function, public pad_pow2(size)
create array size, which is padded to powers of 2
This module defines functions over batches of mesh functions.
subroutine, public zmesh_batch_exchange_points(mesh, aa, forward_map, backward_map)
This functions exchanges points of a mesh according to a certain map. Two possible maps can be given....
subroutine, public dmesh_batch_exchange_points(mesh, aa, forward_map, backward_map)
This functions exchanges points of a mesh according to a certain map. Two possible maps can be given....
This module defines various routines, operating on mesh functions.
subroutine, public zmf_random(mesh, ff, pre_shift, post_shift, seed, normalized)
This subroutine fills a function with random values.
subroutine, public dmf_random(mesh, ff, pre_shift, post_shift, seed, normalized)
This subroutine fills a function with random values.
subroutine, public zmf_normalize(mesh, dim, psi, norm)
Normalize a mesh function psi.
This module defines the meshes, which are used in Octopus.
subroutine, public messages_not_implemented(feature, namespace)
subroutine, public messages_warning(no_lines, all_nodes, namespace)
subroutine, public messages_obsolete_variable(namespace, name, rep)
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 messages_input_error(namespace, var, details, row, column)
subroutine, public messages_experimental(name, namespace)
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
general module for modelmb particles
subroutine, public modelmb_particles_end(this)
subroutine, public modelmb_particles_init(this, namespace, space)
==============================================================
subroutine, public modelmb_particles_copy(modelmb_out, modelmb_in)
subroutine mpi_grp_copy(mpi_grp_out, mpi_grp_in)
type(mpi_comm), parameter, public mpi_comm_undefined
used to indicate a communicator has not been initialized
subroutine mpi_grp_init(grp, comm)
Initialize MPI group instance.
This module handles the communicators for the various parallelization strategies.
logical pure function, public multicomm_strategy_is_parallel(mc, level)
subroutine, public multicomm_all_pairs_copy(apout, apin)
logical pure function, public multicomm_have_slaves(this)
subroutine, public multicomm_create_all_pairs(mpi_grp, ap)
This routine uses the one-factorization (or near-one-factorization of a complete graph to construct a...
logical function, public parse_is_defined(namespace, name)
integer function, public parse_block(namespace, name, blk, check_varinfo_)
subroutine, public profiling_out(label)
Increment out counter and sum up difference between entry and exit time.
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
subroutine, public restart_read(restart, iunit, lines, nlines, ierr)
subroutine, public restart_close(restart, iunit)
Close a file previously opened with restart_open.
integer function, public restart_open(restart, filename, status, position, silent)
Open file "filename" found inside the current restart directory. Depending on the type of restart,...
subroutine, public smear_fill_occupations(this, eigenvalues, occupations, nik, nst)
subroutine, public smear_find_fermi_energy(this, namespace, eigenvalues, occupations, qtot, nik, nst, kweights)
subroutine, public smear_copy(to, from)
integer, parameter, public smear_fixed_occ
subroutine, public smear_init(this, namespace, ispin, fixed_occ, integral_occs, kpoints)
logical pure function, public smear_is_semiconducting(this)
subroutine, public states_set_complex(st)
pure logical function, public states_are_complex(st)
pure logical function, public states_are_real(st)
This module handles spin dimensions of the states and the k-point distribution.
subroutine, public states_elec_dim_copy(dout, din)
subroutine, public states_elec_dim_end(dim)
This module handles groups of electronic batches and their parallel distribution.
subroutine, public states_elec_group_copy(d, group_in, group_out, copy_data, special)
make a copy of a group
subroutine, public states_elec_group_end(this, d)
finalize the local blocks of wave functions and release local arrays
real(real64) function, public states_elec_wfns_memory(st, mesh)
return the memory usage of a states_elec_t object
integer pure function, public states_elec_block_max(st, ib)
return index of last state in block ib
subroutine zstates_elec_generate_random_vector(mesh, st, vector, normalized, reset_seed)
Generate a random vector.
subroutine, public states_elec_null(st)
subroutine, public states_elec_distribute_nodes(st, namespace, mc)
@Brief. Distribute states over the processes for states parallelization
subroutine, public states_elec_fermi(st, namespace, mesh, compute_spin)
calculate the Fermi level for the states in this object
logical function, public state_kpt_is_local(st, ist, ik)
check whether a given state (ist, ik) is on the local node
subroutine, public states_elec_calc_quantities(gr, st, kpoints, nlcc, kinetic_energy_density, paramagnetic_current, density_gradient, density_laplacian, gi_kinetic_energy_density, st_end)
calculated selected quantities
real(real64) function, dimension(1:gr%box%dim) states_elec_calculate_dipole(this, gr)
calculate the expectation value of the dipole moment of electrons
subroutine, public states_elec_densities_init(st, gr)
subroutine, public states_elec_end(st)
finalize the states_elec_t object
subroutine zstates_elec_get_state1(st, mesh, idim, ist, iqn, psi)
Write one component (dim) of a wave function into a state_elec_t object.
subroutine, public states_elec_deallocate_wfns(st)
Deallocates the KS wavefunctions defined within a states_elec_t structure.
subroutine, public states_elec_allocate_wfns(st, mesh, wfs_type, skip, packed)
Allocates the KS wavefunctions defined within a states_elec_t structure.
subroutine zstates_elec_get_points1(st, start_point, end_point, iqn, psi)
Returns the value of all the states for given k-point in the range of points [start_point:end_point].
subroutine states_elec_read_initial_occs(st, namespace, excess_charge, kpoints)
Reads from the input file the initial occupations, if the block "Occupations" is present....
subroutine, public states_elec_init(st, namespace, space, valence_charge, kpoints)
Initialize a new states_elec_t object.
subroutine zstates_elec_set_state1(st, mesh, idim, ist, iqn, psi)
get one dimension of local wave function for given k-point and states index from a states_elec_t obje...
real(real64) function, public states_elec_eigenvalues_sum(st, alt_eig)
function to calculate the eigenvalues sum using occupations as weights
integer, parameter, public par_independent
Method used to generate random states.
subroutine dstates_elec_generate_random_vector(mesh, st, vector, normalized, reset_seed)
Generate a random vector.
subroutine, public occupied_states(st, namespace, ik, n_filled, n_partially_filled, n_half_filled, filled, partially_filled, half_filled)
return information about occupied orbitals in many-body state
subroutine dstates_elec_get_state2(st, mesh, ist, iqn, psi)
Write a wave function into a state_elec_t object.
subroutine zstates_elec_set_state2(st, mesh, ist, iqn, psi)
get local wave function for given k-point and states index from a states_elec_t object
subroutine zstates_elec_get_state2(st, mesh, ist, iqn, psi)
Write a wave function into a state_elec_t object.
subroutine dstates_elec_get_points1(st, start_point, end_point, iqn, psi)
Returns the value of all the states for given k-point in the range of points [start_point:end_point].
integer pure function, public states_elec_block_size(st, ib)
return number of states in block ib
subroutine, public kpoints_distribute(this, mc)
distribute k-points over the nodes in the corresponding communicator
integer, parameter, public par_dependent
subroutine states_elec_pack(st, copy)
pack the batches in this states object
subroutine, public states_elec_choose_kpoints(st, kpoints, namespace)
double up k-points for SPIN_POLARIZED calculations
subroutine zstates_elec_get_state3(st, mesh, iqn, psi)
get local wave functions for given k-point from a states_elec_t object
logical function, public state_is_local(st, ist)
check whether a given state (ist) is on the local node
subroutine dstates_elec_set_state1(st, mesh, idim, ist, iqn, psi)
get one dimension of local wave function for given k-point and states index from a states_elec_t obje...
subroutine, public states_elec_exec_init(st, namespace, mc)
Further initializations.
subroutine dstates_elec_set_state4(st, mesh, psi)
set all local wave functions in a states_elec_t object
subroutine dstates_elec_set_state3(st, mesh, iqn, psi)
set local wave functions for given k-point in a states_elec_t object
subroutine states_elec_write_info(st, namespace)
write information about the states object
subroutine states_elec_init_block(st, mesh, verbose, skip, packed)
Initializes the data components in st that describe how the states are distributed in blocks:
real(real64) function, dimension(1:3) state_spin(mesh, f1)
calculate the spin vector for a spinor wave function f1
subroutine states_elec_kpoints_distribution(st)
Assign the start and end indices for states and kpoints, for "st_kpt_mpi_grp" communicator.
subroutine zstates_elec_get_state4(st, mesh, psi)
get all local wave functions from a states_elec_t object
subroutine, public states_elec_copy(stout, stin, exclude_wfns, exclude_eigenval, special)
make a (selective) copy of a states_elec_t object
subroutine, public states_elec_generate_random(st, mesh, kpoints, ist_start_, ist_end_, ikpt_start_, ikpt_end_, normalized)
randomize states
subroutine dstates_elec_get_state1(st, mesh, idim, ist, iqn, psi)
Write one component (dim) of a wave function into a state_elec_t object.
subroutine dstates_elec_set_state2(st, mesh, ist, iqn, psi)
get local wave function for given k-point and states index from a states_elec_t object
subroutine zstates_elec_set_state4(st, mesh, psi)
set all local wave functions in a states_elec_t object
subroutine states_elec_read_initial_spins(st, namespace)
Reads, if present, the "InitialSpins" block.
subroutine dstates_elec_get_state3(st, mesh, iqn, psi)
get local wave functions for given k-point from a states_elec_t object
subroutine dstates_elec_get_points2(st, start_point, end_point, psi)
Returns the value of all the states in the range of points [start_point:end_point].
integer pure function, public states_elec_block_min(st, ib)
return index of first state in block ib
subroutine states_elec_unpack(st, copy)
unpack the batches in this states object
subroutine, public states_elec_look(restart, nik, dim, nst, ierr)
Reads the 'states' file in the restart directory, and finds out the nik, dim, and nst contained in it...
subroutine dstates_elec_get_state4(st, mesh, psi)
get all local wave functions from a states_elec_t object
subroutine, public states_elec_count_pairs(st, namespace, n_pairs, n_occ, n_unocc, is_included, is_frac_occ)
number of occupied-unoccupied pairs for Casida
subroutine, public states_elec_allocate_current(st, space, mesh)
subroutine, public states_elec_set_zero(st)
expclicitely set all wave functions in the states to zero
subroutine zstates_elec_set_state3(st, mesh, iqn, psi)
set local wave functions for given k-point in a states_elec_t object
subroutine zstates_elec_get_points2(st, start_point, end_point, psi)
Returns the value of all the states in the range of points [start_point:end_point].
type(type_t), public type_float
type(type_t), public type_cmplx
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
This module defines the unit system, used for input and output.
type(unit_t), public unit_megabytes
For large amounts of data (natural code units are bytes)
subroutine, public zwfs_elec_init(this, dim, st_start, st_end, np, ik, special, packed)
initialize an empty wfs_elec_t object
subroutine, public dwfs_elec_init(this, dim, st_start, st_end, np, ik, special, packed)
initialize an empty wfs_elec_t object
subroutine reduce_all(grp)
subroutine print_kpoints_debug
Class defining batches of mesh functions.
Distribution of N instances over mpi_grpsize processes, for the local rank mpi_grprank....
Extension of space that contains the knowledge of the spin dimension.
Description of the grid, containing information on derivatives, stencil, and symmetries.
Describes mesh distribution to nodes.
This is defined even when running serial.
An all-pairs communication schedule for a given group.
Stores all communicators and groups.
abstract class for states
The states_elec_t class contains all electronic wave functions.