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