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