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