Octopus
states_elec.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-2006 M. Marques, A. Castro, A. Rubio, G. Bertsch
2!!
3!! This program is free software; you can redistribute it and/or modify
4!! it under the terms of the GNU General Public License as published by
5!! the Free Software Foundation; either version 2, or (at your option)
6!! any later version.
7!!
8!! This program is distributed in the hope that it will be useful,
9!! but WITHOUT ANY WARRANTY; without even the implied warranty of
10!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
11!! GNU General Public License for more details.
12!!
13!! You should have received a copy of the GNU General Public License
14!! along with this program; if not, write to the Free Software
15!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
16!! 02110-1301, USA.
17!!
18#include "global.h"
19
21 use accel_oct_m
22 use batch_oct_m
26 use batch_oct_m
29 use comm_oct_m
30 use debug_oct_m
34 use global_oct_m
35 use grid_oct_m
36 use io_oct_m
37 use, intrinsic :: iso_fortran_env
40 use loct_oct_m
41 use math_oct_m
42 use mesh_oct_m
47 use mpi_oct_m
50#ifdef HAVE_OPENMP
51 use omp_lib
52#endif
53 use parser_oct_m
57 use smear_oct_m
58 use space_oct_m
62 use types_oct_m
63 use unit_oct_m
67
68 implicit none
69
70 private
71
72 public :: &
103 stress_t, &
105
106
107 ! this type must be moved to stress module but due to circular dependency it is not possible now
108 type stress_t
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
113 real(real64) :: xc_nlcc(3,3) = m_zero
114 real(real64) :: ps_local(3,3) = m_zero
115 real(real64) :: ps_nl(3,3) = m_zero
116 real(real64) :: ion_ion(3,3) = m_zero
117 real(real64) :: vdw(3,3) = m_zero
118 real(real64) :: hubbard(3,3) = m_zero
119
120 real(real64) :: kinetic_sumrule
121 real(real64) :: hartree_sumrule
122 end type stress_t
123
124 ! TODO(Alex) Issue #672 Decouple k-point info from `states_elec_dim_t`
125
133 !
134 type, extends(states_abst_t) :: states_elec_t
135 ! Components are public by default
136 type(states_elec_dim_t) :: d
137 integer :: nst_conv
138
139 logical :: only_userdef_istates
140
141 type(states_elec_group_t) :: group
142 integer :: block_size
146 logical :: pack_states
150
151
152 character(len=1024), allocatable :: user_def_states(:,:,:)
155
156 ! TODO(Alex) Issue #821. Collate current quantities into an object.
157 ! the densities and currents (after all we are doing DFT :)
158 real(real64), allocatable :: rho(:,:)
159 real(real64), allocatable :: rho_core(:)
160
161 real(real64), allocatable :: current(:, :, :)
162 real(real64), allocatable :: current_para(:, :, :)
163 real(real64), allocatable :: current_dia(:, :, :)
164 real(real64), allocatable :: current_mag(:, :, :)
165 real(real64), allocatable :: current_kpt(:,:,:)
166
167 ! TODO(Alex) Issue #673. Create frozen density class and replace in states_elec_t
168 ! It may be required to "freeze" the deepest orbitals during the evolution; the density
169 ! of these orbitals is kept in frozen_rho. It is different from rho_core.
170 real(real64), allocatable :: frozen_rho(:, :)
171 real(real64), allocatable :: frozen_tau(:, :)
172 real(real64), allocatable :: frozen_gdens(:,:,:)
173 real(real64), allocatable :: frozen_ldens(:,:)
174
175 logical :: uniform_occ
176
177 real(real64), allocatable :: eigenval(:,:)
178 logical :: fixed_occ
179 logical :: restart_fixed_occ
180 logical :: restart_reorder_occs
181 real(real64), allocatable :: occ(:,:)
182 real(real64), allocatable :: kweights(:)
183 integer :: nik
184
185 logical :: fixed_spins
187 real(real64), allocatable :: spin(:, :, :)
188
189 real(real64) :: qtot
190 real(real64) :: val_charge
191
192 type(stress_t) :: stress_tensors
193
194 logical :: fromScratch
195 type(smear_t) :: smear
196
197 ! TODO(Alex) Issue #823 Move modelmbparticles out of states_elec_t
198 type(modelmb_particle_t) :: modelmbparticles
199
200 ! TODO(Alex) Issue #824. Package the communicators in a single instance prior to removing
201 ! or consider creating a distributed_t instance for each (as distributed_t contains the an instance of mpi_grp_t)
202 type(mpi_grp_t) :: mpi_grp
203 type(mpi_grp_t) :: dom_st_mpi_grp
204 type(mpi_grp_t) :: st_kpt_mpi_grp
205 type(mpi_grp_t) :: dom_st_kpt_mpi_grp
206 type(blacs_proc_grid_t) :: dom_st_proc_grid
208 type(distributed_t) :: dist
209 logical :: scalapack_compatible
210 logical :: parallel_in_states = .false.
212 ! TODO(Alex) Issue #820. Remove lnst, st_start, st_end and node, as they are all contained within dist
213 integer :: lnst
214 integer :: st_start, st_end
215 integer, allocatable :: node(:)
217 ! TODO(Alex) Issue #824. Either change from data to a method, or package with `st_kpt_mpi_grp`
218 integer, allocatable :: st_kpt_task(:,:)
219
220 type(multicomm_all_pairs_t), private :: ap
221 logical :: symmetrize_density
222 integer :: randomization
223 integer :: orth_method = 0
224
225 real(real64) :: gpu_states_mem
226
227 contains
228 procedure :: nullify => states_elec_null
229 procedure :: write_info => states_elec_write_info
230 procedure :: pack => states_elec_pack
231 procedure :: unpack => states_elec_unpack
233 procedure :: dipole => states_elec_calculate_dipole
235
237 integer, public, parameter :: &
238 par_independent = 1, &
239 par_dependent = 2
240
242 interface states_elec_get_state
245 end interface states_elec_get_state
246
250 end interface states_elec_set_state
251
252 interface states_elec_get_points
255
258 end interface
260contains
261
262 ! TODO(Alex): Issue #826. Rename to something like "states_elec_default_wfs_type", or remove
263 subroutine states_elec_null(st)
264 class(states_elec_t), intent(inout) :: st
268 st%wfs_type = type_float ! By default, calculations use real wavefunctions
269
270 st%packed = .false.
271
273 end subroutine states_elec_null
277 subroutine states_elec_init(st, namespace, space, valence_charge, kpoints, calc_mode_id)
278 type(states_elec_t), target, intent(inout) :: st
279 type(namespace_t), intent(in) :: namespace
280 type(electron_space_t), intent(in) :: space
281 real(real64), intent(in) :: valence_charge
282 type(kpoints_t), intent(in) :: kpoints
283 integer, optional, intent(in) :: calc_mode_id
285 real(real64) :: excess_charge, nempty_percent
286 integer :: nempty, ntot, default
287 integer :: nempty_conv, nempty_conv_default
288 logical :: force, adapt_for_chebyshev
289 real(real64), parameter :: tol = 1e-13_real64
290 integer :: es
291 integer, parameter :: rs_chebyshev = 12
292 integer(int64), parameter :: chebyshev_compatible_modes(3) = &
293 [option__calculationmode__gs, option__calculationmode__go, option__calculationmode__unocc]
294
295 push_sub_with_profile(states_elec_init)
296
297 st%fromScratch = .true. ! this will be reset if restart_read is called
300 ! We get the spin dimension from the electronic space
301 ! TODO: Remove spin space information from states_elec_dim
302 st%d%ispin = space%ispin
304 ! Use of spinors requires complex wavefunctions.
305 if (st%d%ispin == spinors) call states_set_complex(st)
306
307 if (st%d%ispin /= unpolarized .and. kpoints%use_time_reversal) then
308 message(1) = "Time reversal symmetry is only implemented for unpolarized spins."
309 message(2) = "Use KPointsUseTimeReversal = no."
310 call messages_fatal(2, namespace=namespace)
311 end if
312
314 !%Variable ExcessCharge
315 !%Type float
316 !%Default 0.0
317 !%Section States
318 !%Description
319 !% The net charge of the system. A negative value means that we are adding
320 !% electrons, while a positive value means we are taking electrons
321 !% from the system.
322 !%End
323 call parse_variable(namespace, 'ExcessCharge', m_zero, excess_charge)
325 !%Variable TotalStates
326 !%Type integer
327 !%Default 0
328 !%Section States
329 !%Description
330 !% This variable sets the total number of states that Octopus will
331 !% use. This is normally not necessary since by default Octopus
332 !% sets the number of states to the minimum necessary to hold the
333 !% electrons present in the system. (This default behavior is
334 !% obtained by setting <tt>TotalStates</tt> to 0).
335 !%
336 !% If you want to add some unoccupied states, probably it is more convenient to use the variable
337 !% <tt>ExtraStates</tt>.
338 !%End
339 call parse_variable(namespace, 'TotalStates', 0, ntot)
340 if (ntot < 0) then
341 write(message(1), '(a,i5,a)') "Input: '", ntot, "' is not a valid value for TotalStates."
342 call messages_fatal(1, namespace=namespace)
343 end if
344
345 !%Variable ExtraStates
346 !%Type integer
347 !%Default 0
348 !%Section States
349 !%Description
350 !% The number of states is in principle calculated considering the minimum
351 !% numbers of states necessary to hold the electrons present in the system.
352 !% The number of electrons is
353 !% in turn calculated considering the nature of the species supplied in the
354 !% <tt>Species</tt> block, and the value of the <tt>ExcessCharge</tt> variable.
355 !% However, one may command <tt>Octopus</tt> to use more states, which is necessary if one wants to
356 !% use fractional occupational numbers, either fixed from the beginning through
357 !% the <tt>Occupations</tt> block or by prescribing
358 !% an electronic temperature with <tt>Smearing</tt>, or in order to calculate
359 !% excited states (including with <tt>CalculationMode = unocc</tt>).
360 !%End
361 call parse_variable(namespace, 'ExtraStates', 0, nempty)
362 if (nempty < 0) then
363 write(message(1), '(a,i5,a)') "Input: '", nempty, "' is not a valid value for ExtraStates."
364 message(2) = '(0 <= ExtraStates)'
365 call messages_fatal(2, namespace=namespace)
366 end if
367
368 if (ntot > 0 .and. nempty > 0) then
369 message(1) = 'You cannot set TotalStates and ExtraStates at the same time.'
370 call messages_fatal(1, namespace=namespace)
371 end if
373 !%Variable ExtraStatesInPercent
374 !%Type float
375 !%Default 0
376 !%Section States
377 !%Description
378 !% This variable allows to set the number of extra/empty states as percentage of the
379 !% used occupied states. For example, a value 35 for ExtraStatesInPercent would amount
380 !% to ceiling(35/100 * nstates) extra states, where nstates denotes the amount of occupied
381 !% states Octopus is using for the system at hand.
382 !%End
383 call parse_variable(namespace, 'ExtraStatesInPercent', m_zero, nempty_percent)
384 if (nempty_percent < 0) then
385 write(message(1), '(a,f8.6,a)') "Input: '", nempty_percent, &
386 "' should be a percentage value x (where x is parts in hundred) larger or equal 0"
387 call messages_fatal(1, namespace=namespace)
388 end if
389
390 if (nempty > 0 .and. nempty_percent > 0) then
391 message(1) = 'You cannot set ExtraStates and ExtraStatesInPercent at the same time.'
392 call messages_fatal(1, namespace=namespace)
393 end if
394
395 ! Only use extra states for Chebyshev filtering in ground state or geometry optimization runs
396 adapt_for_chebyshev = .false.
397 if (present(calc_mode_id)) then
398 if (any(calc_mode_id == chebyshev_compatible_modes)) then
399 ! this variable is documented in electrons/eigensolver.F90
400 call parse_variable(namespace, 'Eigensolver', rs_chebyshev, es)
401 if (es == rs_chebyshev) then
402 ! reuse this condition below
403 adapt_for_chebyshev = .true.
404 end if
405 end if
406 end if
407 if (adapt_for_chebyshev) then
408 if (nempty == 0 .and. is_close(nempty_percent, m_zero)) then
409 nempty_percent = 15
410 message(1) = 'Chebyshev filtering eigensolver detected. Setting ExtraStatesInPercent = 15'
411 call messages_info(1, namespace=namespace)
412 end if
413 end if
414
415 ! For non-periodic systems this should just return the Gamma point
416 call states_elec_choose_kpoints(st, kpoints, namespace)
417
418 st%val_charge = valence_charge
419
420 st%qtot = -(st%val_charge + excess_charge)
421
422 if (st%qtot < -m_epsilon) then
423 write(message(1),'(a,f12.6,a)') 'Total charge = ', st%qtot, ' < 0'
424 message(2) = 'Check Species and ExcessCharge.'
425 call messages_fatal(2, only_root_writes = .true., namespace=namespace)
426 end if
427
428 select case (st%d%ispin)
429 case (unpolarized)
430 st%d%dim = 1
431 st%nst = nint(st%qtot/2)
432 if (st%nst*2 - st%qtot < -tol) st%nst = st%nst + 1
433 st%d%nspin = 1
434 st%d%spin_channels = 1
435 case (spin_polarized)
436 st%d%dim = 1
437 st%nst = nint(st%qtot/2)
438 if (st%nst*2 - st%qtot < -tol) st%nst = st%nst + 1
439 st%d%nspin = 2
440 st%d%spin_channels = 2
441 case (spinors)
442 st%d%dim = 2
443 st%nst = nint(st%qtot)
444 if (st%nst - st%qtot < -tol) st%nst = st%nst + 1
445 st%d%nspin = 4
446 st%d%spin_channels = 2
447 end select
448
449 if (nempty_percent > 0) then
450 nempty = ceiling(nempty_percent * st%nst / 100)
451 end if
452
453 !%Variable ExtraStatesToConverge
454 !%Type integer
455 !%Default <tt>ExtraStates</tt> (Default 0)
456 !%Section States
457 !%Description
458 !% For <tt>gs</tt> and <tt>unocc</tt> calculations.
459 !% (For the <tt>gs</tt> calculation one needs to set <tt>ConvEigenError=yes</tt>)
460 !% Specifies the number of extra states that will be considered for reaching the convergence.
461 !% The calculation will consider the number off occupied states plus
462 !% <tt>ExtraStatesToConverge</tt> for the convergence criteria.
463 !% By default, all extra states need to be converged (For <tt>gs</tt> calculations only with <tt>ConvEigenError=yes</tt>).
464 !% Thus, together with <tt>ExtraStates</tt>, one can have some more states which will not be
465 !% considered for the convergence criteria, thus making the convergence of the
466 !% unocc calculation faster.
467 !%
468 !% If chebyshev filtering is used, the default is <tt>min(0.8*ExtraStates, ExtraStates - 1)</tt>.
469 !%End
470 nempty_conv_default = nempty
471 if (adapt_for_chebyshev) then
472 ! use 80% of the extra states, but at least do not converge one extra state
473 nempty_conv_default = min(int(0.8*nempty), nempty - 1)
474 end if
475 call parse_variable(namespace, 'ExtraStatesToConverge', nempty_conv_default, nempty_conv)
476 if (nempty_conv < 0) then
477 write(message(1), '(a,i5,a)') "Input: '", nempty_conv, "' is not a valid value for ExtraStatesToConverge."
478 message(2) = '(0 <= ExtraStatesToConverge)'
479 call messages_fatal(2, namespace=namespace)
480 end if
481
482 if (nempty_conv > nempty) then
483 nempty_conv = nempty
484 message(1) = 'You cannot set ExtraStatesToConverge to a higher value than ExtraStates.'
485 message(2) = 'Capping ExtraStatesToConverge to ExtraStates.'
486 call messages_warning(2, namespace=namespace)
487 endif
488
489 if (ntot > 0) then
490 if (ntot < st%nst) then
491 message(1) = 'TotalStates is smaller than the number of states required by the system.'
492 call messages_fatal(1, namespace=namespace)
493 end if
494
495 st%nst = ntot
496 end if
497
498 st%nst_conv = st%nst + nempty_conv
499 st%nst = st%nst + nempty
500 if (st%nst == 0) then
501 message(1) = "Cannot run with number of states = zero."
502 call messages_fatal(1, namespace=namespace)
503 end if
504
505 !%Variable StatesBlockSize
506 !%Type integer
507 !%Section Execution::Optimization
508 !%Description
509 !% Some routines work over blocks of eigenfunctions, which
510 !% generally improves performance at the expense of increased
511 !% memory consumption. This variable selects the size of the
512 !% blocks to be used. If GPUs are used, the default is the
513 !% warp size (32 for NVIDIA, 32 or 64 for AMD);
514 !% otherwise it is 4.
515 !%End
516
517 if (accel_is_enabled()) then
518 ! use the warp size (usually 32, for some AMD GPUs it is 64)
519 default = accel%warp_size
520 else
521 default = 4
522 end if
523
524 if (default > pad_pow2(st%nst)) default = pad_pow2(st%nst)
525
526 assert(default > 0)
527
528 call parse_variable(namespace, 'StatesBlockSize', default, st%block_size)
529 if (st%block_size < 1) then
530 call messages_write("The variable 'StatesBlockSize' must be greater than 0.")
531 call messages_fatal(namespace=namespace)
532 end if
533
534 st%block_size = min(st%block_size, st%nst)
535 conf%target_states_block_size = st%block_size
536
537 safe_allocate(st%eigenval(1:st%nst, 1:st%nik))
538 st%eigenval = huge(st%eigenval)
539
540 ! Periodic systems require complex wavefunctions
541 ! but not if it is Gamma-point only
542 if (.not. kpoints%gamma_only()) then
543 call states_set_complex(st)
544 end if
545
546 !%Variable OnlyUserDefinedInitialStates
547 !%Type logical
548 !%Default no
549 !%Section States
550 !%Description
551 !% If true, then only user-defined states from the block <tt>UserDefinedStates</tt>
552 !% will be used as initial states for a time-propagation. No attempt is made
553 !% to load ground-state orbitals from a previous ground-state run.
554 !%End
555 call parse_variable(namespace, 'OnlyUserDefinedInitialStates', .false., st%only_userdef_istates)
556
557 ! we now allocate some arrays
558 safe_allocate(st%occ (1:st%nst, 1:st%nik))
559 st%occ = m_zero
560 ! allocate space for formula strings that define user-defined states
561 if (parse_is_defined(namespace, 'UserDefinedStates') .or. parse_is_defined(namespace, 'OCTInitialUserdefined') &
562 .or. parse_is_defined(namespace, 'OCTTargetUserdefined')) then
563 safe_allocate(st%user_def_states(1:st%d%dim, 1:st%nst, 1:st%nik))
564 ! initially we mark all 'formulas' as undefined
565 st%user_def_states(1:st%d%dim, 1:st%nst, 1:st%nik) = 'undefined'
566 end if
567
568 if (st%d%ispin == spinors) then
569 safe_allocate(st%spin(1:3, 1:st%nst, 1:st%nik))
570 end if
571
572 !%Variable StatesRandomization
573 !%Type integer
574 !%Default par_independent
575 !%Section States
576 !%Description
577 !% The randomization of states can be done in two ways:
578 !% i) a parallelisation independent way (default), where the random states are identical,
579 !% irrespectively of the number of tasks and
580 !% ii) a parallelisation dependent way, which can prevent linear dependency
581 !% to occur for large systems.
582 !%Option par_independent 1
583 !% Parallelisation-independent randomization of states.
584 !%Option par_dependent 2
585 !% The randomization depends on the number of taks used in the calculation.
586 !%End
587 call parse_variable(namespace, 'StatesRandomization', par_independent, st%randomization)
588
589
590 call states_elec_read_initial_occs(st, namespace, excess_charge, kpoints)
591 call states_elec_read_initial_spins(st, namespace)
592
593 ! This test can only be done here, as smear_init is called inside states_elec_read_initial_occs, and
594 ! only there smear%photodop is set.
595
596 if (st%smear%photodop) then
597 if (nempty == 0) then
598 write(message(1), '(a,i5,a)') "PhotoDoping requires to specify ExtraStates."
599 message(2) = '(0 == ExtraStates)'
600 call messages_fatal(2, namespace=namespace)
601 end if
602 end if
603
604 st%st_start = 1
605 st%st_end = st%nst
606 st%lnst = st%nst
607 safe_allocate(st%node(1:st%nst))
608 st%node(1:st%nst) = 0
609
610 call mpi_grp_init(st%mpi_grp, mpi_comm_undefined)
611 st%parallel_in_states = .false.
612
613 call distributed_nullify(st%d%kpt, st%nik)
614
615 call modelmb_particles_init(st%modelmbparticles, namespace, space, st%nst)
616
617 !%Variable SymmetrizeDensity
618 !%Type logical
619 !%Default no
620 !%Section States
621 !%Description
622 !% When enabled the density is symmetrized. Currently, this can
623 !% only be done for periodic systems. (Experimental.)
624 !%End
625 call parse_variable(namespace, 'SymmetrizeDensity', kpoints%use_symmetries, st%symmetrize_density)
626 call messages_print_var_value('SymmetrizeDensity', st%symmetrize_density, namespace=namespace)
627
628 !%Variable ForceComplex
629 !%Type logical
630 !%Default no
631 !%Section Execution::Debug
632 !%Description
633 !% Normally <tt>Octopus</tt> determines automatically the type necessary
634 !% for the wavefunctions. When set to yes this variable will
635 !% force the use of complex wavefunctions.
636 !%
637 !% Warning: This variable is designed for testing and
638 !% benchmarking and normal users need not use it.
639 !%End
640 call parse_variable(namespace, 'ForceComplex', .false., force)
641
642 if (force) call states_set_complex(st)
643
644 st%packed = .false.
645
646 pop_sub_with_profile(states_elec_init)
647 end subroutine states_elec_init
648
649 ! ---------------------------------------------------------
652 !
653 subroutine states_elec_look(restart, nik, dim, nst, ierr)
654 class(restart_t), intent(in) :: restart
655 integer, intent(out) :: nik
656 integer, intent(out) :: dim
657 integer, intent(out) :: nst
658 integer, intent(out) :: ierr
659
660 character(len=256) :: lines(3)
661 character(len=20) :: char
662 integer :: iunit
663
664 push_sub(states_elec_look)
665
666 ierr = 0
667
668 iunit = restart%open('states')
669 call restart%read(iunit, lines, 3, ierr)
670 if (ierr == 0) then
671 read(lines(1), *) char, nst
672 read(lines(2), *) char, dim
673 read(lines(3), *) char, nik
674 end if
675 call restart%close(iunit)
676
677 pop_sub(states_elec_look)
678 end subroutine states_elec_look
679
680 ! ---------------------------------------------------------
689 !
690 subroutine states_elec_read_initial_occs(st, namespace, excess_charge, kpoints)
691 type(states_elec_t), intent(inout) :: st
692 type(namespace_t), intent(in) :: namespace
693 real(real64), intent(in) :: excess_charge
694 type(kpoints_t), intent(in) :: kpoints
695
696 integer :: ik, ist, ispin, nspin, ncols, nrows, el_per_state, icol, start_pos, spin_n
697 type(block_t) :: blk
698 real(real64) :: rr, charge
699 logical :: integral_occs, unoccupied_states
700 real(real64), allocatable :: read_occs(:, :)
701 real(real64) :: charge_in_block
702
704
705 !%Variable RestartFixedOccupations
706 !%Type logical
707 !%Default yes
708 !%Section States
709 !%Description
710 !% Setting this variable will make the restart proceed as
711 !% if the occupations from the previous calculation had been set via the <tt>Occupations</tt> block,
712 !% <i>i.e.</i> fixed. Otherwise, occupations will be determined by smearing.
713 !%End
714 call parse_variable(namespace, 'RestartFixedOccupations', .true., st%restart_fixed_occ)
715 ! we will turn on st%fixed_occ if restart_read is ever called
716
717 !%Variable Occupations
718 !%Type block
719 !%Section States
720 !%Description
721 !% The occupation numbers of the orbitals can be fixed through the use of this
722 !% variable. For example:
723 !%
724 !% <tt>%Occupations
725 !% <br>&nbsp;&nbsp;2 | 2 | 2 | 2 | 2
726 !% <br>%</tt>
727 !%
728 !% would fix the occupations of the five states to 2. There can be
729 !% at most as many columns as states in the calculation. If there are fewer columns
730 !% than states, then the code will assume that the user is indicating the occupations
731 !% of the uppermost states where all lower states have full occupation (i.e. 2 for spin-unpolarized
732 !% calculations, 1 otherwise) and all higher states have zero occupation. The first column
733 !% will be taken to refer to the lowest state such that the occupations would be consistent
734 !% with the correct total charge. For example, if there are 8 electrons and 10 states (from
735 !% <tt>ExtraStates = 6</tt>), then an abbreviated specification
736 !%
737 !% <tt>%Occupations
738 !% <br>&nbsp;&nbsp;1 | 0 | 1
739 !% <br>%</tt>
740 !%
741 !% would be equivalent to a full specification
742 !%
743 !% <tt>%Occupations
744 !% <br>&nbsp;&nbsp;2 | 2 | 2 | 1 | 0 | 1 | 0 | 0 | 0 | 0
745 !% <br>%</tt>
746 !%
747 !% This is an example of use for constrained density-functional theory,
748 !% crudely emulating a HOMO->LUMO+1 optical excitation.
749 !% The number of rows should be equal
750 !% to the number of k-points times the number of spins. For example, for a finite system
751 !% with <tt>SpinComponents == spin_polarized</tt>,
752 !% this block should contain two lines, one for each spin channel.
753 !% All rows must have the same number of columns.
754 !%
755 !% The <tt>Occupations</tt> block is useful for the ground state of highly symmetric
756 !% small systems (like an open-shell atom), to fix the occupation numbers
757 !% of degenerate states in order to help <tt>octopus</tt> to converge. This is to
758 !% be used in conjuction with <tt>ExtraStates</tt>. For example, to calculate the
759 !% carbon atom, one would do:
760 !%
761 !% <tt>ExtraStates = 2
762 !% <br>%Occupations
763 !% <br>&nbsp;&nbsp;2 | 2/3 | 2/3 | 2/3
764 !% <br>%</tt>
765 !%
766 !% If you want the calculation to be spin-polarized (which makes more sense), you could do:
767 !%
768 !% <tt>ExtraStates = 2
769 !% <br>%Occupations
770 !% <br>&nbsp;&nbsp; 2/3 | 2/3 | 2/3
771 !% <br>&nbsp;&nbsp; 0 | 0 | 0
772 !% <br>%</tt>
773 !%
774 !% Note that in this case the first state is absent, the code will calculate four states
775 !% (two because there are four electrons, plus two because <tt>ExtraStates</tt> = 2), and since
776 !% it finds only three columns, it will occupy the first state with one electron for each
777 !% of the spin options.
778 !%
779 !% If the sum of occupations is not equal to the total charge set by <tt>ExcessCharge</tt>,
780 !% an error message is printed.
781 !% If <tt>FromScratch = no</tt> and <tt>RestartFixedOccupations = yes</tt>,
782 !% this block will be ignored.
783 !%End
784
785 integral_occs = .true.
786
787 occ_fix: if (parse_block(namespace, 'Occupations', blk) == 0) then
788 ! read in occupations
789 st%fixed_occ = .true.
790
791 ncols = parse_block_cols(blk, 0)
792 if (ncols > st%nst) then
793 call messages_input_error(namespace, "Occupations", "Too many columns in block Occupations.")
794 end if
795
796 nrows = parse_block_n(blk)
797 if (nrows /= st%nik) then
798 call messages_input_error(namespace, "Occupations", "Wrong number of rows in block Occupations.")
799 end if
800
801 do ik = 1, st%nik - 1
802 if (parse_block_cols(blk, ik) /= ncols) then
803 call messages_input_error(namespace, "Occupations", &
804 "All rows in block Occupations must have the same number of columns.")
805 end if
806 end do
807
808 ! Now we fill all the "missing" states with the maximum occupation.
809 if (st%d%ispin == unpolarized) then
810 el_per_state = 2
811 else
812 el_per_state = 1
813 end if
814
815 safe_allocate(read_occs(1:ncols, 1:st%nik))
816
817 charge_in_block = m_zero
818 do ik = 1, st%nik
819 do icol = 1, ncols
820 call parse_block_float(blk, ik - 1, icol - 1, read_occs(icol, ik))
821 charge_in_block = charge_in_block + read_occs(icol, ik) * st%kweights(ik)
822 end do
823 end do
824
825 spin_n = 2
826 select case (st%d%ispin)
827 case (unpolarized)
828 spin_n = 2
829 case (spin_polarized)
830 spin_n = 2
831 case (spinors)
832 spin_n = 1
833 end select
834
835 start_pos = nint((st%qtot - charge_in_block)/spin_n)
836
837 if (start_pos + ncols > st%nst) then
838 message(1) = "To balance charge, the first column in block Occupations is taken to refer to state"
839 write(message(2),'(a,i6,a)') "number ", start_pos, " but there are too many columns for the number of states."
840 write(message(3),'(a,i6,a)') "Solution: set ExtraStates = ", start_pos + ncols - st%nst
841 call messages_fatal(3, namespace=namespace)
842 end if
843
844 do ik = 1, st%nik
845 do ist = 1, start_pos
846 st%occ(ist, ik) = el_per_state
847 end do
848 end do
849
850 do ik = 1, st%nik
851 do ist = start_pos + 1, start_pos + ncols
852 st%occ(ist, ik) = read_occs(ist - start_pos, ik)
853 integral_occs = integral_occs .and. &
854 abs((st%occ(ist, ik) - el_per_state) * st%occ(ist, ik)) <= m_epsilon
855 end do
856 end do
857
858 do ik = 1, st%nik
859 do ist = start_pos + ncols + 1, st%nst
860 st%occ(ist, ik) = m_zero
861 end do
862 end do
863
864 call parse_block_end(blk)
865
866 safe_deallocate_a(read_occs)
867
868 else
869 st%fixed_occ = .false.
870 integral_occs = .false.
871
872 ! first guess for occupation...paramagnetic configuration
873 rr = m_one
874 if (st%d%ispin == unpolarized) rr = m_two
875
876 st%occ = m_zero
877 st%qtot = -(st%val_charge + excess_charge)
878
879 nspin = 1
880 if (st%d%nspin == 2) nspin = 2
881
882 do ik = 1, st%nik, nspin
883 charge = m_zero
884 do ist = 1, st%nst
885 do ispin = ik, ik + nspin - 1
886 st%occ(ist, ispin) = min(rr, -(st%val_charge + excess_charge) - charge)
887 charge = charge + st%occ(ist, ispin)
888 end do
889 end do
890 end do
891
892 end if occ_fix
893
894 !%Variable RestartReorderOccs
895 !%Type logical
896 !%Default no
897 !%Section States
898 !%Description
899 !% Consider doing a ground-state calculation, and then restarting with new occupations set
900 !% with the <tt>Occupations</tt> block, in an attempt to populate the orbitals of the original
901 !% calculation. However, the eigenvalues may reorder as the density changes, in which case the
902 !% occupations will now be referring to different orbitals. Setting this variable to yes will
903 !% try to solve this issue when the restart data is being read, by reordering the occupations
904 !% according to the order of the expectation values of the restart wavefunctions.
905 !%End
906 if (st%fixed_occ) then
907 call parse_variable(namespace, 'RestartReorderOccs', .false., st%restart_reorder_occs)
908 else
909 st%restart_reorder_occs = .false.
910 end if
911
912 call smear_init(st%smear, namespace, st%d%ispin, st%fixed_occ, integral_occs, kpoints)
913
914 unoccupied_states = (st%d%ispin /= spinors .and. st%nst*2 > st%qtot) .or. (st%d%ispin == spinors .and. st%nst > st%qtot)
915
916 if (.not. smear_is_semiconducting(st%smear) .and. .not. st%smear%method == smear_fixed_occ) then
917 if (.not. unoccupied_states) then
918 call messages_write('Smearing needs unoccupied states (via ExtraStates or TotalStates) to be useful.')
919 call messages_warning(namespace=namespace)
920 end if
921 end if
922
923 ! sanity check
924 charge = m_zero
925 do ist = 1, st%nst
926 charge = charge + sum(st%occ(ist, 1:st%nik) * st%kweights(1:st%nik))
927 end do
928 if (abs(charge - st%qtot) > 1e-6_real64) then
929 message(1) = "Initial occupations do not integrate to total charge."
930 write(message(2), '(6x,f12.6,a,f12.6)') charge, ' != ', st%qtot
931 call messages_fatal(2, only_root_writes = .true., namespace=namespace)
932 end if
933
934 st%uniform_occ = smear_is_semiconducting(st%smear) .and. .not. unoccupied_states
935
937 end subroutine states_elec_read_initial_occs
938
939
940 ! ---------------------------------------------------------
947 !
948 subroutine states_elec_read_initial_spins(st, namespace)
949 type(states_elec_t), intent(inout) :: st
950 type(namespace_t), intent(in) :: namespace
951
952 integer :: i, j, nrows
953 type(block_t) :: blk
954
956
957 st%fixed_spins = .false.
958 if (st%d%ispin /= spinors) then
960 return
961 end if
962
963 !%Variable InitialSpins
964 !%Type block
965 !%Section States
966 !%Description
967 !% The spin character of the initial random guesses for the spinors can
968 !% be fixed by making use of this block. Note that this will not "fix" the
969 !% the spins during the calculation (this cannot be done in spinors mode, in
970 !% being able to change the spins is why the spinors mode exists in the first
971 !% place).
972 !%
973 !% This block is meaningless and ignored if the run is not in spinors mode
974 !% (<tt>SpinComponents = spinors</tt>).
975 !%
976 !% The structure of the block is very simple: each column contains the desired
977 !% <math>\left< S_x \right>, \left< S_y \right>, \left< S_z \right> </math> for each spinor.
978 !% If the calculation is for a periodic system
979 !% and there is more than one <i>k</i>-point, the spins of all the <i>k</i>-points are
980 !% the same.
981 !%
982 !% For example, if we have two spinors, and we want one in the <math>S_x</math> "down" state,
983 !% and another one in the <math>S_x</math> "up" state:
984 !%
985 !% <tt>%InitialSpins
986 !% <br>&nbsp;&nbsp;&nbsp; 0.5 | 0.0 | 0.0
987 !% <br>&nbsp;&nbsp; -0.5 | 0.0 | 0.0
988 !% <br>%</tt>
989 !%
990 !% WARNING: if the calculation is for a system described by pseudopotentials (as
991 !% opposed to user-defined potentials or model systems), this option is
992 !% meaningless since the random spinors are overwritten by the atomic orbitals.
993 !%
994 !% This constraint must be fulfilled:
995 !% <br><math> \left< S_x \right>^2 + \left< S_y \right>^2 + \left< S_z \right>^2 = \frac{1}{4} </math>
996 !%End
997 spin_fix: if (parse_block(namespace, 'InitialSpins', blk) == 0) then
998 nrows = parse_block_n(blk)
999 if (nrows < st%nst) then
1000 message(1) = "Please specify one row for each state in InitialSpins, also for extra states."
1001 call messages_fatal(1, namespace=namespace)
1002 end if
1003 do i = 1, st%nst
1004 do j = 1, 3
1005 call parse_block_float(blk, i-1, j-1, st%spin(j, i, 1))
1006 end do
1007 if (abs(sum(st%spin(1:3, i, 1)**2) - m_fourth) > 1.0e-6_real64) call messages_input_error(namespace, 'InitialSpins')
1008 end do
1009 call parse_block_end(blk)
1010 st%fixed_spins = .true.
1011 do i = 2, st%nik
1012 st%spin(:, :, i) = st%spin(:, :, 1)
1013 end do
1014 end if spin_fix
1015
1017 end subroutine states_elec_read_initial_spins
1018
1019
1020 ! ---------------------------------------------------------
1022 !
1023 subroutine states_elec_allocate_wfns(st, mesh, wfs_type, skip, packed)
1024 type(states_elec_t), intent(inout) :: st
1025 class(mesh_t), intent(in) :: mesh
1026 type(type_t), optional, intent(in) :: wfs_type
1027 logical, optional, intent(in) :: skip(:)
1028 logical, optional, intent(in) :: packed
1029
1031
1032 if (present(wfs_type)) then
1033 assert(wfs_type == type_float .or. wfs_type == type_cmplx)
1034 st%wfs_type = wfs_type
1035 end if
1036
1037 call states_elec_init_block(st, mesh, skip = skip, packed=packed)
1038 call states_elec_set_zero(st)
1039
1041 end subroutine states_elec_allocate_wfns
1042
1043 !---------------------------------------------------------------------
1060 subroutine states_elec_init_block(st, mesh, verbose, skip, packed)
1061 type(states_elec_t), intent(inout) :: st
1062 type(mesh_t), intent(in) :: mesh
1063 logical, optional, intent(in) :: verbose
1064 logical, optional, intent(in) :: skip(:)
1065 logical, optional, intent(in) :: packed
1066
1067 integer :: ib, iqn, ist, istmin, istmax
1068 logical :: same_node, verbose_, packed_
1069 integer, allocatable :: bstart(:), bend(:)
1070
1071 push_sub(states_elec_init_block)
1072
1073 safe_allocate(bstart(1:st%nst))
1074 safe_allocate(bend(1:st%nst))
1075 safe_allocate(st%group%iblock(1:st%nst))
1076
1077 st%group%iblock = 0
1078
1079 verbose_ = optional_default(verbose, .true.)
1080 packed_ = optional_default(packed, .false.)
1081
1082 !In case we have a list of state to skip, we do not allocate them
1083 istmin = 1
1084 if (present(skip)) then
1085 do ist = 1, st%nst
1086 if (.not. skip(ist)) then
1087 istmin = ist
1088 exit
1089 end if
1090 end do
1091 end if
1092
1093 istmax = st%nst
1094 if (present(skip)) then
1095 do ist = st%nst, istmin, -1
1096 if (.not. skip(ist)) then
1097 istmax = ist
1098 exit
1099 end if
1100 end do
1101 end if
1102
1103 if (present(skip) .and. verbose_) then
1104 call messages_write('Info: Allocating states from ')
1105 call messages_write(istmin, fmt = 'i8')
1106 call messages_write(' to ')
1107 call messages_write(istmax, fmt = 'i8')
1108 call messages_info()
1109 end if
1110
1111 ! count and assign blocks
1112 ib = 0
1113 st%group%nblocks = 0
1114 bstart(1) = istmin
1115 do ist = istmin, istmax
1116 ib = ib + 1
1117
1118 st%group%iblock(ist) = st%group%nblocks + 1
1119
1120 same_node = .true.
1121 if (st%parallel_in_states .and. ist /= istmax) then
1122 ! We have to avoid that states that are in different nodes end
1123 ! up in the same block
1124 same_node = (st%node(ist + 1) == st%node(ist))
1125 end if
1126
1127 if (ib == st%block_size .or. ist == istmax .or. .not. same_node) then
1128 ib = 0
1129 st%group%nblocks = st%group%nblocks + 1
1130 bend(st%group%nblocks) = ist
1131 if (ist /= istmax) bstart(st%group%nblocks + 1) = ist + 1
1132 end if
1133 end do
1134
1135 safe_allocate(st%group%psib(1:st%group%nblocks, st%d%kpt%start:st%d%kpt%end))
1136
1137 safe_allocate(st%group%block_is_local(1:st%group%nblocks, 1:st%nik))
1138 st%group%block_is_local = .false.
1139 st%group%block_start = -1
1140 st%group%block_end = -2 ! this will make that loops block_start:block_end do not run if not initialized
1141
1142 do ib = 1, st%group%nblocks
1143 if (bstart(ib) >= st%st_start .and. bend(ib) <= st%st_end) then
1144 if (st%group%block_start == -1) st%group%block_start = ib
1145 st%group%block_end = ib
1146 do iqn = st%d%kpt%start, st%d%kpt%end
1147 st%group%block_is_local(ib, iqn) = .true.
1148
1149 if (states_are_real(st)) then
1150 call dwfs_elec_init(st%group%psib(ib, iqn), st%d%dim, bstart(ib), bend(ib), mesh%np_part, iqn, &
1151 special=.true., packed=packed_)
1152 else
1153 call zwfs_elec_init(st%group%psib(ib, iqn), st%d%dim, bstart(ib), bend(ib), mesh%np_part, iqn, &
1154 special=.true., packed=packed_)
1155 end if
1156
1157 end do
1158 end if
1159 end do
1160
1161 safe_allocate(st%group%block_range(1:st%group%nblocks, 1:2))
1162 safe_allocate(st%group%block_size(1:st%group%nblocks))
1163
1164 st%group%block_range(1:st%group%nblocks, 1) = bstart(1:st%group%nblocks)
1165 st%group%block_range(1:st%group%nblocks, 2) = bend(1:st%group%nblocks)
1166 st%group%block_size(1:st%group%nblocks) = bend(1:st%group%nblocks) - bstart(1:st%group%nblocks) + 1
1167
1168 st%group%block_initialized = .true.
1169
1170 safe_allocate(st%group%block_node(1:st%group%nblocks, 1:st%nik))
1171 st%group%block_node = 0
1172
1173 assert(allocated(st%node))
1174 assert(all(st%node >= 0) .and. all(st%node < st%mpi_grp%size))
1175
1176 do iqn = st%d%kpt%start, st%d%kpt%end
1177 do ib = st%group%block_start, st%group%block_end
1178 st%group%block_node(ib, iqn) = st%st_kpt_mpi_grp%rank
1179 end do
1180 end do
1181
1182 call comm_allreduce(st%st_kpt_mpi_grp, st%group%block_node)
1183
1184 if (verbose_) then
1185 call messages_write('Info: Blocks of states')
1186 call messages_info()
1187 do ib = 1, st%group%nblocks
1188 call messages_write(' Block ')
1189 call messages_write(ib, fmt = 'i8')
1190 call messages_write(' contains ')
1191 call messages_write(st%group%block_size(ib), fmt = 'i8')
1192 call messages_write(' states')
1193 if (st%group%block_size(ib) > 0) then
1194 call messages_write(':')
1195 call messages_write(st%group%block_range(ib, 1), fmt = 'i8')
1196 call messages_write(' - ')
1197 call messages_write(st%group%block_range(ib, 2), fmt = 'i8')
1198 end if
1199 call messages_info()
1200 end do
1201 end if
1202
1203!!$!!!!DEBUG
1204!!$ ! some debug output that I will keep here for the moment
1205!!$ if (mpi_world%is_root()) then
1206!!$ print*, "NST ", st%nst
1207!!$ print*, "BLOCKSIZE ", st%block_size
1208!!$ print*, "NBLOCKS ", st%group%nblocks
1209!!$
1210!!$ print*, "==============="
1211!!$ do ist = 1, st%nst
1212!!$ print*, st%node(ist), ist, st%group%iblock(ist)
1213!!$ end do
1214!!$ print*, "==============="
1215!!$
1216!!$ do ib = 1, st%group%nblocks
1217!!$ print*, ib, bstart(ib), bend(ib)
1218!!$ end do
1219!!$
1220!!$ end if
1221!!$!!!!ENDOFDEBUG
1222
1223 safe_deallocate_a(bstart)
1224 safe_deallocate_a(bend)
1225 pop_sub(states_elec_init_block)
1226 end subroutine states_elec_init_block
1227
1228
1229 ! ---------------------------------------------------------
1231 subroutine states_elec_deallocate_wfns(st)
1232 type(states_elec_t), intent(inout) :: st
1233
1235
1236 call states_elec_group_end(st%group, st%d)
1237
1239 end subroutine states_elec_deallocate_wfns
1240
1241
1242 ! ---------------------------------------------------------
1243 subroutine states_elec_densities_init(st, gr)
1244 type(states_elec_t), target, intent(inout) :: st
1245 type(grid_t), intent(in) :: gr
1246
1247 real(real64) :: fsize
1248
1250
1251 safe_allocate(st%rho(1:gr%np_part, 1:st%d%nspin))
1252 st%rho = m_zero
1253
1254 fsize = gr%np_part*8.0_real64*st%block_size
1255
1256 call messages_write('Info: states-block size = ')
1257 call messages_write(fsize, fmt = '(f10.1)', align_left = .true., units = unit_megabytes, print_units = .true.)
1258 call messages_info()
1259
1261 end subroutine states_elec_densities_init
1262
1263 !---------------------------------------------------------------------
1264 subroutine states_elec_allocate_current(st, space, mesh)
1265 type(states_elec_t), intent(inout) :: st
1266 class(space_t), intent(in) :: space
1267 class(mesh_t), intent(in) :: mesh
1268
1270
1271 if (.not. allocated(st%current)) then
1272 safe_allocate(st%current(1:mesh%np_part, 1:space%dim, 1:st%d%nspin))
1273 st%current = m_zero
1274 end if
1275
1276 if (.not. allocated(st%current_para)) then
1277 safe_allocate(st%current_para(1:mesh%np_part, 1:space%dim, 1:st%d%nspin))
1278 st%current_para = m_zero
1279 end if
1280
1281 if (.not. allocated(st%current_dia)) then
1282 safe_allocate(st%current_dia(1:mesh%np_part, 1:space%dim, 1:st%d%nspin))
1283 st%current_dia= m_zero
1284 end if
1285
1286 if (.not. allocated(st%current_mag)) then
1287 safe_allocate(st%current_mag(1:mesh%np_part, 1:space%dim, 1:st%d%nspin))
1288 st%current_mag= m_zero
1289 end if
1290
1291 if (.not. allocated(st%current_kpt)) then
1292 safe_allocate(st%current_kpt(1:mesh%np,1:space%dim,st%d%kpt%start:st%d%kpt%end))
1293 st%current_kpt = m_zero
1294 end if
1295
1297 end subroutine states_elec_allocate_current
1298
1299 !---------------------------------------------------------------------
1307 subroutine states_elec_exec_init(st, namespace, mc)
1308 type(states_elec_t), intent(inout) :: st
1309 type(namespace_t), intent(in) :: namespace
1310 type(multicomm_t), intent(in) :: mc
1311
1312 integer :: default
1313
1314 push_sub(states_elec_exec_init)
1315
1316 !%Variable StatesPack
1317 !%Type logical
1318 !%Section Execution::Optimization
1319 !%Description
1320 !% When set to yes, states are stored in packed mode, which improves
1321 !% performance considerably. Not all parts of the code will profit from
1322 !% this, but should nevertheless work regardless of how the states are
1323 !% stored.
1324 !%
1325 !% If GPUs are used and this variable is set to yes, Octopus
1326 !% will store the wave-functions in device (GPU) memory. If
1327 !% there is not enough memory to store all the wave-functions,
1328 !% execution will stop with an error.
1329 !%
1330 !% See also the related <tt>HamiltonianApplyPacked</tt> variable.
1331 !%
1332 !% The default is yes.
1333 !%End
1334
1335 call parse_variable(namespace, 'StatesPack', .true., st%pack_states)
1336
1337 call messages_print_var_value('StatesPack', st%pack_states, namespace=namespace)
1339 call messages_obsolete_variable(namespace, 'StatesMirror')
1340
1341 !%Variable StatesOrthogonalization
1342 !%Type integer
1343 !%Section SCF::Eigensolver
1344 !%Description
1345 !% The full orthogonalization method used by some
1346 !% eigensolvers. The default is <tt>cholesky_serial</tt>, except with state
1347 !% parallelization, the default is <tt>cholesky_parallel</tt>.
1348 !%Option cholesky_serial 1
1349 !% Cholesky decomposition implemented using
1350 !% BLAS/LAPACK. Can be used with domain parallelization but not
1351 !% state parallelization.
1352 !%Option cholesky_parallel 2
1353 !% Cholesky decomposition implemented using
1354 !% ScaLAPACK. Compatible with states parallelization.
1355 !%Option cgs 3
1356 !% Classical Gram-Schmidt (CGS) orthogonalization.
1357 !% Can be used with domain parallelization but not state parallelization.
1358 !% The algorithm is defined in Giraud et al., Computers and Mathematics with Applications 50, 1069 (2005).
1359 !%Option mgs 4
1360 !% Modified Gram-Schmidt (MGS) orthogonalization.
1361 !% Can be used with domain parallelization but not state parallelization.
1362 !% The algorithm is defined in Giraud et al., Computers and Mathematics with Applications 50, 1069 (2005).
1363 !%Option drcgs 5
1364 !% Classical Gram-Schmidt orthogonalization with double-step reorthogonalization.
1365 !% Can be used with domain parallelization but not state parallelization.
1366 !% The algorithm is taken from Giraud et al., Computers and Mathematics with Applications 50, 1069 (2005).
1367 !% According to this reference, this is much more precise than CGS or MGS algorithms.
1368 !% The MGS version seems not to improve much the stability and would require more communications over the domains.
1369 !%End
1370
1371 default = option__statesorthogonalization__cholesky_serial
1372#ifdef HAVE_SCALAPACK
1374 default = option__statesorthogonalization__cholesky_parallel
1375 end if
1376#endif
1377
1378 call parse_variable(namespace, 'StatesOrthogonalization', default, st%orth_method)
1379
1380 if (.not. varinfo_valid_option('StatesOrthogonalization', st%orth_method)) then
1381 call messages_input_error(namespace, 'StatesOrthogonalization')
1382 end if
1383 call messages_print_var_option('StatesOrthogonalization', st%orth_method, namespace=namespace)
1384
1385 !%Variable StatesDeviceMemory
1386 !%Type float
1387 !%Section Execution::Optimization
1388 !%Default -512
1389 !%Description
1390 !% This variable selects the amount of device memory that
1391 !% will be used by Octopus to store the states.
1392 !%
1393 !% A positive number smaller than 1 indicates a fraction of the total
1394 !% device memory. A number larger than one indicates an absolute
1395 !% amount of memory in megabytes. A negative number indicates an
1396 !% amount of memory in megabytes that would be subtracted from
1397 !% the total device memory.
1398 !%End
1399 call parse_variable(namespace, 'StatesDeviceMemory', -512.0_real64, st%gpu_states_mem)
1400 call messages_obsolete_variable(namespace, 'StatesCLDeviceMemory', 'StatesDeviceMemory')
1401
1403 end subroutine states_elec_exec_init
1404
1405
1406 ! ---------------------------------------------------------
1408 !
1409 subroutine states_elec_copy(stout, stin, exclude_wfns, exclude_eigenval, special)
1410 type(states_elec_t), target, intent(inout) :: stout
1411 type(states_elec_t), intent(in) :: stin
1412 logical, optional, intent(in) :: exclude_wfns
1413 logical, optional, intent(in) :: exclude_eigenval
1414 logical, optional, intent(in) :: special
1415
1416 logical :: exclude_wfns_
1417
1418 push_sub(states_elec_copy)
1419
1420 exclude_wfns_ = optional_default(exclude_wfns, .false.)
1421
1422 call states_elec_null(stout)
1423
1424 call states_elec_dim_copy(stout%d, stin%d)
1425 safe_allocate_source_a(stout%kweights, stin%kweights)
1426 stout%nik = stin%nik
1427
1428 call modelmb_particles_copy(stout%modelmbparticles, stin%modelmbparticles)
1429
1430 stout%wfs_type = stin%wfs_type
1431 stout%nst = stin%nst
1432 stout%block_size = stin%block_size
1433 stout%orth_method = stin%orth_method
1434
1435 stout%gpu_states_mem = stin%gpu_states_mem
1436 stout%pack_states = stin%pack_states
1437
1438 stout%only_userdef_istates = stin%only_userdef_istates
1439
1440 if (.not. exclude_wfns_) then
1441 safe_allocate_source_a(stout%rho, stin%rho)
1442 end if
1443
1444 stout%uniform_occ = stin%uniform_occ
1445
1446 if (.not. optional_default(exclude_eigenval, .false.)) then
1447 safe_allocate_source_a(stout%eigenval, stin%eigenval)
1448 safe_allocate_source_a(stout%occ, stin%occ)
1449 safe_allocate_source_a(stout%spin, stin%spin)
1450 end if
1451
1452 ! the call to init_block is done at the end of this subroutine
1453 ! it allocates iblock, psib, block_is_local
1454 stout%group%nblocks = stin%group%nblocks
1455
1456 safe_allocate_source_a(stout%user_def_states, stin%user_def_states)
1457
1458 safe_allocate_source_a(stout%current, stin%current)
1459 safe_allocate_source_a(stout%current_kpt, stin%current_kpt)
1460 safe_allocate_source_a(stout%rho_core, stin%rho_core)
1461 safe_allocate_source_a(stout%frozen_rho, stin%frozen_rho)
1462 safe_allocate_source_a(stout%frozen_tau, stin%frozen_tau)
1463 safe_allocate_source_a(stout%frozen_gdens, stin%frozen_gdens)
1464 safe_allocate_source_a(stout%frozen_ldens, stin%frozen_ldens)
1465
1466 stout%fixed_occ = stin%fixed_occ
1467 stout%restart_fixed_occ = stin%restart_fixed_occ
1468
1469 stout%fixed_spins = stin%fixed_spins
1470
1471 stout%qtot = stin%qtot
1472 stout%val_charge = stin%val_charge
1473
1474 call smear_copy(stout%smear, stin%smear)
1475
1476 stout%parallel_in_states = stin%parallel_in_states
1477 call mpi_grp_copy(stout%mpi_grp, stin%mpi_grp)
1478 call mpi_grp_copy(stout%dom_st_kpt_mpi_grp, stin%dom_st_kpt_mpi_grp)
1479 call mpi_grp_copy(stout%st_kpt_mpi_grp, stin%st_kpt_mpi_grp)
1480 call mpi_grp_copy(stout%dom_st_mpi_grp, stin%dom_st_mpi_grp)
1481 safe_allocate_source_a(stout%node, stin%node)
1482 safe_allocate_source_a(stout%st_kpt_task, stin%st_kpt_task)
1483
1484#ifdef HAVE_SCALAPACK
1485 call blacs_proc_grid_copy(stin%dom_st_proc_grid, stout%dom_st_proc_grid)
1486#endif
1487
1488 call distributed_copy(stin%dist, stout%dist)
1489
1490 stout%scalapack_compatible = stin%scalapack_compatible
1491
1492 stout%lnst = stin%lnst
1493 stout%st_start = stin%st_start
1494 stout%st_end = stin%st_end
1495
1496 if (stin%parallel_in_states) call multicomm_all_pairs_copy(stout%ap, stin%ap)
1497
1498 stout%symmetrize_density = stin%symmetrize_density
1499
1500 if (.not. exclude_wfns_) call states_elec_group_copy(stin%d, stin%group, stout%group, special=special)
1501
1502 stout%packed = stin%packed
1503
1504 stout%randomization = stin%randomization
1505
1506 pop_sub(states_elec_copy)
1507 end subroutine states_elec_copy
1508
1509
1510 ! ---------------------------------------------------------
1512 !
1513 subroutine states_elec_end(st)
1514 type(states_elec_t), intent(inout) :: st
1515
1516 push_sub(states_elec_end)
1517
1518 call states_elec_dim_end(st%d)
1519
1520 call modelmb_particles_end(st%modelmbparticles)
1521
1522 ! this deallocates dpsi, zpsi, psib, iblock, iblock
1524
1525 safe_deallocate_a(st%user_def_states)
1526
1527 safe_deallocate_a(st%rho)
1528 safe_deallocate_a(st%eigenval)
1529
1530 safe_deallocate_a(st%current)
1531 safe_deallocate_a(st%current_para)
1532 safe_deallocate_a(st%current_dia)
1533 safe_deallocate_a(st%current_mag)
1534 safe_deallocate_a(st%current_kpt)
1535 safe_deallocate_a(st%rho_core)
1536 safe_deallocate_a(st%frozen_rho)
1537 safe_deallocate_a(st%frozen_tau)
1538 safe_deallocate_a(st%frozen_gdens)
1539 safe_deallocate_a(st%frozen_ldens)
1540 safe_deallocate_a(st%occ)
1541 safe_deallocate_a(st%spin)
1542 safe_deallocate_a(st%kweights)
1543
1544
1545#ifdef HAVE_SCALAPACK
1546 call blacs_proc_grid_end(st%dom_st_proc_grid)
1547#endif
1548
1549 call distributed_end(st%dist)
1550
1551 safe_deallocate_a(st%node)
1552 safe_deallocate_a(st%st_kpt_task)
1553
1554 if (st%parallel_in_states) then
1555 safe_deallocate_a(st%ap%schedule)
1556 end if
1557
1558 pop_sub(states_elec_end)
1559 end subroutine states_elec_end
1560
1561
1562 ! TODO(Alex) Issue #684. Abstract duplicate code in states_elec_generate_random, to get it to
1563 ! a point where it can be refactored.
1565 subroutine states_elec_generate_random(st, mesh, kpoints, ist_start_, ist_end_, ikpt_start_, ikpt_end_, normalized)
1566 type(states_elec_t), intent(inout) :: st
1567 class(mesh_t), intent(in) :: mesh
1568 type(kpoints_t), intent(in) :: kpoints
1569 integer, optional, intent(in) :: ist_start_
1570 integer, optional, intent(in) :: ist_end_
1571 integer, optional, intent(in) :: ikpt_start_
1572 integer, optional, intent(in) :: ikpt_end_
1573 logical, optional, intent(in) :: normalized
1575
1576 integer :: ist, ik, id, ist_start, ist_end, jst, ikpt_start, ikpt_end
1577 complex(real64) :: alpha, beta
1578 real(real64), allocatable :: dpsi(:, :)
1579 complex(real64), allocatable :: zpsi(:, :), zpsi2(:)
1580 integer :: ikpoint, ip
1581 type(batch_t) :: ffb
1582
1583 logical :: normalized_
1584
1585 normalized_ = optional_default(normalized, .true.)
1586
1588
1589 ist_start = optional_default(ist_start_, 1)
1590 ist_end = optional_default(ist_end_, st%nst)
1591 ikpt_start = optional_default(ikpt_start_, 1)
1592 ikpt_end = optional_default(ikpt_end_, st%nik)
1593
1594 safe_allocate(dpsi(1:mesh%np, 1:st%d%dim))
1595 if (states_are_complex(st)) then
1596 safe_allocate(zpsi(1:mesh%np, 1:st%d%dim))
1597 end if
1598
1599 select case (st%d%ispin)
1601
1602 do ik = ikpt_start, ikpt_end
1603 ikpoint = st%d%get_kpoint_index(ik)
1604 do ist = ist_start, ist_end
1605 if (states_are_real(st).or.kpoints_point_is_gamma(kpoints, ikpoint)) then
1606 if (st%randomization == par_independent) then
1607 call dmf_random(mesh, dpsi(:, 1), &
1608 pre_shift = mesh%pv%xlocal-1, &
1609 post_shift = mesh%pv%np_global - mesh%pv%xlocal - mesh%np + 1, &
1610 normalized = normalized)
1611 ! Ensures that the grid points are properly distributed in the domain parallel case
1612 if(mesh%parallel_in_domains) then
1613 call batch_init(ffb, dpsi(:,1))
1614 call dmesh_batch_exchange_points(mesh, ffb, backward_map = .true.)
1615 call ffb%end()
1616 end if
1617 else
1618 call dmf_random(mesh, dpsi(:, 1), normalized = normalized)
1619 end if
1620 if (.not. state_kpt_is_local(st, ist, ik)) cycle
1621 if (states_are_complex(st)) then !Gamma point
1622 do ip = 1, mesh%np
1623 zpsi(ip,1) = cmplx(dpsi(ip,1), m_zero, real64)
1624 end do
1625 call states_elec_set_state(st, mesh, ist, ik, zpsi)
1626 else
1627 call states_elec_set_state(st, mesh, ist, ik, dpsi)
1628 end if
1629 else
1630 if (st%randomization == par_independent) then
1631 call zmf_random(mesh, zpsi(:, 1), &
1632 pre_shift = mesh%pv%xlocal-1, &
1633 post_shift = mesh%pv%np_global - mesh%pv%xlocal - mesh%np + 1, &
1634 normalized = normalized)
1635 ! Ensures that the grid points are properly distributed in the domain parallel case
1636 if(mesh%parallel_in_domains) then
1637 call batch_init(ffb, zpsi(:,1))
1638 call zmesh_batch_exchange_points(mesh, ffb, backward_map = .true.)
1639 call ffb%end()
1640 end if
1641 else
1642 call zmf_random(mesh, zpsi(:, 1), normalized = normalized)
1643 end if
1644 if (.not. state_kpt_is_local(st, ist, ik)) cycle
1645 call states_elec_set_state(st, mesh, ist, ik, zpsi)
1646 end if
1647 end do
1648 end do
1649
1650 case (spinors)
1651
1652 assert(states_are_complex(st))
1653
1654 if (st%fixed_spins) then
1655
1656 do ik = ikpt_start, ikpt_end
1657 ikpoint = st%d%get_kpoint_index(ik)
1658 do ist = ist_start, ist_end
1659 if (kpoints_point_is_gamma(kpoints, ikpoint)) then
1660 if (st%randomization == par_independent) then
1661 call dmf_random(mesh, dpsi(:, 1), &
1662 pre_shift = mesh%pv%xlocal-1, &
1663 post_shift = mesh%pv%np_global - mesh%pv%xlocal - mesh%np + 1, &
1664 normalized = normalized)
1665 ! Ensures that the grid points are properly distributed in the domain parallel case
1666 if(mesh%parallel_in_domains) then
1667 call batch_init(ffb, dpsi(:,1))
1668 call dmesh_batch_exchange_points(mesh, ffb, backward_map = .true.)
1669 call ffb%end()
1670 end if
1671 else
1672 call dmf_random(mesh, dpsi(:, 1), normalized = normalized)
1673 if (.not. state_kpt_is_local(st, ist, ik)) cycle
1674 end if
1675 do ip = 1, mesh%np
1676 zpsi(ip,1) = cmplx(dpsi(ip,1), m_zero, real64)
1677 end do
1678 call states_elec_set_state(st, mesh, ist, ik, zpsi)
1679 else
1680 if (st%randomization == par_independent) then
1681 call zmf_random(mesh, zpsi(:, 1), &
1682 pre_shift = mesh%pv%xlocal-1, &
1683 post_shift = mesh%pv%np_global - mesh%pv%xlocal - mesh%np + 1, &
1684 normalized = normalized)
1685 ! Ensures that the grid points are properly distributed in the domain parallel case
1686 if(mesh%parallel_in_domains) then
1687 call batch_init(ffb, zpsi(:,1))
1688 call zmesh_batch_exchange_points(mesh, ffb, backward_map = .true.)
1689 call ffb%end()
1690 end if
1691 else
1692 call zmf_random(mesh, zpsi(:, 1), normalized = normalized)
1693 if (.not. state_kpt_is_local(st, ist, ik)) cycle
1694 end if
1695 end if
1696 if (.not. state_kpt_is_local(st, ist, ik)) cycle
1697 ! In this case, the spinors are made of a spatial part times a vector [alpha beta]^T in
1698 ! spin space (i.e., same spatial part for each spin component). So (alpha, beta)
1699 ! determines the spin values. The values of (alpha, beta) can be be obtained
1700 ! with simple formulae from <Sx>, <Sy>, <Sz>.
1701 !
1702 ! Note that here we orthonormalize the orbital part. This ensures that the spinors
1703 ! are untouched later in the general orthonormalization, and therefore the spin values
1704 ! of each spinor remain the same.
1705 safe_allocate(zpsi2(1:mesh%np))
1706 do jst = ist_start, ist - 1
1707 call states_elec_get_state(st, mesh, 1, jst, ik, zpsi2)
1708 zpsi(1:mesh%np, 1) = zpsi(1:mesh%np, 1) - zmf_dotp(mesh, zpsi(:, 1), zpsi2)*zpsi2(1:mesh%np)
1709 end do
1710 safe_deallocate_a(zpsi2)
1711
1712 call zmf_normalize(mesh, 1, zpsi)
1713 zpsi(1:mesh%np, 2) = zpsi(1:mesh%np, 1)
1714
1715 alpha = cmplx(sqrt(m_half + st%spin(3, ist, ik)), m_zero, real64)
1716 beta = cmplx(sqrt(m_one - abs(alpha)**2), m_zero, real64)
1717 if (abs(alpha) > m_zero) then
1718 beta = cmplx(st%spin(1, ist, ik) / abs(alpha), st%spin(2, ist, ik) / abs(alpha), real64)
1719 end if
1720 zpsi(1:mesh%np, 1) = alpha*zpsi(1:mesh%np, 1)
1721 zpsi(1:mesh%np, 2) = beta*zpsi(1:mesh%np, 2)
1722 call states_elec_set_state(st, mesh, ist, ik, zpsi)
1723 end do
1724 end do
1725 else
1726 do ik = ikpt_start, ikpt_end
1727 do ist = ist_start, ist_end
1728 do id = 1, st%d%dim
1729 if (st%randomization == par_independent) then
1730 call zmf_random(mesh, zpsi(:, id), &
1731 pre_shift = mesh%pv%xlocal-1, &
1732 post_shift = mesh%pv%np_global - mesh%pv%xlocal - mesh%np + 1, &
1733 normalized = .false.)
1734 ! Ensures that the grid points are properly distributed in the domain parallel case
1735 if(mesh%parallel_in_domains) then
1736 call batch_init(ffb, zpsi(:, id))
1737 call zmesh_batch_exchange_points(mesh, ffb, backward_map = .true.)
1738 call ffb%end()
1739 end if
1740 else
1741 call zmf_random(mesh, zpsi(:, id), normalized = .false.)
1742 end if
1743 end do
1744 ! We need to generate the wave functions even if not using them in order to be consistent
1745 ! with the random seed in parallel runs.
1746 if (.not. state_kpt_is_local(st, ist, ik)) cycle
1747 ! Note that mf_random normalizes each spin channel independently to 1.
1748 ! Therefore we need to renormalize the spinor:
1749 if (normalized_) call zmf_normalize(mesh, st%d%dim, zpsi)
1750 call states_elec_set_state(st, mesh, ist, ik, zpsi)
1751 end do
1752 end do
1753 end if
1754
1755 end select
1756
1757 safe_deallocate_a(dpsi)
1758 safe_deallocate_a(zpsi)
1759
1761 end subroutine states_elec_generate_random
1762
1763 ! ---------------------------------------------------------
1765 !
1766 subroutine states_elec_fermi(st, namespace, mesh, compute_spin)
1767 type(states_elec_t), intent(inout) :: st
1768 type(namespace_t), intent(in) :: namespace
1769 class(mesh_t), intent(in) :: mesh
1770 logical, optional, intent(in) :: compute_spin
1771
1773 integer :: ist, ik
1774 real(real64) :: charge
1775 complex(real64), allocatable :: zpsi(:, :)
1776
1777 push_sub(states_elec_fermi)
1778
1779 call smear_find_fermi_energy(st%smear, namespace, st%eigenval, st%occ, st%qtot, &
1780 st%nik, st%nst, st%kweights)
1781
1782 call smear_fill_occupations(st%smear, st%eigenval, st%occ, st%nik, st%nst)
1783
1784 ! check if everything is OK
1785 charge = m_zero
1786 do ist = 1, st%nst
1787 charge = charge + sum(st%occ(ist, 1:st%nik) * st%kweights(1:st%nik))
1788 end do
1789 if (abs(charge-st%qtot) > 1e-6_real64) then
1790 message(1) = 'Occupations do not integrate to total charge.'
1791 write(message(2), '(6x,f12.8,a,f12.8)') charge, ' != ', st%qtot
1792 call messages_warning(2, namespace=namespace)
1793 if (charge < m_epsilon) then
1794 message(1) = "There don't seem to be any electrons at all!"
1795 call messages_fatal(1, namespace=namespace)
1796 end if
1797 end if
1798
1799 if (st%d%ispin == spinors .and. optional_default(compute_spin,.true.)) then
1800 assert(states_are_complex(st))
1801
1802 st%spin(:,:,:) = m_zero
1803
1804 safe_allocate(zpsi(1:mesh%np, st%d%dim))
1805 do ik = st%d%kpt%start, st%d%kpt%end
1806 do ist = st%st_start, st%st_end
1807 call states_elec_get_state(st, mesh, ist, ik, zpsi)
1808 st%spin(1:3, ist, ik) = state_spin(mesh, zpsi)
1809 end do
1810 end do
1811 safe_deallocate_a(zpsi)
1812
1813 if (st%parallel_in_states .or. st%d%kpt%parallel) then
1814 call comm_allreduce(st%st_kpt_mpi_grp, st%spin)
1815 end if
1816
1817 end if
1818
1819 pop_sub(states_elec_fermi)
1820 end subroutine states_elec_fermi
1821
1822
1823 ! ---------------------------------------------------------
1825 !
1826 function states_elec_eigenvalues_sum(st, alt_eig) result(tot)
1827 type(states_elec_t), intent(in) :: st
1828 real(real64), optional, intent(in) :: alt_eig(st%st_start:, st%d%kpt%start:)
1830 real(real64) :: tot
1831
1832 integer :: ik
1833
1835
1836 tot = m_zero
1837 do ik = st%d%kpt%start, st%d%kpt%end
1838 if (present(alt_eig)) then
1839 tot = tot + st%kweights(ik) * sum(st%occ(st%st_start:st%st_end, ik) * &
1840 alt_eig(st%st_start:st%st_end, ik))
1841 else
1842 tot = tot + st%kweights(ik) * sum(st%occ(st%st_start:st%st_end, ik) * &
1843 st%eigenval(st%st_start:st%st_end, ik))
1844 end if
1845 end do
1846
1847 if (st%parallel_in_states .or. st%d%kpt%parallel) call comm_allreduce(st%st_kpt_mpi_grp, tot)
1848
1850 end function states_elec_eigenvalues_sum
1851
1852
1854 subroutine states_elec_distribute_nodes(st, namespace, mc)
1855 type(states_elec_t), intent(inout) :: st
1856 type(namespace_t), intent(in) :: namespace
1857 type(multicomm_t), intent(in) :: mc
1858
1859 logical :: default_scalapack_compatible
1860
1862
1863 ! TODO(Alex) Issue #820. This is superflous. These defaults are set in initialisation of
1864 ! states, and in the state distribution instance
1865 ! Defaults.
1866 st%node(:) = 0
1867 st%st_start = 1
1868 st%st_end = st%nst
1869 st%lnst = st%nst
1870 st%parallel_in_states = .false.
1871
1872 call mpi_grp_init(st%mpi_grp, mc%group_comm(p_strategy_states))
1873 call mpi_grp_init(st%dom_st_kpt_mpi_grp, mc%dom_st_kpt_comm)
1874 call mpi_grp_init(st%dom_st_mpi_grp, mc%dom_st_comm)
1875 call mpi_grp_init(st%st_kpt_mpi_grp, mc%st_kpt_comm)
1876
1877 default_scalapack_compatible = calc_mode_par%scalapack_compat() .and. .not. st%d%kpt%parallel
1878
1879 !%Variable ScaLAPACKCompatible
1880 !%Type logical
1881 !%Section Execution::Parallelization
1882 !%Description
1883 !% Whether to use a layout for states parallelization which is compatible with ScaLAPACK.
1884 !% The default is yes for <tt>CalculationMode = gs, unocc, go</tt> without k-point parallelization,
1885 !% and no otherwise. (Setting to other than default is experimental.)
1886 !% The value must be yes if any ScaLAPACK routines are called in the course of the run;
1887 !% it must be set by hand for <tt>td</tt> with <tt>TDDynamics = bo</tt>.
1888 !% This variable has no effect unless you are using states parallelization and have linked ScaLAPACK.
1889 !% Note: currently, use of ScaLAPACK is not compatible with task parallelization (<i>i.e.</i> slaves).
1890 !%End
1891 call parse_variable(namespace, 'ScaLAPACKCompatible', default_scalapack_compatible, st%scalapack_compatible)
1892
1893#ifdef HAVE_SCALAPACK
1894 if (default_scalapack_compatible .neqv. st%scalapack_compatible) then
1895 call messages_experimental('Setting ScaLAPACKCompatible to other than default', namespace=namespace)
1896 end if
1897
1898 if (st%scalapack_compatible) then
1899 if (multicomm_have_slaves(mc)) then
1900 call messages_not_implemented("ScaLAPACK usage with task parallelization (slaves)", namespace=namespace)
1901 end if
1902 call blacs_proc_grid_init(st%dom_st_proc_grid, st%dom_st_mpi_grp)
1903 end if
1904#else
1905 st%scalapack_compatible = .false.
1906#endif
1907
1909
1910#ifdef HAVE_MPI
1911 call multicomm_create_all_pairs(st%mpi_grp, st%ap)
1912#endif
1913
1914 if (st%nst < st%mpi_grp%size) then
1915 message(1) = "Have more processors than necessary"
1916 write(message(2),'(i4,a,i4,a)') st%mpi_grp%size, " processors and ", st%nst, " states."
1917 call messages_fatal(2, namespace=namespace)
1918 end if
1919
1920 call distributed_init(st%dist, st%nst, st%mpi_grp%comm, "states", scalapack_compat = st%scalapack_compatible)
1922 st%parallel_in_states = st%dist%parallel
1923
1924 ! TODO(Alex) Issue #820. Remove lnst, st_start, st_end and node, as they are all contained within dist
1925 st%st_start = st%dist%start
1926 st%st_end = st%dist%end
1927 st%lnst = st%dist%nlocal
1928 st%node(1:st%nst) = st%dist%node(1:st%nst)
1929
1930 end if
1931
1933
1935 end subroutine states_elec_distribute_nodes
1936
1937
1943 !
1944 subroutine states_elec_calc_quantities(gr, st, kpoints, nlcc, &
1945 kinetic_energy_density, paramagnetic_current, density_gradient, density_laplacian, &
1946 gi_kinetic_energy_density, st_end)
1947 type(grid_t), intent(in) :: gr
1948 type(states_elec_t), intent(in) :: st
1949 type(kpoints_t), intent(in) :: kpoints
1950 logical, intent(in) :: nlcc
1951 real(real64), contiguous, optional, target, intent(out) :: kinetic_energy_density(:,:)
1952 real(real64), contiguous, optional, target, intent(out) :: paramagnetic_current(:,:,:)
1953 real(real64), contiguous, optional, intent(out) :: density_gradient(:,:,:)
1954 real(real64), contiguous, optional, intent(out) :: density_laplacian(:,:)
1955 real(real64), contiguous, optional, intent(out) :: gi_kinetic_energy_density(:,:)
1956 integer, optional, intent(in) :: st_end
1957
1958 real(real64), pointer, contiguous :: jp(:, :, :)
1959 real(real64), pointer, contiguous :: tau(:, :)
1960 complex(real64), allocatable :: wf_psi(:,:), gwf_psi(:,:,:), wf_psi_conj(:,:), lwf_psi(:,:)
1961 real(real64), allocatable :: abs_wf_psi(:,:), abs_gwf_psi(:,:)
1962 complex(real64), allocatable :: psi_gpsi(:,:)
1963 complex(real64) :: c_tmp
1964 integer :: is, ik, ist, i_dim, st_dim, ii, st_end_
1965 real(real64) :: ww, kpoint(gr%der%dim)
1966 logical :: something_to_do
1967
1968 call profiling_in("STATES_CALC_QUANTITIES")
1969
1971
1972 st_end_ = min(st%st_end, optional_default(st_end, st%st_end))
1973
1974 something_to_do = present(kinetic_energy_density) .or. present(gi_kinetic_energy_density) .or. &
1975 present(paramagnetic_current) .or. present(density_gradient) .or. present(density_laplacian)
1976 assert(something_to_do)
1977
1978 safe_allocate( wf_psi(1:gr%np_part, 1:st%d%dim))
1979 safe_allocate( wf_psi_conj(1:gr%np_part, 1:st%d%dim))
1980 safe_allocate(gwf_psi(1:gr%np, 1:gr%der%dim, 1:st%d%dim))
1981 safe_allocate(abs_wf_psi(1:gr%np, 1:st%d%dim))
1982 safe_allocate(abs_gwf_psi(1:gr%np, 1:st%d%dim))
1983 safe_allocate(psi_gpsi(1:gr%np, 1:st%d%dim))
1984 if (present(density_laplacian)) then
1985 safe_allocate(lwf_psi(1:gr%np, 1:st%d%dim))
1986 end if
1987
1988 nullify(tau)
1989 if (present(kinetic_energy_density)) tau => kinetic_energy_density
1990
1991 nullify(jp)
1992 if (present(paramagnetic_current)) jp => paramagnetic_current
1993
1994 ! for the gauge-invariant kinetic energy density we need the
1995 ! current and the kinetic energy density
1996 if (present(gi_kinetic_energy_density)) then
1997 if (.not. present(paramagnetic_current) .and. states_are_complex(st)) then
1998 safe_allocate(jp(1:gr%np, 1:gr%der%dim, 1:st%d%nspin))
1999 end if
2000 if (.not. present(kinetic_energy_density)) then
2001 safe_allocate(tau(1:gr%np, 1:st%d%nspin))
2002 end if
2003 end if
2004
2005 if (associated(tau)) tau = m_zero
2006 if (associated(jp)) jp = m_zero
2007 if (present(density_gradient)) density_gradient(:,:,:) = m_zero
2008 if (present(density_laplacian)) density_laplacian(:,:) = m_zero
2009 if (present(gi_kinetic_energy_density)) gi_kinetic_energy_density = m_zero
2010
2011 do ik = st%d%kpt%start, st%d%kpt%end
2012
2013 kpoint(1:gr%der%dim) = kpoints%get_point(st%d%get_kpoint_index(ik))
2014 is = st%d%get_spin_index(ik)
2015
2016 do ist = st%st_start, st_end_
2017 ww = st%kweights(ik)*st%occ(ist, ik)
2018 if (abs(ww) <= m_epsilon) cycle
2019
2020 ! all calculations will be done with complex wavefunctions
2021 call states_elec_get_state(st, gr, ist, ik, wf_psi)
2022
2023 do st_dim = 1, st%d%dim
2024 call boundaries_set(gr%der%boundaries, gr, wf_psi(:, st_dim))
2025 end do
2026
2027 ! calculate gradient of the wavefunction
2028 do st_dim = 1, st%d%dim
2029 call zderivatives_grad(gr%der, wf_psi(:,st_dim), gwf_psi(:,:,st_dim), set_bc = .false.)
2030 end do
2031
2032 ! calculate the Laplacian of the wavefunction
2033 if (present(density_laplacian)) then
2034 do st_dim = 1, st%d%dim
2035 call zderivatives_lapl(gr%der, wf_psi(:,st_dim), lwf_psi(:,st_dim), ghost_update = .false., set_bc = .false.)
2036 end do
2037 end if
2038
2039 ! We precompute some quantites, to avoid to compute it many times
2040 wf_psi_conj(1:gr%np, 1:st%d%dim) = conjg(wf_psi(1:gr%np,1:st%d%dim))
2041 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)
2042
2043 if (present(density_laplacian)) then
2044 density_laplacian(1:gr%np, is) = density_laplacian(1:gr%np, is) + &
2045 ww * m_two*real(wf_psi_conj(1:gr%np, 1) * lwf_psi(1:gr%np, 1), real64)
2046 if (st%d%ispin == spinors) then
2047 density_laplacian(1:gr%np, 2) = density_laplacian(1:gr%np, 2) + &
2048 ww * m_two*real(wf_psi_conj(1:gr%np, 2) * lwf_psi(1:gr%np, 2), real64)
2049 !$omp parallel do private(c_tmp)
2050 do ii = 1, gr%np
2051 c_tmp = ww*(lwf_psi(ii, 1) * wf_psi_conj(ii, 2) + wf_psi(ii, 1) * conjg(lwf_psi(ii, 2)))
2052 density_laplacian(ii, 3) = density_laplacian(ii, 3) + real(c_tmp, real64)
2053 density_laplacian(ii, 4) = density_laplacian(ii, 4) + aimag(c_tmp)
2054 end do
2055 end if
2056 end if
2057
2058 if (associated(tau)) then
2059 tau(1:gr%np, is) = tau(1:gr%np, is) &
2060 + ww * sum(kpoint(1:gr%der%dim)**2) * abs_wf_psi(1:gr%np, 1)
2061 if (st%d%ispin == spinors) then
2062 tau(1:gr%np, 2) = tau(1:gr%np, 2) &
2063 + ww * sum(kpoint(1:gr%der%dim)**2) * abs_wf_psi(1:gr%np, 2)
2064
2065 !$omp parallel do private(c_tmp)
2066 do ii = 1, gr%np
2067 c_tmp = ww * sum(kpoint(1:gr%der%dim)**2) * wf_psi(ii, 1) * wf_psi_conj(ii, 2)
2068 tau(ii, 3) = tau(ii, 3) + real(c_tmp, real64)
2069 tau(ii, 4) = tau(ii, 4) + aimag(c_tmp)
2070 end do
2071 end if
2072 end if
2073
2074 do i_dim = 1, gr%der%dim
2075
2076 ! We precompute some quantites, to avoid to compute them many times
2077 psi_gpsi(1:gr%np, 1:st%d%dim) = wf_psi_conj(1:gr%np, 1:st%d%dim) &
2078 * gwf_psi(1:gr%np,i_dim,1:st%d%dim)
2079 abs_gwf_psi(1:gr%np, 1:st%d%dim) = real(conjg(gwf_psi(1:gr%np, i_dim, 1:st%d%dim)) &
2080 * gwf_psi(1:gr%np, i_dim, 1:st%d%dim), real64)
2081
2082 if (present(density_gradient)) then
2083 density_gradient(1:gr%np, i_dim, is) = density_gradient(1:gr%np, i_dim, is) &
2084 + ww * m_two * real(psi_gpsi(1:gr%np, 1), real64)
2085 if (st%d%ispin == spinors) then
2086 density_gradient(1:gr%np, i_dim, 2) = density_gradient(1:gr%np, i_dim, 2) &
2087 + ww * m_two*real(psi_gpsi(1:gr%np, 2), real64)
2088 !$omp parallel do private(c_tmp)
2089 do ii = 1, gr%np
2090 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)))
2091 density_gradient(ii, i_dim, 3) = density_gradient(ii, i_dim, 3) + real(c_tmp, real64)
2092 density_gradient(ii, i_dim, 4) = density_gradient(ii, i_dim, 4) + aimag(c_tmp)
2093 end do
2094 end if
2095 end if
2096
2097 if (present(density_laplacian)) then
2098 call lalg_axpy(gr%np, ww*m_two, abs_gwf_psi(:,1), density_laplacian(:,is))
2099 if (st%d%ispin == spinors) then
2100 call lalg_axpy(gr%np, ww*m_two, abs_gwf_psi(:,2), density_laplacian(:,2))
2101 !$omp parallel do private(c_tmp)
2102 do ii = 1, gr%np
2103 c_tmp = m_two * ww * gwf_psi(ii, i_dim, 1) * conjg(gwf_psi(ii, i_dim, 2))
2104 density_laplacian(ii, 3) = density_laplacian(ii, 3) + real(c_tmp, real64)
2105 density_laplacian(ii, 4) = density_laplacian(ii, 4) + aimag(c_tmp)
2106 end do
2107 end if
2108 end if
2109
2110 ! the expression for the paramagnetic current with spinors is
2111 ! j = ( jp(1) jp(3) + i jp(4) )
2112 ! (-jp(3) + i jp(4) jp(2) )
2113 if (associated(jp)) then
2114 if (.not.(states_are_real(st))) then
2115 jp(1:gr%np, i_dim, is) = jp(1:gr%np, i_dim, is) + &
2116 ww*(aimag(psi_gpsi(1:gr%np, 1)) - abs_wf_psi(1:gr%np, 1)*kpoint(i_dim))
2117 if (st%d%ispin == spinors) then
2118 jp(1:gr%np, i_dim, 2) = jp(1:gr%np, i_dim, 2) + &
2119 ww*( aimag(psi_gpsi(1:gr%np, 2)) - abs_wf_psi(1:gr%np, 2)*kpoint(i_dim))
2120 !$omp parallel do private(c_tmp)
2121 do ii = 1, gr%np
2122 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)) &
2123 - m_two * m_zi*wf_psi(ii, 1)*wf_psi_conj(ii, 2)*kpoint(i_dim))
2124 jp(ii, i_dim, 3) = jp(ii, i_dim, 3) + real(c_tmp, real64)
2125 jp(ii, i_dim, 4) = jp(ii, i_dim, 4) + aimag(c_tmp)
2126 end do
2127 end if
2128 end if
2129 end if
2130
2131 ! the expression for the paramagnetic current with spinors is
2132 ! t = ( tau(1) tau(3) + i tau(4) )
2133 ! ( tau(3) - i tau(4) tau(2) )
2134 if (associated(tau)) then
2135 tau(1:gr%np, is) = tau(1:gr%np, is) + ww*(abs_gwf_psi(1:gr%np,1) &
2136 - m_two*aimag(psi_gpsi(1:gr%np, 1))*kpoint(i_dim))
2137 if (st%d%ispin == spinors) then
2138 tau(1:gr%np, 2) = tau(1:gr%np, 2) + ww*(abs_gwf_psi(1:gr%np, 2) &
2139 - m_two*aimag(psi_gpsi(1:gr%np, 2))*kpoint(i_dim))
2140 !$omp parallel do private(c_tmp)
2141 do ii = 1, gr%np
2142 c_tmp = ww * ( gwf_psi(ii, i_dim, 1)*conjg(gwf_psi(ii, i_dim, 2)) &
2143 + m_zi * (gwf_psi(ii,i_dim,1)*wf_psi_conj(ii, 2) &
2144 - wf_psi(ii, 1)*conjg(gwf_psi(ii,i_dim,2)))*kpoint(i_dim))
2145 tau(ii, 3) = tau(ii, 3) + real(c_tmp, real64)
2146 tau(ii, 4) = tau(ii, 4) + aimag(c_tmp)
2147 end do
2148 end if
2149 end if
2150
2151 end do
2152
2153 end do
2154 end do
2155
2156 safe_deallocate_a(wf_psi_conj)
2157 safe_deallocate_a(abs_wf_psi)
2158 safe_deallocate_a(abs_gwf_psi)
2159 safe_deallocate_a(psi_gpsi)
2160
2161 if (.not. present(gi_kinetic_energy_density)) then
2162 if (.not. present(paramagnetic_current)) then
2163 safe_deallocate_p(jp)
2164 end if
2165 if (.not. present(kinetic_energy_density)) then
2166 safe_deallocate_p(tau)
2167 end if
2168 end if
2169
2170 if (st%parallel_in_states .or. st%d%kpt%parallel) call reduce_all(st%st_kpt_mpi_grp)
2171
2172 ! We have to symmetrize everything as they are calculated from the
2173 ! wavefunctions.
2174 ! This must be done before compute the gauge-invariant kinetic energy density
2175 if (st%symmetrize_density) then
2176 do is = 1, st%d%nspin
2177 if (associated(tau)) then
2178 call dgrid_symmetrize_scalar_field(gr, tau(:, is), suppress_warning = .true.)
2179 end if
2180
2181 if (present(density_laplacian)) then
2182 call dgrid_symmetrize_scalar_field(gr, density_laplacian(:, is), suppress_warning = .true.)
2183 end if
2184
2185 if (associated(jp)) then
2186 call dgrid_symmetrize_vector_field(gr, jp(:, :, is), suppress_warning = .true.)
2187 end if
2188
2189 if (present(density_gradient)) then
2190 call dgrid_symmetrize_vector_field(gr, density_gradient(:, :, is), suppress_warning = .true.)
2191 end if
2192 end do
2193 end if
2194
2195
2196 if (allocated(st%rho_core) .and. nlcc .and. (present(density_laplacian) .or. present(density_gradient))) then
2197 do ii = 1, gr%np
2198 wf_psi(ii, 1) = st%rho_core(ii)/st%d%spin_channels
2199 end do
2200
2201 call boundaries_set(gr%der%boundaries, gr, wf_psi(:, 1))
2202
2203 if (present(density_gradient)) then
2204 ! calculate gradient of the NLCC
2205 call zderivatives_grad(gr%der, wf_psi(:,1), gwf_psi(:,:,1), set_bc = .false.)
2206 do is = 1, st%d%spin_channels
2207 density_gradient(1:gr%np, 1:gr%der%dim, is) = density_gradient(1:gr%np, 1:gr%der%dim, is) + &
2208 real(gwf_psi(1:gr%np, 1:gr%der%dim,1))
2209 end do
2210 end if
2211
2212 ! calculate the Laplacian of the wavefunction
2213 if (present(density_laplacian)) then
2214 call zderivatives_lapl(gr%der, wf_psi(:,1), lwf_psi(:,1), set_bc = .false.)
2215
2216 do is = 1, st%d%spin_channels
2217 density_laplacian(1:gr%np, is) = density_laplacian(1:gr%np, is) + real(lwf_psi(1:gr%np, 1))
2218 end do
2219 end if
2220 end if
2221
2222 !If we freeze some of the orbitals, we need to had the contributions here
2223 !Only in the case we are not computing it
2224 if (allocated(st%frozen_tau) .and. .not. present(st_end)) then
2225 call lalg_axpy(gr%np, st%d%nspin, m_one, st%frozen_tau, tau)
2226 end if
2227 if (allocated(st%frozen_gdens) .and. .not. present(st_end)) then
2228 do is = 1, st%d%nspin
2229 call lalg_axpy(gr%np, gr%der%dim, m_one, st%frozen_gdens(:,:,is), density_gradient(:,:,is))
2230 end do
2231 end if
2232 if (allocated(st%frozen_tau) .and. .not. present(st_end)) then
2233 call lalg_axpy(gr%np, st%d%nspin, m_one, st%frozen_ldens, density_laplacian)
2234 end if
2235
2236 safe_deallocate_a(wf_psi)
2237 safe_deallocate_a(gwf_psi)
2238 safe_deallocate_a(lwf_psi)
2239
2240
2241 !We compute the gauge-invariant kinetic energy density
2242 if (present(gi_kinetic_energy_density)) then
2243 do is = 1, st%d%nspin
2244 assert(associated(tau))
2245 gi_kinetic_energy_density(1:gr%np, is) = tau(1:gr%np, is)
2246 end do
2247 if (states_are_complex(st)) then
2248 assert(associated(jp))
2249 if (st%d%ispin /= spinors) then
2250 do is = 1, st%d%nspin
2251 !$omp parallel do
2252 do ii = 1, gr%np
2253 if (st%rho(ii, is) < 1.0e-7_real64) cycle
2254 gi_kinetic_energy_density(ii, is) = &
2255 gi_kinetic_energy_density(ii, is) - sum(jp(ii,1:gr%der%dim, is)**2)/st%rho(ii, is)
2256 end do
2257 end do
2258 else ! Note that this is only the U(1) part of the gauge-invariant KED
2259 !$omp parallel do
2260 do ii = 1, gr%np
2261 gi_kinetic_energy_density(ii, 1) = gi_kinetic_energy_density(ii, 1) &
2262 - sum(jp(ii,1:gr%der%dim, 1)**2)/(safe_tol(st%rho(ii, 1),m_epsilon))
2263 gi_kinetic_energy_density(ii, 2) = gi_kinetic_energy_density(ii, 2) &
2264 - sum(jp(ii,1:gr%der%dim, 2)**2)/(safe_tol(st%rho(ii, 2),m_epsilon))
2265 gi_kinetic_energy_density(ii, 3) = &
2266 gi_kinetic_energy_density(ii, 3) - sum(jp(ii,1:gr%der%dim, 3)**2 + jp(ii,1:gr%der%dim, 4)**2) &
2267 /(safe_tol((st%rho(ii, 3)**2 + st%rho(ii, 4)**2), m_epsilon))*st%rho(ii, 3)
2268 gi_kinetic_energy_density(ii, 4) = &
2269 gi_kinetic_energy_density(ii, 4) + sum(jp(ii,1:gr%der%dim, 3)**2 + jp(ii,1:gr%der%dim, 4)**2) &
2270 /(safe_tol((st%rho(ii, 3)**2 + st%rho(ii, 4)**2), m_epsilon))*st%rho(ii, 4)
2271 end do
2272 end if
2273 end if
2274 end if
2275
2276 if (.not. present(kinetic_energy_density)) then
2277 safe_deallocate_p(tau)
2278 end if
2279 if (.not. present(paramagnetic_current)) then
2280 safe_deallocate_p(jp)
2281 end if
2282
2283
2285
2286 call profiling_out("STATES_CALC_QUANTITIES")
2287
2288 contains
2289
2290 subroutine reduce_all(grp)
2291 type(mpi_grp_t), intent(in) :: grp
2292
2294
2295 if (associated(tau)) call comm_allreduce(grp, tau, dim = (/gr%np, st%d%nspin/))
2296
2297 if (present(density_laplacian)) call comm_allreduce(grp, density_laplacian, dim = (/gr%np, st%d%nspin/))
2298
2299 do is = 1, st%d%nspin
2300 if (associated(jp)) call comm_allreduce(grp, jp(:, :, is), dim = (/gr%np, gr%der%dim/))
2301
2302 if (present(density_gradient)) then
2303 call comm_allreduce(grp, density_gradient(:, :, is), dim = (/gr%np, gr%der%dim/))
2304 end if
2305 end do
2306
2308 end subroutine reduce_all
2309
2310 end subroutine states_elec_calc_quantities
2311
2312
2313 ! ---------------------------------------------------------
2315 !
2316 function state_spin(mesh, f1) result(spin)
2317 type(mesh_t), intent(in) :: mesh
2318 complex(real64), intent(in) :: f1(:, :)
2319 real(real64) :: spin(1:3)
2320
2321 complex(real64) :: z
2322
2323 push_sub(state_spin)
2324
2325 z = zmf_dotp(mesh, f1(:, 1) , f1(:, 2))
2326
2327 spin(1) = m_two*real(z, real64)
2328 spin(2) = m_two*aimag(z)
2329 spin(3) = zmf_nrm2(mesh, f1(:, 1))**2 - zmf_nrm2(mesh, f1(:, 2))**2
2330 spin = m_half*spin ! spin is half the sigma matrix.
2331
2332 pop_sub(state_spin)
2333 end function state_spin
2334
2335 ! ---------------------------------------------------------
2337 !
2338 logical function state_is_local(st, ist)
2339 type(states_elec_t), intent(in) :: st
2340 integer, intent(in) :: ist
2341
2342 push_sub(state_is_local)
2343
2344 state_is_local = ist >= st%st_start.and.ist <= st%st_end
2345
2346 pop_sub(state_is_local)
2347 end function state_is_local
2348
2349 ! ---------------------------------------------------------
2351 !
2352 logical function state_kpt_is_local(st, ist, ik)
2353 type(states_elec_t), intent(in) :: st
2354 integer, intent(in) :: ist
2355 integer, intent(in) :: ik
2356
2357 push_sub(state_kpt_is_local)
2358
2359 state_kpt_is_local = ist >= st%st_start .and. ist <= st%st_end .and. &
2360 ik >= st%d%kpt%start .and. ik <= st%d%kpt%end
2361
2362 pop_sub(state_kpt_is_local)
2363 end function state_kpt_is_local
2364
2365
2366 ! ---------------------------------------------------------
2368 real(real64) function states_elec_wfns_memory(st, mesh) result(memory)
2369 type(states_elec_t), intent(in) :: st
2370 class(mesh_t), intent(in) :: mesh
2371
2372 push_sub(states_elec_wfns_memory)
2373 memory = m_zero
2374
2375 ! orbitals
2376 memory = memory + real64*real(mesh%np_part_global, real64) *st%d%dim*real(st%nst, real64) *st%d%kpt%nglobal
2377
2379 end function states_elec_wfns_memory
2380
2381 ! ---------------------------------------------------------
2383 !
2384 subroutine states_elec_pack(st, copy)
2385 class(states_elec_t), intent(inout) :: st
2386 logical, optional, intent(in) :: copy
2387
2388 integer :: iqn, ib
2389 integer(int64) :: max_mem, mem
2390
2391 push_sub(states_elec_pack)
2392
2393 ! nothing to do, already packed
2394 if (st%packed) then
2395 pop_sub(states_elec_pack)
2396 return
2397 end if
2398
2399 st%packed = .true.
2400
2401 if (accel_is_enabled()) then
2402 max_mem = accel_global_memory_size()
2403
2404 if (st%gpu_states_mem > m_one) then
2405 max_mem = int(st%gpu_states_mem, int64)*(1024_8)**2
2406 else if (st%gpu_states_mem < 0.0_real64) then
2407 max_mem = max_mem + int(st%gpu_states_mem, int64)*(1024_8)**2
2408 else
2409 max_mem = int(st%gpu_states_mem*real(max_mem, real64) , int64)
2410 end if
2411 else
2412 max_mem = huge(max_mem)
2413 end if
2414
2415 mem = 0
2416 qnloop: do iqn = st%d%kpt%start, st%d%kpt%end
2417 do ib = st%group%block_start, st%group%block_end
2418
2419 mem = mem + st%group%psib(ib, iqn)%pack_total_size()
2420
2421 if (mem > max_mem) then
2422 call messages_write('Not enough CL device memory to store all states simultaneously.', new_line = .true.)
2423 call messages_write('Only ')
2424 call messages_write(ib - st%group%block_start)
2425 call messages_write(' of ')
2426 call messages_write(st%group%block_end - st%group%block_start + 1)
2427 call messages_write(' blocks will be stored in device memory.', new_line = .true.)
2428 call messages_warning()
2429 exit qnloop
2430 end if
2431
2432 call st%group%psib(ib, iqn)%do_pack(copy)
2433 end do
2434 end do qnloop
2436 pop_sub(states_elec_pack)
2437 end subroutine states_elec_pack
2438
2439 ! ------------------------------------------------------------
2441 !
2442 subroutine states_elec_unpack(st, copy)
2443 class(states_elec_t), intent(inout) :: st
2444 logical, optional, intent(in) :: copy
2445
2446 integer :: iqn, ib
2447
2448 push_sub(states_elec_unpack)
2449
2450 if (st%packed) then
2451 st%packed = .false.
2452
2453 do iqn = st%d%kpt%start, st%d%kpt%end
2454 do ib = st%group%block_start, st%group%block_end
2455 if (st%group%psib(ib, iqn)%is_packed()) call st%group%psib(ib, iqn)%do_unpack(copy)
2456 end do
2457 end do
2458 end if
2459
2460 pop_sub(states_elec_unpack)
2461 end subroutine states_elec_unpack
2462
2463 ! -----------------------------------------------------------
2465 !
2466 subroutine states_elec_write_info(st, namespace)
2467 class(states_elec_t), intent(in) :: st
2468 type(namespace_t), intent(in) :: namespace
2469
2470 push_sub(states_elec_write_info)
2471
2472 call messages_print_with_emphasis(msg="States", namespace=namespace)
2473
2474 write(message(1), '(a,f12.3)') 'Total electronic charge = ', st%qtot
2475 write(message(2), '(a,i8)') 'Number of states = ', st%nst
2476 write(message(3), '(a,i8)') 'States block-size = ', st%block_size
2477 call messages_info(3, namespace=namespace)
2478
2479 call messages_print_with_emphasis(namespace=namespace)
2480
2481 pop_sub(states_elec_write_info)
2482 end subroutine states_elec_write_info
2483
2485 subroutine states_elec_set_zero(st)
2486 class(states_elec_t), intent(inout) :: st
2487
2488 integer :: iqn, ib
2489
2490 push_sub(states_elec_set_zero)
2491
2492 do iqn = st%d%kpt%start, st%d%kpt%end
2493 do ib = st%group%block_start, st%group%block_end
2494 call batch_set_zero(st%group%psib(ib, iqn))
2495 end do
2496 end do
2497
2498 pop_sub(states_elec_set_zero)
2499 end subroutine states_elec_set_zero
2500
2501 ! ------------------------------------------------------------
2503 !
2504 integer pure function states_elec_block_min(st, ib) result(range)
2505 type(states_elec_t), intent(in) :: st
2506 integer, intent(in) :: ib
2507
2508 range = st%group%block_range(ib, 1)
2509 end function states_elec_block_min
2510
2511 ! ------------------------------------------------------------
2513 !
2514 integer pure function states_elec_block_max(st, ib) result(range)
2515 type(states_elec_t), intent(in) :: st
2516 integer, intent(in) :: ib
2517
2518 range = st%group%block_range(ib, 2)
2519 end function states_elec_block_max
2520
2521 ! ------------------------------------------------------------
2523 !
2524 integer pure function states_elec_block_size(st, ib) result(size)
2525 type(states_elec_t), intent(in) :: st
2526 integer, intent(in) :: ib
2527
2528 size = st%group%block_size(ib)
2529 end function states_elec_block_size
2530
2531 ! ---------------------------------------------------------
2533 !
2534 subroutine states_elec_count_pairs(st, namespace, n_pairs, n_occ, n_unocc, is_included, is_frac_occ)
2535 type(states_elec_t), intent(in) :: st
2536 type(namespace_t), intent(in) :: namespace
2537 integer, intent(out) :: n_pairs
2538 integer, intent(out) :: n_occ(:)
2539 integer, intent(out) :: n_unocc(:)
2540 logical, allocatable, intent(out) :: is_included(:,:,:)
2542 logical, intent(out) :: is_frac_occ
2543
2544 integer :: ik, ist, ast, n_filled, n_partially_filled, n_half_filled
2545 character(len=80) :: nst_string, default, wfn_list
2546 real(real64) :: energy_window
2547
2548 push_sub(states_elec_count_pairs)
2550 is_frac_occ = .false.
2551 do ik = 1, st%nik
2552 call occupied_states(st, namespace, ik, n_filled, n_partially_filled, n_half_filled)
2553 if (n_partially_filled > 0 .or. n_half_filled > 0) is_frac_occ = .true.
2554 n_occ(ik) = n_filled + n_partially_filled + n_half_filled
2555 n_unocc(ik) = st%nst - n_filled
2556 ! when we implement occupations, partially occupied levels need to be counted as both occ and unocc.
2557 end do
2558
2559 !%Variable CasidaKSEnergyWindow
2560 !%Type float
2561 !%Section Linear Response::Casida
2562 !%Description
2563 !% An alternative to <tt>CasidaKohnShamStates</tt> for specifying which occupied-unoccupied
2564 !% transitions will be used: all those whose eigenvalue differences are less than this
2565 !% number will be included. If a value less than 0 is supplied, this criterion will not be used.
2566 !%End
2567
2568 call parse_variable(namespace, 'CasidaKSEnergyWindow', -m_one, energy_window, units_inp%energy)
2569
2570 !%Variable CasidaKohnShamStates
2571 !%Type string
2572 !%Section Linear Response::Casida
2573 !%Default all states
2574 !%Description
2575 !% The calculation of the excitation spectrum of a system in the Casida frequency-domain
2576 !% formulation of linear-response time-dependent density functional theory (TDDFT)
2577 !% implies the use of a basis set of occupied/unoccupied Kohn-Sham orbitals. This
2578 !% basis set should, in principle, include all pairs formed by all occupied states,
2579 !% and an infinite number of unoccupied states. In practice, one has to truncate this
2580 !% basis set, selecting a number of occupied and unoccupied states that will form the
2581 !% pairs. These states are specified with this variable. If there are, say, 15 occupied
2582 !% states, and one sets this variable to the value "10-18", this means that occupied
2583 !% states from 10 to 15, and unoccupied states from 16 to 18 will be considered.
2584 !%
2585 !% This variable is a string in list form, <i>i.e.</i> expressions such as "1,2-5,8-15" are
2586 !% valid. You should include a non-zero number of unoccupied states and a non-zero number
2587 !% of occupied states.
2588 !%End
2589
2590 n_pairs = 0
2591 safe_allocate(is_included(maxval(n_occ), minval(n_occ) + 1:st%nst , st%nik))
2592 is_included(:,:,:) = .false.
2593
2594 if (energy_window < m_zero) then
2595 write(nst_string,'(i6)') st%nst
2596 write(default,'(a,a)') "1-", trim(adjustl(nst_string))
2597 call parse_variable(namespace, 'CasidaKohnShamStates', default, wfn_list)
2598
2599 write(message(1),'(a,a)') "Info: States that form the basis: ", trim(wfn_list)
2600 call messages_info(1, namespace=namespace)
2601
2602 ! count pairs
2603 n_pairs = 0
2604 do ik = 1, st%nik
2605 do ast = n_occ(ik) + 1, st%nst
2606 if (loct_isinstringlist(ast, wfn_list)) then
2607 do ist = 1, n_occ(ik)
2608 if (loct_isinstringlist(ist, wfn_list)) then
2609 n_pairs = n_pairs + 1
2610 is_included(ist, ast, ik) = .true.
2611 end if
2612 end do
2613 end if
2614 end do
2615 end do
2616
2617 else ! using CasidaKSEnergyWindow
2618
2619 write(message(1),'(a,f12.6,a)') "Info: including transitions with energy < ", &
2620 units_from_atomic(units_out%energy, energy_window), trim(units_abbrev(units_out%energy))
2621 call messages_info(1, namespace=namespace)
2622
2623 ! count pairs
2624 n_pairs = 0
2625 do ik = 1, st%nik
2626 do ast = n_occ(ik) + 1, st%nst
2627 do ist = 1, n_occ(ik)
2628 if (st%eigenval(ast, ik) - st%eigenval(ist, ik) < energy_window) then
2629 n_pairs = n_pairs + 1
2630 is_included(ist, ast, ik) = .true.
2631 end if
2632 end do
2633 end do
2634 end do
2635
2636 end if
2637
2639 end subroutine states_elec_count_pairs
2640
2641
2657 !
2658 subroutine occupied_states(st, namespace, ik, n_filled, n_partially_filled, n_half_filled, &
2659 filled, partially_filled, half_filled)
2660 type(states_elec_t), intent(in) :: st
2661 type(namespace_t), intent(in) :: namespace
2662 integer, intent(in) :: ik
2663 integer, intent(out) :: n_filled, n_partially_filled, n_half_filled
2664 integer, optional, intent(out) :: filled(:), partially_filled(:), half_filled(:)
2665
2666 integer :: ist
2667
2668 push_sub(occupied_states)
2669
2670 if (present(filled)) filled(:) = 0
2671 if (present(partially_filled)) partially_filled(:) = 0
2672 if (present(half_filled)) half_filled(:) = 0
2673 n_filled = 0
2674 n_partially_filled = 0
2675 n_half_filled = 0
2676
2677 select case (st%d%ispin)
2678 case (unpolarized)
2679 do ist = 1, st%nst
2680 if (abs(st%occ(ist, ik) - m_two) < m_min_occ) then
2681 n_filled = n_filled + 1
2682 if (present(filled)) filled(n_filled) = ist
2683 else if (abs(st%occ(ist, ik) - m_one) < m_min_occ) then
2684 n_half_filled = n_half_filled + 1
2685 if (present(half_filled)) half_filled(n_half_filled) = ist
2686 else if (st%occ(ist, ik) > m_min_occ) then
2687 n_partially_filled = n_partially_filled + 1
2688 if (present(partially_filled)) partially_filled(n_partially_filled) = ist
2689 else if (abs(st%occ(ist, ik)) > m_min_occ) then
2690 write(message(1),*) 'Internal error in occupied_states: Illegal occupation value ', st%occ(ist, ik)
2691 call messages_fatal(1, namespace=namespace)
2692 end if
2693 end do
2694 case (spin_polarized, spinors)
2695 do ist = 1, st%nst
2696 if (abs(st%occ(ist, ik)-m_one) < m_min_occ) then
2697 n_filled = n_filled + 1
2698 if (present(filled)) filled(n_filled) = ist
2699 else if (st%occ(ist, ik) > m_min_occ) then
2700 n_partially_filled = n_partially_filled + 1
2701 if (present(partially_filled)) partially_filled(n_partially_filled) = ist
2702 else if (abs(st%occ(ist, ik)) > m_min_occ) then
2703 write(message(1),*) 'Internal error in occupied_states: Illegal occupation value ', st%occ(ist, ik)
2704 call messages_fatal(1, namespace=namespace)
2705 end if
2706 end do
2707 end select
2708
2709 pop_sub(occupied_states)
2710 end subroutine occupied_states
2711
2712
2714 subroutine kpoints_distribute(this, mc)
2715 type(states_elec_t), intent(inout) :: this
2716 type(multicomm_t), intent(in) :: mc
2717
2718 push_sub(kpoints_distribute)
2719 call distributed_init(this%d%kpt, this%nik, mc%group_comm(p_strategy_kpoints), "k-points")
2720
2721 pop_sub(kpoints_distribute)
2722 end subroutine kpoints_distribute
2723
2724
2725 ! TODO(Alex) Issue #824. Consider converting this to a function to returns `st_kpt_task`
2726 ! as this is only called in a couple of places, or package with the `st_kpt_mpi_grp` when split
2727 ! from st instance
2730 type(states_elec_t), intent(inout) :: st
2731
2733
2734 if (.not. allocated(st%st_kpt_task)) then
2735 safe_allocate(st%st_kpt_task(0:st%st_kpt_mpi_grp%size-1, 1:4))
2736 end if
2737
2738 st%st_kpt_task(0:st%st_kpt_mpi_grp%size-1, :) = 0
2739 st%st_kpt_task(st%st_kpt_mpi_grp%rank, :) = [st%st_start, st%st_end, st%d%kpt%start, st%d%kpt%end]
2740
2741 if (st%parallel_in_states .or. st%d%kpt%parallel) then
2742 call comm_allreduce(st%st_kpt_mpi_grp, st%st_kpt_task)
2743 end if
2744
2747
2748 ! ---------------------------------------------------------
2750 !
2751 subroutine states_elec_choose_kpoints(st, kpoints, namespace)
2752 type(states_elec_t), target, intent(inout) :: st
2753 type(kpoints_t), intent(in) :: kpoints
2754 type(namespace_t), intent(in) :: namespace
2755
2756 integer :: ik, iq
2757 type(states_elec_dim_t), pointer :: dd
2758
2759
2761
2762 dd => st%d
2763
2764 st%nik = kpoints_number(kpoints)
2765
2766 if (dd%ispin == spin_polarized) st%nik = 2*st%nik
2767
2768 safe_allocate(st%kweights(1:st%nik))
2769
2770 do iq = 1, st%nik
2771 ik = dd%get_kpoint_index(iq)
2772 st%kweights(iq) = kpoints%get_weight(ik)
2773 end do
2774
2775 if (debug%info) call print_kpoints_debug
2777
2778 contains
2779 subroutine print_kpoints_debug
2780 integer :: iunit
2781
2783
2784 call io_mkdir('debug/', namespace)
2785 iunit = io_open('debug/kpoints', namespace, action = 'write')
2786 call kpoints%write_info(iunit=iunit, absolute_coordinates = .true.)
2787 call io_close(iunit)
2788
2790 end subroutine print_kpoints_debug
2791
2792 end subroutine states_elec_choose_kpoints
2793
2798 !
2799 function states_elec_calculate_dipole(this, gr) result(dipole)
2800 class(states_elec_t), intent(in) :: this
2801 class(mesh_t), intent(in) :: gr
2802 real(real64) :: dipole(1:gr%box%dim)
2803
2804 integer :: ispin
2805 real(real64) :: e_dip(1:gr%box%dim, this%d%spin_channels)
2806
2808
2809 do ispin = 1, this%d%spin_channels
2810 call dmf_dipole(gr, this%rho(:, ispin), e_dip(:, ispin))
2811 end do
2813 dipole(:) = sum(e_dip(:, 1:this%d%spin_channels), 2) ! dipole moment <mu_el> = \sum_i -e <x_i>
2814
2816 end function states_elec_calculate_dipole
2817
2818
2819#include "undef.F90"
2820#include "real.F90"
2821#include "states_elec_inc.F90"
2822
2823#include "undef.F90"
2824#include "complex.F90"
2825#include "states_elec_inc.F90"
2826#include "undef.F90"
2827
2828end module states_elec_oct_m
2829
2830
2831!! Local Variables:
2832!! mode: f90
2833!! coding: utf-8
2834!! End:
initialize a batch with existing memory
Definition: batch.F90:277
constant times a vector plus a vector
Definition: lalg_basic.F90:173
Prints out to iunit a message in the form: ["InputVariable" = value] where "InputVariable" is given b...
Definition: messages.F90:182
pure logical function, public accel_is_enabled()
Definition: accel.F90:392
type(accel_t), public accel
Definition: accel.F90:248
This module implements batches of mesh functions.
Definition: batch.F90:135
This module implements common operations on batches of mesh functions.
Definition: batch_ops.F90:118
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.
Definition: boundaries.F90:124
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
Definition: global.F90:193
real(real64), parameter, public m_zero
Definition: global.F90:191
real(real64), parameter, public m_fourth
Definition: global.F90:200
complex(real64), parameter, public m_zi
Definition: global.F90:205
real(real64), parameter, public m_epsilon
Definition: global.F90:207
type(conf_t), public conf
Global instance of Octopus configuration.
Definition: global.F90:181
real(real64), parameter, public m_half
Definition: global.F90:197
real(real64), parameter, public m_one
Definition: global.F90:192
This module implements the underlying real-space grid.
Definition: grid.F90:119
subroutine, public dgrid_symmetrize_vector_field(gr, field, suppress_warning)
Definition: grid.F90:698
subroutine, public dgrid_symmetrize_scalar_field(gr, field, suppress_warning)
Definition: grid.F90:672
Definition: io.F90:116
logical pure function, public kpoints_point_is_gamma(this, ik)
Definition: kpoints.F90:1593
System information (time, memory, sysname)
Definition: loct.F90:117
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:117
integer pure function, public pad_pow2(size)
create array size, which is padded to powers of 2
Definition: math.F90:892
This module defines functions over batches of mesh functions.
Definition: mesh_batch.F90:118
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.
Definition: mesh.F90:120
subroutine, public messages_not_implemented(feature, namespace)
Definition: messages.F90:1091
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:525
subroutine, public messages_obsolete_variable(namespace, name, rep)
Definition: messages.F90:1023
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:162
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:410
subroutine, public messages_input_error(namespace, var, details, row, column)
Definition: messages.F90:691
subroutine, public messages_experimental(name, namespace)
Definition: messages.F90:1063
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:594
general module for modelmb particles
subroutine, public modelmb_particles_end(this)
subroutine, public modelmb_particles_init(this, namespace, space, nst)
==============================================================
subroutine, public modelmb_particles_copy(modelmb_out, modelmb_in)
subroutine mpi_grp_copy(mpi_grp_out, mpi_grp_in)
MPI_THREAD_FUNNELED allows for calls to MPI from an OMP region if the thread is the team master.
Definition: mpi.F90:389
type(mpi_comm), parameter, public mpi_comm_undefined
used to indicate a communicator has not been initialized
Definition: mpi.F90:138
subroutine mpi_grp_init(grp, comm)
Initialize MPI group instance.
Definition: mpi.F90:347
This module handles the communicators for the various parallelization strategies.
Definition: multicomm.F90:147
logical pure function, public multicomm_strategy_is_parallel(mc, level)
Definition: multicomm.F90:719
subroutine, public multicomm_all_pairs_copy(apout, apin)
Definition: multicomm.F90:248
logical pure function, public multicomm_have_slaves(this)
Definition: multicomm.F90:845
logical function, public parse_is_defined(namespace, name)
Definition: parser.F90:455
integer function, public parse_block(namespace, name, blk, check_varinfo_)
Definition: parser.F90:615
subroutine, public profiling_out(label)
Increment out counter and sum up difference between entry and exit time.
Definition: profiling.F90:631
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
Definition: profiling.F90:554
subroutine, public smear_fill_occupations(this, eigenvalues, occupations, nik, nst)
Definition: smear.F90:534
subroutine, public smear_find_fermi_energy(this, namespace, eigenvalues, occupations, qtot, nik, nst, kweights)
Definition: smear.F90:346
subroutine, public smear_copy(to, from)
Definition: smear.F90:325
integer, parameter, public smear_fixed_occ
Definition: smear.F90:173
subroutine, public smear_init(this, namespace, ispin, fixed_occ, integral_occs, kpoints)
Definition: smear.F90:187
logical pure function, public smear_is_semiconducting(this)
Definition: smear.F90:838
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 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_init(st, namespace, space, valence_charge, kpoints, calc_mode_id)
Initialize a new 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)
Explicitly 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
Definition: types.F90:135
type(type_t), public type_cmplx
Definition: types.F90:136
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:134
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
Definition: wfs_elec.F90:547
subroutine, public dwfs_elec_init(this, dim, st_start, st_end, np, ik, special, packed)
initialize an empty wfs_elec_t object
Definition: wfs_elec.F90:403
subroutine reduce_all(grp)
subroutine print_kpoints_debug
Class defining batches of mesh functions.
Definition: batch.F90:161
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.
Definition: grid.F90:171
Describes mesh distribution to nodes.
Definition: mesh.F90:187
This is defined even when running serial.
Definition: mpi.F90:144
An all-pairs communication schedule for a given group.
Definition: multicomm.F90:237
Stores all communicators and groups.
Definition: multicomm.F90:208
abstract class for states
The states_elec_t class contains all electronic wave functions.
int true(void)