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