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