Octopus
eigensolver.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
19#include "global.h"
20
22 use batch_oct_m
24 use debug_oct_m
32 use global_oct_m
33 use grid_oct_m
37 use, intrinsic :: iso_fortran_env
41 use loct_oct_m
42 use mesh_oct_m
46 use mpi_oct_m
50 use parser_oct_m
54 use smear_oct_m
55 use space_oct_m
62 use types_oct_m
63 use unit_oct_m
66 use xc_oct_m
67
68 implicit none
69
70 private
71 public :: &
75
76 type eigensolver_t
77 private
78 integer, public :: es_type
79
80 real(real64), public :: tolerance
81 integer, public :: es_maxiter
82
83 real(real64) :: tau
84 real(real64) :: tau0
85 real(real64) :: time
86 integer :: it_propagator
87 type(potential_interpolation_t) :: vks_old
88 real(real64), allocatable :: normalization_energies(:, :)
89 ! !! evolution eigensolver
90 real(real64), allocatable :: normalization_energies_prev(:, :)
92 logical :: variable_timestep
93
95 real(real64), allocatable, public :: diff(:, :)
96 integer, public :: matvec
97 integer, allocatable, public :: converged(:)
98
100 type(preconditioner_t), public :: pre
101
103 type(subspace_t) :: sdiag
104
105 integer :: rmmdiis_minimization_iter
106
107 logical, public :: folded_spectrum
108
109 ! cg options
110 logical, public :: orthogonalize_to_all
111 integer, public :: conjugate_direction
112 real(real64), public :: energy_change_threshold
113
114 ! Chebyshev filtering options
115 type(eigen_chebyshev_t), public :: cheby_params
117 type(exponential_t) :: exponential_operator
118 contains
119 procedure :: run => eigensolver_run
120 procedure :: set_lower_bound_is_known => eigensolver_set_lower_bound_is_known
121 end type eigensolver_t
122
123
124 integer, public, parameter :: &
125 RS_CG = 5, &
126 rs_evo = 9, &
127 rs_rmmdiis = 10, &
128 rs_chebyshev = 12
129
130contains
131
132 ! ---------------------------------------------------------
133 subroutine eigensolver_init(eigens, namespace, gr, st, hm, mc, space, deactivate_oracle)
134 type(eigensolver_t), intent(out) :: eigens
135 type(namespace_t), intent(in) :: namespace
136 type(grid_t), intent(in) :: gr
137 type(states_elec_t), intent(in) :: st
138 type(hamiltonian_elec_t), intent(in) :: hm
139 type(multicomm_t), intent(in) :: mc
140 class(space_t), intent(in) :: space
141 logical, optional, intent(in) :: deactivate_oracle
142
143 integer :: default_iter
144 real(real64) :: default_tol
145 real(real64) :: mem
146
147 push_sub(eigensolver_init)
148
149 !%Variable Eigensolver
150 !%Type integer
151 !%Section SCF::Eigensolver
152 !%Description
153 !% Which eigensolver to use to obtain the lowest eigenvalues and
154 !% eigenfunctions of the Kohn-Sham Hamiltonian. The default is
155 !% conjugate gradients (<tt>cg</tt>), except that when parallelization in states is
156 !% enabled, the default is <tt>chebyshev_filter</tt>.
157 !%Option cg 5
158 !% Conjugate-gradients algorithm.
159 !%Option evolution 9
160 !% (Experimental) Propagation in imaginary time.
161 !% This eigensolver uses a time propagation in imaginary time (see, e.g., Aichinger, M. & Krotscheck, E.,
162 !% Computational Materials Science 34, 188–212 (2005), 10.1016/j.commatsci.2004.11.002).
163 !% The timestep is set using <tt>EigensolverImaginaryTime</tt>. The convergence typically depends
164 !% on the total imaginary time in a propagation which means that it should be faster with larger
165 !% timesteps. However, simulations typically diverge above a certain timestep that depends on the system.
166 !% In general, it is advisable to use an exponential method that always converges, especially for large
167 !% timesteps, i.e., set <tt>TDExponentialMethod</tt> to <tt>lanczos</tt> or <tt>chebyshev</tt>.
168 !% The propagation method can be set with <tt>ImaginaryTimePropagator</tt> and variable timestepping
169 !% can be enabled with <tt>ImaginaryTimeVariableTimestep</tt>.
170 !% Moreover, this method usually converges quite slowly, so you might need to increase <tt>MaximumIter</tt>.
171 !% It is incompatible with smearing for the occupations.
172 !%Option rmmdiis 10
173 !% Residual minimization scheme, direct inversion in the
174 !% iterative subspace eigensolver, based on the implementation of
175 !% Kresse and Furthm&uuml;ller [<i>Phys. Rev. B</i> <b>54</b>, 11169
176 !% (1996)]. This eigensolver requires almost no orthogonalization
177 !% so it can be considerably faster than the other options for
178 !% large systems. To improve its performance a large number of <tt>ExtraStates</tt>
179 !% are required (around 10-20% of the number of occupied states).
180 !% Note: with <tt>unocc</tt>, you will need to stop the calculation
181 !% by hand, since the highest states will probably never converge.
182 !% Usage with more than one block of states per node is experimental, unfortunately.
183 !%Option chebyshev_filter 12
184 !% Chebyshev polynomials are used to construct a spectral filter that is iteratively applied to a trial subspace,
185 !% amplifying a subset of the lowest eigenvalues of the Hamiltonian, which effectively isolates the invariant
186 !% subspace associated with the these states. The filtered subspace is projected to define a reduced eigenvalue problem,
187 !% which is then solved using dense diagonalisation. Finally, the eigenstates of the original Hamiltonian are recovered
188 !% via subspace rotation. For further details, see [Zhou et. al.](http://dx.doi.org/10.1016/j.jcp.2014.06.056)
189 !%End
190 call parse_variable(namespace, 'Eigensolver', rs_chebyshev, eigens%es_type)
192 if(st%parallel_in_states .and. .not. eigensolver_parallel_in_states(eigens)) then
193 message(1) = "The selected eigensolver is not parallel in states."
194 message(2) = "Please use the rmmdiis or Chebyshev filtering eigensolvers."
195 call messages_fatal(2, namespace=namespace)
196 end if
197
198 call messages_obsolete_variable(namespace, 'EigensolverVerbose')
199 call messages_obsolete_variable(namespace, 'EigensolverSubspaceDiag', 'SubspaceDiagonalization')
201 default_iter = 25
203 select case(eigens%es_type)
204 case(rs_cg)
205 !%Variable CGOrthogonalizeAll
206 !%Type logical
207 !%Default no
208 !%Section SCF::Eigensolver
209 !%Description
210 !% Used by the cg solver only.
211 !% During the cg iterations, the current band can be orthogonalized
212 !% against all other bands or only against the lower bands. Orthogonalizing
213 !% against all other bands can improve convergence properties, whereas
214 !% orthogonalizing against lower bands needs less operations.
215 !% Moreover, orthogonalizing against all bands can make converging
216 !% the highest band or unoccupied bands more difficult.
217 !%End
218 call parse_variable(namespace, 'CGOrthogonalizeAll', .false., eigens%orthogonalize_to_all)
220 !%Variable CGDirection
221 !%Type integer
222 !%Section SCF::Eigensolver
223 !%Description
224 !% Used by the cg solver only.
225 !% The conjugate direction is updated using a certain coefficient to the previous
226 !% direction. This coeffiction can be computed in different ways. The default is
227 !% to use Fletcher-Reeves (FR), an alternative is Polak-Ribiere (PR).
228 !%Option fletcher 1
229 !% The coefficient for Fletcher-Reeves consists of the current norm of the
230 !% steepest descent vector divided by that of the previous iteration.
231 !%Option polak 2
232 !% For the Polak-Ribiere scheme, a product of the current with the previous
233 !% steepest descent vector is subtracted in the nominator.
234 !%End
235 call parse_variable(namespace, 'CGDirection', option__cgdirection__fletcher, eigens%conjugate_direction)
236
237 !%Variable CGEnergyChangeThreshold
238 !%Type float
239 !%Section SCF::Eigensolver
240 !%Default 0.1
241 !%Description
242 !% Used by the cg solver only.
243 !% For each band, the CG iterations are stopped when the change in energy is smaller than the
244 !% change in the first iteration multiplied by this factor. This limits the number of CG
245 !% iterations for each band, while still showing good convergence for the SCF cycle. The criterion
246 !% is discussed in Sec. V.B.6 of Payne et al. (1992), Rev. Mod. Phys. 64, 4.
247 !% The default value is 0.1, which is usually a good choice for LDA and GGA potentials. If you
248 !% are solving the OEP equation, you might want to set this value to 1e-3 or smaller. In general,
249 !% smaller values might help if you experience convergence problems.
250 !% For very small convergence tolerances, choose 0 to disable this criterion.
251 !%End
252 call parse_variable(namespace, 'CGEnergyChangeThreshold', 0.1_real64, eigens%energy_change_threshold)
253
254 case(rs_evo)
255 call messages_experimental("imaginary-time evolution eigensolver")
256
257 !%Variable EigensolverImaginaryTime
258 !%Type float
259 !%Default 0.1
260 !%Section SCF::Eigensolver
261 !%Description
262 !% The imaginary-time step that is used in the imaginary-time evolution
263 !% method (<tt>Eigensolver = evolution</tt>) to obtain the lowest eigenvalues/eigenvectors.
264 !% It must satisfy <tt>EigensolverImaginaryTime > 0</tt>.
265 !% Increasing this value can make the propagation faster, but could lead to unstable propagations.
266 !%End
267 call parse_variable(namespace, 'EigensolverImaginaryTime', 0.1_real64, eigens%tau)
268 if(eigens%tau <= m_zero) call messages_input_error(namespace, 'EigensolverImaginaryTime')
269 ! save initial timestep
270 eigens%tau0 = eigens%tau
271
272 call exponential_init(eigens%exponential_operator, namespace)
273
274 if(st%smear%method /= smear_semiconductor .and. st%smear%method /= smear_fixed_occ) then
275 message(1) = "Smearing of occupations is incompatible with imaginary time evolution."
276 call messages_fatal(1, namespace=namespace)
277 end if
278
279 !%Variable ImaginaryTimePropagator
280 !%Type integer
281 !%Section SCF::Eigensolver
282 !%Default it_forward_euler
283 !%Description
284 !% The propagation method for imaginary-time evolution.
285 !%Option it_forward_euler 1
286 !% Approximate the evolution operator by one forward Euler step with the exponential:
287 !% <math>
288 !% U_{\rm FE}(t+\delta t, t) = \exp \left( -\tau H_{t}\right)\,.
289 !% </math>
290 !% This scheme is first order.
291 !%Option it_expmid 2
292 !% Similar to forward Euler, but evaluate the Hamiltonian at the half-timestep by
293 !% extrapolation:
294 !% <math>
295 !% U_{\rm EM}(t+\delta t, t) = \exp \left( -\tau H_{t+\delta t/2}\right)\,.
296 !% </math>
297 !% This scheme is second order.
298 !%End
299 call parse_variable(namespace, 'ImaginaryTimePropagator', it_forward_euler, eigens%it_propagator)
300 ! initialize potential interpolation
301 if (eigens%it_propagator == it_expmid) then
302 call hm%ks_pot%init_interpolation(eigens%vks_old)
303 call hm%ks_pot%run_zero_iter(eigens%vks_old)
304 end if
305
306 !%Variable ImaginaryTimeVariableTimestep
307 !%Type logical
308 !%Default no
309 !%Section SCF::Eigensolver
310 !%Description
311 !% Enable variable timesteps for the imaginary time propagation.
312 !%End
313 call parse_variable(namespace, 'ImaginaryTimeVariableTimestep', .false., eigens%variable_timestep)
314
315 eigens%time = m_zero
316
317 safe_allocate(eigens%normalization_energies(1:st%nst, 1:st%nik))
318 safe_allocate(eigens%normalization_energies_prev(1:st%nst, 1:st%nik))
319
320 case(rs_rmmdiis)
321 default_iter = 5
322
323 !%Variable EigensolverMinimizationIter
324 !%Type integer
325 !%Default 0
326 !%Section SCF::Eigensolver
327 !%Description
328 !% During the first iterations, the RMMDIIS eigensolver requires
329 !% some steepest-descent minimizations to improve
330 !% convergence. This variable determines the number of those
331 !% minimizations.
332 !%End
333
334 call parse_variable(namespace, 'EigensolverMinimizationIter', 0, eigens%rmmdiis_minimization_iter)
335
336 if(gr%use_curvilinear) call messages_experimental("RMMDIIS eigensolver for curvilinear coordinates")
337
338 case(rs_chebyshev)
339 !%Variable ChebyshevFilterLanczosOrder
340 !%Type integer
341 !%Default -1
342 !%Section SCF::Eigensolver
343 !%Description
344 !% Used by the Chebyshev filter only.
345 !% The number of Lanczos iterations used to construct the tridiagonal matrix,
346 !% from which the upper bound of H is estimated.
347 !% A value in the range 4 <= ChebyshevFilterLanczosOrder <= 10 is reasonable.
348 !% Values greater than 30 will raise an assertion.
349 !%
350 !% The default value, -1, indicate that the order is selected automatically.
351 !% In this case, the trial vector is also reused to make the procedure faster.
352 !%End
353 call parse_variable(namespace, 'ChebyshevFilterLanczosOrder', default_chebyshev_params%n_lanczos, &
354 eigens%cheby_params%n_lanczos)
355
356 ! TODO Rewrite this, having changed the behaviour
357 !%Variable ChebyshevFilterDegree
358 !%Type integer
359 !%Default 25
360 !%Section SCF::Eigensolver
361 !%Description
362 !% Used by the Chebyshev filter only.
363 !% The degree of the Chebyshev polynomial used to dampen the interval of eigenstates one is not interested in.
364 !% If used in conjunction with "OptimizeChebyshevFilterDegree", which is the default, "ChebyshevFilterDegree" defines the
365 !% the maximum Chebyshev polynomial degree Octopus will consider when determining an approximate, optimal degree.
366 !%
367 !% If not used in conjunction with "OptimizeChebyshevFilterDegree", this defines the polynomial degree used for all SCF steps.
368 !% The larger the degree, the larger the separation between the subspaces, which will reduce the number of SCF iterations
369 !% required to reach convergence. However, ChebyshevFilterDegree also directly correlates with the number of matrix-vector
370 !% products performed on the Hamiltonian per SCF step, and so larger values result in slower SCF cycles.
371 !% A value in the range 8 <= ChebyshevFilterDegree <= 20 is reasonable.
372 !%End
373 call parse_variable(namespace, 'ChebyshevFilterDegree', default_chebyshev_params%degree, eigens%cheby_params%degree)
374
375 !%Variable OptimizeChebyshevFilterDegree
376 !%Type logical
377 !%Default yes
378 !%Section SCF::Eigensolver
379 !%Description
380 !% Used by the Chebyshev filter only.
381 !% Octopus generates a best-estimate for the degree of the Chebyshev polynomial used to filter the subspace.
382 !% A separate estimate is generated for each state block, per SCF iteration. Note that if running parallelism
383 !% over states, the block/batch containing the largest eigenstates will converge more slowly and will
384 !% therefore use a larger degree relative to all other batches. One should therefore avoid setting "ChebyshevFilterDegree"
385 !% to an excessive value (for example > 50). For more details regarding how the degree is estimated, one can
386 !% refer to Section 5.4 in "Computer Physics Communications 187 (2015) 98–105"
387 !% (http://dx.doi.org/10.1016/j.cpc.2014.10.015).
388 !%
389 !% This is deactivated by default for spinors as this is not found to be advantageous in test systems.
390 !%End
391 if (st%d%ispin == spinors .or. optional_default(deactivate_oracle, .false.)) then
392 default_chebyshev_params%optimize_degree = .false.
393 end if
394 call parse_variable(namespace, 'OptimizeChebyshevFilterDegree', default_chebyshev_params%optimize_degree, &
395 eigens%cheby_params%optimize_degree)
396
397 !%Variable ChebyshevFilterBoundMixing
398 !%Type float
399 !%Default 0.5
400 !%Section SCF::Eigensolver
401 !%Description
402 !% In the first application of the filter, for the first SCF step, the initial lower bound estimate
403 !% is defined as a linear combination of the smallest and largest eigenvalues of the Hamiltonian.
404 !% The bound mixing determines the proportion of linear mixing, beta:
405 !% $b_{lower} = beta min{\lambda} + (\beta - 1) max{\lambda}$
406 !% of the smallest and largest eigenvalues.
407 !%End
408 call parse_variable(namespace, 'ChebyshevFilterBoundMixing', default_chebyshev_params%bound_mixing, &
409 eigens%cheby_params%bound_mixing)
410
411 !%Variable ChebyshevFilterNIter
412 !%Type integer
413 !%Default 5
414 !%Section SCF::Eigensolver
415 !%Description
416 !% The max number of iterations in the first SCF step of the Chebyshev solver. In practice this
417 !% does not need to exceed 10, as the eigenstates will vary alot between the first and second
418 !% SCF steps.
419 !%End
420 call parse_variable(namespace, 'ChebyshevFilterNIter', default_chebyshev_params%n_iter, &
421 eigens%cheby_params%n_iter)
422
423 case default
424 call messages_input_error(namespace, 'Eigensolver')
425 end select
426
427 call messages_print_with_emphasis(msg='Eigensolver', namespace=namespace)
428
429 call messages_print_var_option("Eigensolver", eigens%es_type, namespace=namespace)
430
431 call messages_obsolete_variable(namespace, 'EigensolverInitTolerance', 'EigensolverTolerance')
432 call messages_obsolete_variable(namespace, 'EigensolverFinalTolerance', 'EigensolverTolerance')
433 call messages_obsolete_variable(namespace, 'EigensolverFinalToleranceIteration')
434
435 ! this is an internal option that makes the solver use the
436 ! folded operator (H-shift)^2 to converge first eigenvalues around
437 ! the values of shift
438 ! c.f. L. W. Wang and A. Zunger
439 ! JCP 100, 2394 (1994); doi: http://dx.doi.org/10.1063/1.466486
440 eigens%folded_spectrum = .false.
441
442 !%Variable EigensolverTolerance
443 !%Type float
444 !%Section SCF::Eigensolver
445 !%Description
446 !% This is the tolerance for the eigenvectors. The default is 1e-7.
447 !% However, if <tt>ConvRelDens</tt> is defined, the default is set to a tenth of its value.
448 !%End
449 default_tol = 1e-7_real64
450 if (parse_is_defined(namespace, 'ConvRelDens')) then
451 call parse_variable(namespace, 'ConvRelDens', 1e-6_real64, default_tol)
452 default_tol = default_tol / 10.0_real64
453 end if
454 call parse_variable(namespace, 'EigensolverTolerance', default_tol, eigens%tolerance)
455
456 !%Variable EigensolverMaxIter
457 !%Type integer
458 !%Section SCF::Eigensolver
459 !%Description
460 !% Determines the maximum number of iterations that the
461 !% eigensolver will perform if the desired tolerance is not
462 !% achieved. The default is 25 iterations for all eigensolvers
463 !% except for <tt>rmdiis</tt>, which performs only 5 iterations.
464 !% Increasing this value for <tt>rmdiis</tt> increases the convergence speed,
465 !% at the cost of an increased memory footprint.
466 !%End
467 call parse_variable(namespace, 'EigensolverMaxIter', default_iter, eigens%es_maxiter)
468 if(eigens%es_maxiter < 1) call messages_input_error(namespace, 'EigensolverMaxIter')
469
470 if(eigens%es_maxiter > default_iter) then
471 call messages_write('You have specified a large number of eigensolver iterations (')
472 call messages_write(eigens%es_maxiter)
473 call messages_write(').', new_line = .true.)
474 call messages_write('This is not a good idea as it might slow down convergence, even for', new_line = .true.)
475 call messages_write('independent particles, as subspace diagonalization will not be used', new_line = .true.)
476 call messages_write('often enough.')
477 call messages_warning(namespace=namespace)
478 end if
479
480 if (any(eigens%es_type == (/rs_cg, rs_rmmdiis/))) then
481 call preconditioner_init(eigens%pre, namespace, gr, mc, space)
482 end if
483
484 safe_allocate(eigens%diff(1:st%nst, 1:st%nik))
485 eigens%diff(1:st%nst, 1:st%nik) = 0
486
487 safe_allocate(eigens%converged(1:st%nik))
488 eigens%converged(1:st%nik) = 0
489 eigens%matvec = 0
490
491 call eigens%sdiag%init(namespace, st)
492
493 ! print memory requirements
494 select case(eigens%es_type)
495 case(rs_rmmdiis)
496 call messages_write('Info: The rmmdiis eigensolver requires ')
497 mem = (m_two*eigens%es_maxiter - m_one)*st%block_size*real(gr%np_part, real64)
498 if(states_are_real(st)) then
499 mem = mem*8.0_real64
500 else
501 mem = mem*16.0_real64
502 end if
503 call messages_write(mem, units = unit_megabytes, fmt = '(f9.1)')
504 call messages_write(' of additional')
505 call messages_new_line()
506 call messages_write(' memory. This amount can be reduced by decreasing the value')
507 call messages_new_line()
508 call messages_write(' of the variable StatesBlockSize (currently set to ')
509 call messages_write(st%block_size)
510 call messages_write(').')
511 call messages_info()
512 end select
513
514 call messages_print_with_emphasis(namespace=namespace)
515
516 pop_sub(eigensolver_init)
517 end subroutine eigensolver_init
518
519
520 ! ---------------------------------------------------------
521 subroutine eigensolver_end(eigens)
522 type(eigensolver_t), intent(inout) :: eigens
523
524 push_sub(eigensolver_end)
525
526 select case(eigens%es_type)
527 case(rs_cg, rs_rmmdiis)
528 call preconditioner_end(eigens%pre)
529 end select
530
531 safe_deallocate_a(eigens%converged)
532 safe_deallocate_a(eigens%diff)
533
534 safe_deallocate_a(eigens%normalization_energies)
535 safe_deallocate_a(eigens%normalization_energies_prev)
536
537 pop_sub(eigensolver_end)
538 end subroutine eigensolver_end
539
540
541 ! ---------------------------------------------------------
542 subroutine eigensolver_run(eigens, namespace, gr, st, hm, space, ext_partners, iter, conv, nstconv)
543 class(eigensolver_t), intent(inout) :: eigens
544 type(namespace_t), intent(in) :: namespace
545 type(grid_t), intent(in) :: gr
546 type(states_elec_t), intent(inout) :: st
547 type(hamiltonian_elec_t), intent(inout) :: hm
548 class(space_t), intent(in) :: space
549 type(partner_list_t), intent(in) :: ext_partners
550 integer, intent(in) :: iter
551 logical, optional, intent(out) :: conv
552 integer, optional, intent(in) :: nstconv
553 ! !< the convergence criteria
554
555 integer :: nstconv_
556#ifdef HAVE_MPI
557 logical :: conv_reduced
558 integer :: ist, outcount, lmatvec
559 real(real64), allocatable :: ldiff(:), leigenval(:)
560 real(real64), allocatable :: ldiff_out(:), leigenval_out(:)
561 integer, allocatable :: lconv(:)
562#endif
563
564 call profiling_in("EIGEN_SOLVER")
565 push_sub(eigensolver_run)
566
567 if(present(conv)) conv = .false.
568 if(present(nstconv)) then
569 nstconv_ = nstconv
570 else
571 nstconv_ = st%nst
572 end if
573
574 eigens%matvec = 0
575
576 if(mpi_world%is_root() .and. eigensolver_has_progress_bar(eigens) .and. .not. debug%info) then
577 call loct_progress_bar(-1, st%lnst*st%d%kpt%nlocal)
578 end if
579
580 if (states_are_real(st)) then
581 call deigensolver_run(eigens, namespace, gr, st, hm, space, ext_partners, iter)
582 else
583 call zeigensolver_run(eigens, namespace, gr, st, hm, space, ext_partners, iter)
584 end if
585
586 if(mpi_world%is_root() .and. eigensolver_has_progress_bar(eigens) .and. .not. debug%info) then
587 write(stdout, '(1x)')
588 end if
589
590 if(present(conv)) conv = all(eigens%converged(st%d%kpt%start:st%d%kpt%end) >= nstconv_)
591
592#ifdef HAVE_MPI
593 if (st%d%kpt%parallel) then
594 if (present(conv)) then
595 call st%d%kpt%mpi_grp%allreduce(conv, conv_reduced, 1, mpi_logical, mpi_land)
596 conv = conv_reduced
597 end if
598
599 lmatvec = eigens%matvec
600 call st%d%kpt%mpi_grp%allreduce(lmatvec, eigens%matvec, 1, mpi_integer, mpi_sum)
601
602 safe_allocate(lconv(1:st%d%kpt%nlocal))
603 lconv(1:st%d%kpt%nlocal) = eigens%converged(st%d%kpt%start:st%d%kpt%end)
604 call lmpi_gen_allgatherv(st%d%kpt%nlocal, lconv, outcount, eigens%converged, st%d%kpt%mpi_grp)
605 assert(outcount == st%nik)
606 safe_deallocate_a(lconv)
607
608 ! every node needs to know all eigenvalues (and diff)
609 safe_allocate(ldiff(1:st%d%kpt%nlocal))
610 safe_allocate(leigenval(1:st%d%kpt%nlocal))
611 safe_allocate(ldiff_out(1:st%nik))
612 safe_allocate(leigenval_out(1:st%nik))
613 do ist = st%st_start, st%st_end
614 ldiff(1:st%d%kpt%nlocal) = eigens%diff(ist, st%d%kpt%start:st%d%kpt%end)
615 leigenval(1:st%d%kpt%nlocal) = st%eigenval(ist, st%d%kpt%start:st%d%kpt%end)
616 call lmpi_gen_allgatherv(st%d%kpt%nlocal, ldiff, outcount, ldiff_out, st%d%kpt%mpi_grp)
617 eigens%diff(ist, :) = ldiff_out
618 assert(outcount == st%nik)
619 call lmpi_gen_allgatherv(st%d%kpt%nlocal, leigenval, outcount, leigenval_out, st%d%kpt%mpi_grp)
620 st%eigenval(ist, :) = leigenval_out
621 assert(outcount == st%nik)
622 end do
623 safe_deallocate_a(ldiff)
624 safe_deallocate_a(ldiff_out)
625 safe_deallocate_a(leigenval)
626 safe_deallocate_a(leigenval_out)
627 end if
628#endif
629
630 pop_sub(eigensolver_run)
631 call profiling_out("EIGEN_SOLVER")
632 end subroutine eigensolver_run
633
634
635 ! ---------------------------------------------------------
636 logical function eigensolver_parallel_in_states(this) result(par_stat)
637 type(eigensolver_t), intent(in) :: this
638
640
641 par_stat = .false.
642
643 select case(this%es_type)
645 par_stat = .true.
646 end select
647
650
651
652 ! ---------------------------------------------------------
653 logical function eigensolver_has_progress_bar(this) result(has)
654 type(eigensolver_t), intent(in) :: this
655
657
658 has = .false.
659
660 select case(this%es_type)
661 case(rs_rmmdiis, rs_cg)
662 has = .true.
663 end select
664
667
669 pure subroutine eigensolver_set_lower_bound_is_known(this, known_lower_bound)
670 class(eigensolver_t), intent(inout) :: this
671 logical, intent(in) :: known_lower_bound
672
673 this%cheby_params%known_lower_bound = known_lower_bound
675
676
677#include "undef.F90"
678#include "real.F90"
679#include "eigensolver_inc.F90"
680
681#include "undef.F90"
682#include "complex.F90"
683#include "eigensolver_inc.F90"
684
685end module eigensolver_oct_m
687
688!! Local Variables:
689!! mode: f90
690!! coding: utf-8
691!! End:
This module implements batches of mesh functions.
Definition: batch.F90:135
This module implements common operations on batches of mesh functions.
Definition: batch_ops.F90:118
type(debug_t), save, public debug
Definition: debug.F90:158
This module calculates the derivatives (gradients, Laplacians, etc.) of a function.
type(eigen_chebyshev_t), public default_chebyshev_params
Default Chebyshev input parameters Arguments 1 and 2 taken from 10.1016/j.jcp.2006....
integer, parameter, public it_expmid
integer, parameter, public it_forward_euler
subroutine zeigensolver_run(eigens, namespace, mesh, st, hm, space, ext_partners, iter)
integer, parameter, public rs_evo
subroutine, public eigensolver_init(eigens, namespace, gr, st, hm, mc, space, deactivate_oracle)
subroutine deigensolver_run(eigens, namespace, mesh, st, hm, space, ext_partners, iter)
subroutine eigensolver_run(eigens, namespace, gr, st, hm, space, ext_partners, iter, conv, nstconv)
logical function eigensolver_parallel_in_states(this)
integer, parameter, public rs_chebyshev
logical function eigensolver_has_progress_bar(this)
subroutine, public eigensolver_end(eigens)
integer, parameter, public rs_rmmdiis
pure subroutine eigensolver_set_lower_bound_is_known(this, known_lower_bound)
Set the flag lower_bound_is_known.
integer, parameter, public spinors
subroutine, public exponential_init(te, namespace, full_batch)
real(real64), parameter, public m_two
Definition: global.F90:193
real(real64), parameter, public m_zero
Definition: global.F90:191
real(real64), parameter, public m_one
Definition: global.F90:192
This module implements the underlying real-space grid.
Definition: grid.F90:119
This module defines classes and functions for interaction partners.
A module to handle KS potential, without the external potential.
System information (time, memory, sysname)
Definition: loct.F90:117
subroutine, public loct_progress_bar(a, maxcount)
A wrapper around the progress bar, such that it can be silenced without needing to dress the call wit...
Definition: loct.F90:276
This module defines functions over batches of mesh functions.
Definition: mesh_batch.F90:118
This module defines various routines, operating on mesh functions.
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:120
subroutine, public messages_print_with_emphasis(msg, iunit, namespace)
Definition: messages.F90:898
character(len=512), private msg
Definition: messages.F90:167
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:525
subroutine, public messages_obsolete_variable(namespace, name, rep)
Definition: messages.F90:1023
subroutine, public messages_new_line()
Definition: messages.F90:1112
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:162
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:410
subroutine, public messages_input_error(namespace, var, details, row, column)
Definition: messages.F90:691
subroutine, public messages_experimental(name, namespace)
Definition: messages.F90:1063
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:594
This module contains some common usage patterns of MPI routines.
Definition: mpi_lib.F90:117
type(mpi_grp_t), public mpi_world
Definition: mpi.F90:272
This module handles the communicators for the various parallelization strategies.
Definition: multicomm.F90:147
logical function, public parse_is_defined(namespace, name)
Definition: parser.F90:455
subroutine, public preconditioner_end(this)
subroutine, public preconditioner_init(this, namespace, gr, mc, space)
subroutine, public profiling_out(label)
Increment out counter and sum up difference between entry and exit time.
Definition: profiling.F90:631
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
Definition: profiling.F90:554
integer, parameter, public smear_semiconductor
Definition: smear.F90:173
integer, parameter, public smear_fixed_occ
Definition: smear.F90:173
pure logical function, public states_are_real(st)
This module handles spin dimensions of the states and the k-point distribution.
This module provides routines for communicating states when using states parallelization.
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:134
This module defines the unit system, used for input and output.
type(unit_t), public unit_megabytes
For large amounts of data (natural code units are bytes)
Definition: xc.F90:116
Description of the grid, containing information on derivatives, stencil, and symmetries.
Definition: grid.F90:171
The states_elec_t class contains all electronic wave functions.
int true(void)