Octopus
exponential.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-2006 M. Marques, A. Castro, A. Rubio, G. Bertsch
2!! Copyright (C) 2019 M. Oliveira
3!!
4!! This program is free software; you can redistribute it and/or modify
5!! it under the terms of the GNU General Public License as published by
6!! the Free Software Foundation; either version 2, or (at your option)
7!! any later version.
8!!
9!! This program is distributed in the hope that it will be useful,
10!! but WITHOUT ANY WARRANTY; without even the implied warranty of
11!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12!! GNU General Public License for more details.
13!!
14!! You should have received a copy of the GNU General Public License
15!! along with this program; if not, write to the Free Software
16!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
17!! 02110-1301, USA.
18!!
19
20#include "global.h"
21
23 use accel_oct_m
24 use batch_oct_m
26 use blas_oct_m
28 use debug_oct_m
29 use global_oct_m
33 use, intrinsic :: iso_fortran_env
37 use math_oct_m
38 use mesh_oct_m
43 use parser_oct_m
47 use types_oct_m
49 use xc_oct_m
50
51 implicit none
52
53 private
54 public :: &
60
61 integer, public, parameter :: &
62 EXP_LANCZOS = 2, &
63 exp_taylor = 3, &
65
66 type exponential_t
67 private
68 integer, public :: exp_method
69 real(real64) :: lanczos_tol
70 real(real64) :: chebyshev_tol
71 integer, public :: exp_order
72 integer :: arnoldi_gs
73 logical, public :: full_batch = .false.
74 contains
75 procedure :: apply_batch => exponential_apply_batch
76 procedure :: apply_single => exponential_apply_single
77 procedure :: apply_phi_batch => exponential_apply_phi_batch
78 end type exponential_t
79
80contains
81
82 ! ---------------------------------------------------------
83 subroutine exponential_init(te, namespace, full_batch)
84 type(exponential_t), intent(out) :: te
85 type(namespace_t), intent(in) :: namespace
86 logical, optional, intent(in) :: full_batch
87
88 push_sub(exponential_init)
89
90 !%Variable TDExponentialMethod
91 !%Type integer
92 !%Default taylor
93 !%Section Time-Dependent::Propagation
94 !%Description
95 !% Method used to numerically calculate the exponential of the Hamiltonian,
96 !% a core part of the full algorithm used to approximate the evolution
97 !% operator, specified through the variable <tt>TDPropagator</tt>.
98 !% In the case of using the Magnus method, described below, the action of the exponential
99 !% of the Magnus operator is also calculated through the algorithm specified
100 !% by this variable.
101 !%Option lanczos 2
102 !% Allows for larger time-steps.
103 !% However, the larger the time-step, the longer the computational time per time-step.
104 !% In certain cases, if the time-step is too large, the code will emit a warning
105 !% whenever it considers that the evolution may not be properly proceeding --
106 !% the Lanczos process did not converge. The method consists in a Krylov
107 !% subspace approximation of the action of the exponential
108 !% (see M. Hochbruck and C. Lubich, <i>SIAM J. Numer. Anal.</i> <b>34</b>, 1911 (1997) for details).
109 !% The performance of the method is controlled by the tolerance (controlled by <tt>TDLanczosTol</tt>).
110 !% The smaller the tolerance, the more precisely the exponential
111 !% is calculated, but also the larger the dimension of the Arnoldi
112 !% subspace. If the maximum dimension (currently 200) is not enough to meet the criterion,
113 !% the above-mentioned warning is emitted.
114 !% Be aware that the larger the required dimension of the Krylov subspace, the larger
115 !% the memory required for this method. So if you run out of memory, try to reduce
116 !% the time step.
117 !%Option taylor 3
118 !% This method amounts to a straightforward application of the definition of
119 !% the exponential of an operator, in terms of its Taylor expansion.
120 !%
121 !% <math>\exp_{\rm STD} (-i\delta t H) = \sum_{i=0}^{k} {(-i\delta t)^i\over{i!}} H^i.</math>
122 !%
123 !% The order <i>k</i> is determined by variable <tt>TDExpOrder</tt>.
124 !% Some numerical considerations from <a href=http://www.phys.washington.edu/~bertsch/num3.ps>
125 !% Jeff Giansiracusa and George F. Bertsch</a>
126 !% suggest the 4th order as especially suitable and stable.
127 !%Option chebyshev 4
128 !% In principle, the Chebyshev expansion
129 !% of the exponential represents it more accurately than the canonical or standard expansion.
130 !% <tt>TDChebyshevTol</tt> determines the tolerance to which the expansion is computed.
131 !%
132 !% There exists a closed analytic form for the coefficients of the exponential in terms
133 !% of Chebyshev polynomials:
134 !%
135 !% <math>\exp_{\rm CHEB} \left( -i\delta t H \right) = \sum_{k=0}^{\infty} (2-\delta_{k0})(-i)^{k}J_k(\delta t) T_k(H),</math>
136 !%
137 !% where <math>J_k</math> are the Bessel functions of the first kind, and H has to be previously
138 !% scaled to <math>[-1,1]</math>.
139 !% See H. Tal-Ezer and R. Kosloff, <i>J. Chem. Phys.</i> <b>81</b>,
140 !% 3967 (1984); R. Kosloff, <i>Annu. Rev. Phys. Chem.</i> <b>45</b>, 145 (1994);
141 !% C. W. Clenshaw, <i>MTAC</i> <b>9</b>, 118 (1955).
142 !%End
143 call parse_variable(namespace, 'TDExponentialMethod', exp_taylor, te%exp_method)
144
145 select case (te%exp_method)
146 case (exp_taylor)
147 case (exp_chebyshev)
148 !%Variable TDChebyshevTol
149 !%Type float
150 !%Default 1e-10
151 !%Section Time-Dependent::Propagation
152 !%Description
153 !% An internal tolerance variable for the Chebyshev method. The smaller, the more
154 !% precisely the exponential is calculated and the more iterations are needed, i.e.,
155 !% it becomes more expensive. The expansion is terminated once the error estimate
156 !% is below this tolerance.
157 !%End
158 call parse_variable(namespace, 'TDChebyshevTol', 1e-10_real64, te%chebyshev_tol)
159 if (te%chebyshev_tol <= m_zero) call messages_input_error(namespace, 'TDChebyshevTol')
160 case (exp_lanczos)
161 !%Variable TDLanczosTol
162 !%Type float
163 !%Default 1e-5
164 !%Section Time-Dependent::Propagation
165 !%Description
166 !% An internal tolerance variable for the Lanczos method. The smaller, the more
167 !% precisely the exponential is calculated, and also the bigger the dimension
168 !% of the Krylov subspace needed to perform the algorithm. One should carefully
169 !% make sure that this value is not too big, or else the evolution will be
170 !% wrong.
171 !%End
172 call parse_variable(namespace, 'TDLanczosTol', 1e-5_real64, te%lanczos_tol)
173 if (te%lanczos_tol <= m_zero) call messages_input_error(namespace, 'TDLanczosTol')
174
175 case default
176 call messages_input_error(namespace, 'TDExponentialMethod')
177 end select
178 call messages_print_var_option('TDExponentialMethod', te%exp_method, namespace=namespace)
179
180 if (te%exp_method == exp_taylor) then
181 !%Variable TDExpOrder
182 !%Type integer
183 !%Default 4
184 !%Section Time-Dependent::Propagation
185 !%Description
186 !% For <tt>TDExponentialMethod</tt> = <tt>taylor</tt>,
187 !% the order to which the exponential is expanded.
188 !%End
189 call parse_variable(namespace, 'TDExpOrder', default__tdexporder, te%exp_order)
190 if (te%exp_order < 2) call messages_input_error(namespace, 'TDExpOrder')
191
192 end if
193
194 te%arnoldi_gs = option__arnoldiorthogonalization__cgs
195 if (te%exp_method == exp_lanczos) then
196 !%Variable ArnoldiOrthogonalization
197 !%Type integer
198 !%Section Time-Dependent::Propagation
199 !%Description
200 !% The orthogonalization method used for the Arnoldi procedure.
201 !% Only for TDExponentialMethod = lanczos.
202 !%Option cgs 3
203 !% Classical Gram-Schmidt (CGS) orthogonalization.
204 !% The algorithm is defined in Giraud et al., Computers and Mathematics with Applications 50, 1069 (2005).
205 !%Option drcgs 5
206 !% Classical Gram-Schmidt orthogonalization with double-step reorthogonalization.
207 !% The algorithm is taken from Giraud et al., Computers and Mathematics with Applications 50, 1069 (2005).
208 !% According to this reference, this is much more precise than CGS or MGS algorithms.
209 !%End
210 call parse_variable(namespace, 'ArnoldiOrthogonalization', option__arnoldiorthogonalization__cgs, &
211 te%arnoldi_gs)
212 end if
213
214 ! do lanczos expansion for full batch?
215 te%full_batch = optional_default(full_batch, te%full_batch)
216
217 pop_sub(exponential_init)
218 end subroutine exponential_init
219
220 ! ---------------------------------------------------------
221 subroutine exponential_copy(teo, tei)
222 type(exponential_t), intent(inout) :: teo
223 type(exponential_t), intent(in) :: tei
224
225 push_sub(exponential_copy)
226
227 teo%exp_method = tei%exp_method
228 teo%lanczos_tol = tei%lanczos_tol
229 teo%exp_order = tei%exp_order
230 teo%arnoldi_gs = tei%arnoldi_gs
231
232 pop_sub(exponential_copy)
233 end subroutine exponential_copy
234
235 ! ---------------------------------------------------------
237 subroutine exponential_apply_single(te, namespace, mesh, hm, zpsi, ist, ik, deltat, vmagnus, imag_time)
238 class(exponential_t), intent(inout) :: te
239 type(namespace_t), intent(in) :: namespace
240 class(mesh_t), intent(in) :: mesh
241 type(hamiltonian_elec_t), intent(inout) :: hm
242 integer, intent(in) :: ist
243 integer, intent(in) :: ik
244 complex(real64), contiguous, intent(inout) :: zpsi(:, :)
245 real(real64), intent(in) :: deltat
246 real(real64), optional, intent(in) :: vmagnus(mesh%np, hm%d%nspin, 2)
247 logical, optional, intent(in) :: imag_time
248
249 type(wfs_elec_t) :: psib, inh_psib
250 complex(real64), allocatable :: zpsi_inh(:, :)
251
253
254 !We apply the phase only to np points, and the phase for the np+1 to np_part points
255 !will be treated as a phase correction in the Hamiltonian
256 if (hm%phase%is_allocated()) then
257 call hm%phase%apply_to_single(zpsi, mesh%np, hm%d%dim, ik, .false.)
258 end if
259
260 call wfs_elec_init(psib, hm%d%dim, ist, ist, zpsi, ik)
261
262 if (hamiltonian_elec_inh_term(hm)) then
263 safe_allocate(zpsi_inh(1:mesh%np_part, 1:hm%d%dim))
264 call states_elec_get_state(hm%inh_st, mesh, ist, ik, zpsi_inh(:, :))
265 call wfs_elec_init(inh_psib, hm%d%dim, ist, ist, zpsi_inh, ik)
266 call te%apply_batch(namespace, mesh, hm, psib, deltat, &
267 vmagnus=vmagnus, imag_time=imag_time, inh_psib=inh_psib)
268 call inh_psib%end()
269 safe_deallocate_a(zpsi_inh)
270 else
271 call te%apply_batch(namespace, mesh, hm, psib, deltat, &
272 vmagnus=vmagnus, imag_time=imag_time)
273 end if
274
275 call psib%end()
276
277 if (hm%phase%is_allocated()) then
278 call hm%phase%apply_to_single(zpsi, mesh%np, hm%d%dim, ik, .true.)
279 end if
280
282 end subroutine exponential_apply_single
283
284 ! ---------------------------------------------------------
301 ! ---------------------------------------------------------
302 subroutine exponential_apply_batch(te, namespace, mesh, hm, psib, deltat, psib2, deltat2, vmagnus, imag_time, inh_psib)
303 class(exponential_t), intent(inout) :: te
304 type(namespace_t), intent(in) :: namespace
305 class(mesh_t), intent(in) :: mesh
306 class(hamiltonian_abst_t), intent(inout) :: hm
307 class(batch_t), intent(inout) :: psib
308 real(real64), intent(in) :: deltat
309 class(batch_t), optional, intent(inout) :: psib2
310 real(real64), optional, intent(in) :: deltat2
311 real(real64), optional, intent(in) :: vmagnus(:,:,:) !(mesh%np, hm%d%nspin, 2)
312 logical, optional, intent(in) :: imag_time
313 class(batch_t), optional, intent(inout) :: inh_psib
315 complex(real64) :: deltat_, deltat2_
316 class(chebyshev_function_t), pointer :: chebyshev_function
317 logical :: imag_time_
318
320 call profiling_in("EXPONENTIAL_BATCH")
321
322 assert(psib%type() == type_cmplx)
323
324 assert(present(psib2) .eqv. present(deltat2))
325 if (present(inh_psib)) then
326 assert(inh_psib%nst == psib%nst)
327 end if
328
329 if (present(vmagnus)) then
330 select type(hm)
331 class is (hamiltonian_elec_t)
332 assert(size(vmagnus, 1) >= mesh%np)
333 assert(size(vmagnus, 2) == hm%d%nspin)
334 assert(size(vmagnus, 3) == 2)
335 class default
336 write(message(1), '(a)') 'Magnus operators only implemented for electrons at the moment'
337 call messages_fatal(1, namespace=namespace)
338 end select
339 end if
340
341 deltat2_ = cmplx(optional_default(deltat2, m_zero), m_zero, real64)
342
343 imag_time_ = optional_default(imag_time, .false.)
344 if (imag_time_) then
345 deltat_ = -m_zi*deltat
346 if (present(deltat2)) deltat2_ = m_zi*deltat2
347 else
348 deltat_ = cmplx(deltat, m_zero, real64)
349 if (present(deltat2)) deltat2_ = cmplx(deltat2, m_zero, real64)
350 end if
351
352 if (.not. hm%is_hermitian() .and. te%exp_method == exp_chebyshev) then
353 write(message(1), '(a)') 'The Chebyshev expansion cannot be used for non-Hermitian operators.'
354 write(message(2), '(a)') 'Please use the Lanczos exponentiation scheme ("TDExponentialMethod = lanczos")'
355 write(message(3), '(a)') 'or the Taylor expansion ("TDExponentialMethod = taylor") method.'
356 call messages_fatal(3, namespace=namespace)
357 end if
358
359 select case (te%exp_method)
360 case (exp_taylor)
361 ! Note that delttat2_ is only initialized if deltat2 is present.
362 if (present(deltat2)) then
363 call exponential_taylor_series_batch(te, namespace, mesh, hm, psib, deltat_, &
364 psib2, deltat2_, vmagnus)
365 else
366 call exponential_taylor_series_batch(te, namespace, mesh, hm, psib, deltat_, &
367 vmagnus=vmagnus)
368 end if
369 if (present(inh_psib)) then
370 if (present(deltat2)) then
371 call exponential_taylor_series_batch(te, namespace, mesh, hm, psib, deltat_, &
372 psib2, deltat2_, vmagnus, inh_psib)
373 else
374 call exponential_taylor_series_batch(te, namespace, mesh, hm, psib, deltat_, &
375 vmagnus=vmagnus, inh_psib=inh_psib)
376 end if
377 end if
378
379 case (exp_lanczos)
380 if (present(psib2)) call psib%copy_data_to(mesh%np, psib2)
381 call exponential_lanczos_batch(te, namespace, mesh, hm, psib, deltat_, vmagnus)
382 if (present(inh_psib)) then
383 call exponential_lanczos_batch(te, namespace, mesh, hm, psib, deltat_, vmagnus, inh_psib)
384 end if
385 if (present(psib2)) then
386 call exponential_lanczos_batch(te, namespace, mesh, hm, psib2, deltat2_, vmagnus)
387 if (present(inh_psib)) then
388 call exponential_lanczos_batch(te, namespace, mesh, hm, psib2, deltat2_, vmagnus, inh_psib)
389 end if
390 end if
391
392 case (exp_chebyshev)
393 if (present(inh_psib)) then
394 write(message(1), '(a)') 'Chebyshev exponential ("TDExponentialMethod = chebyshev")'
395 write(message(2), '(a)') 'with inhomogeneous term is not implemented'
396 call messages_fatal(2, namespace=namespace)
397 end if
398 ! initialize classes for computing coefficients
399 if (imag_time_) then
400 chebyshev_function => chebyshev_exp_imagtime_t(hm%spectral_half_span, hm%spectral_middle_point, deltat)
401 else
402 chebyshev_function => chebyshev_exp_t(hm%spectral_half_span, hm%spectral_middle_point, deltat)
403 end if
404 if (present(psib2)) call psib%copy_data_to(mesh%np, psib2)
405 call exponential_cheby_batch(te, namespace, mesh, hm, psib, deltat, chebyshev_function, vmagnus)
406 if (present(psib2)) then
407 call exponential_cheby_batch(te, namespace, mesh, hm, psib2, deltat2, chebyshev_function, vmagnus)
408 end if
409 deallocate(chebyshev_function)
410
411 end select
412
413 call profiling_out("EXPONENTIAL_BATCH")
415 end subroutine exponential_apply_batch
416
417 ! ---------------------------------------------------------
418 subroutine exponential_taylor_series_batch(te, namespace, mesh, hm, psib, deltat, psib2, deltat2, vmagnus, inh_psib, phik_shift)
419 type(exponential_t), intent(inout) :: te
420 type(namespace_t), intent(in) :: namespace
421 class(mesh_t), intent(in) :: mesh
422 class(hamiltonian_abst_t), intent(inout) :: hm
423 class(batch_t), intent(inout) :: psib
424 complex(real64), intent(in) :: deltat
425 class(batch_t), optional, intent(inout) :: psib2
426 complex(real64), optional, intent(in) :: deltat2
427 real(real64), optional, intent(in) :: vmagnus(:,:,:) !(mesh%np, hm%d%nspin, 2)
428 class(batch_t), optional, intent(inout) :: inh_psib
429 integer, optional, intent(in) :: phik_shift
430
431 complex(real64) :: zfact, zfact2
432 integer :: iter, denom, phik_shift_
433 logical :: zfact_is_real
434 class(batch_t), allocatable :: psi1b, hpsi1b
435
437 call profiling_in("EXP_TAYLOR_BATCH")
438
439 call psib%clone_to(psi1b)
440 call psib%clone_to(hpsi1b)
441
442 zfact = m_z1
443 zfact2 = m_z1
444 zfact_is_real = abs(deltat-real(deltat, real64)) < m_epsilon
445
446 if (present(psib2)) call psib%copy_data_to(mesh%np, psib2)
447
448 if (present(inh_psib)) then
449 zfact = zfact*deltat
450 call batch_axpy(mesh%np, real(zfact, real64), inh_psib, psib)
451
452 if (present(psib2)) then
453 zfact2 = zfact2*deltat2
454 call batch_axpy(mesh%np, real(zfact2, real64), inh_psib, psib2)
455 end if
456 end if
457
458 ! shift the denominator by this shift for the phi_k functions
459 phik_shift_ = optional_default(phik_shift, 0)
460
461 do iter = 1, te%exp_order
462 denom = iter+phik_shift_
463 if (present(inh_psib)) denom = denom + 1
464 zfact = zfact*(-m_zi*deltat)/denom
465 if (present(deltat2)) zfact2 = zfact2*(-m_zi*deltat2)/denom
466 zfact_is_real = .not. zfact_is_real
467 ! FIXME: need a test here for runaway exponential, e.g. for too large dt.
468 ! in runaway case the problem is really hard to trace back: the positions
469 ! go haywire on the first step of dynamics (often NaN) and with debugging options
470 ! the code stops in ZAXPY below without saying why.
471
472 if (iter /= 1) then
473 call operate_batch(hm, namespace, mesh, psi1b, hpsi1b, vmagnus)
474 else
475 if (present(inh_psib)) then
476 call operate_batch(hm, namespace, mesh, inh_psib, hpsi1b, vmagnus)
477 else
478 call operate_batch(hm, namespace, mesh, psib, hpsi1b, vmagnus)
479 end if
480 end if
481
482 if (zfact_is_real) then
483 call batch_axpy(mesh%np, real(zfact, real64), hpsi1b, psib)
484 if (present(psib2)) call batch_axpy(mesh%np, real(zfact2, real64), hpsi1b, psib2)
485 else
486 call batch_axpy(mesh%np, zfact, hpsi1b, psib)
487 if (present(psib2)) call batch_axpy(mesh%np, zfact2, hpsi1b, psib2)
488 end if
489
490 if (iter /= te%exp_order) call hpsi1b%copy_data_to(mesh%np, psi1b)
491
492 end do
493
494 call psi1b%end()
495 call hpsi1b%end()
496 safe_deallocate_a(psi1b)
497 safe_deallocate_a(hpsi1b)
498
499 call profiling_out("EXP_TAYLOR_BATCH")
502
503 subroutine exponential_lanczos_batch(te, namespace, mesh, hm, psib, deltat, vmagnus, inh_psib)
504 type(exponential_t), intent(inout) :: te
505 type(namespace_t), intent(in) :: namespace
506 class(mesh_t), intent(in) :: mesh
507 class(hamiltonian_abst_t), intent(inout) :: hm
508 class(batch_t), intent(inout) :: psib
509 complex(real64), intent(in) :: deltat
510 real(real64), optional, intent(in) :: vmagnus(:,:,:) !(mesh%np, hm%d%nspin, 2)
511 class(batch_t), optional, intent(in) :: inh_psib
512
513 class(batch_t), allocatable :: tmpb
514
516
517 if (present(inh_psib)) then
518 call inh_psib%clone_to(tmpb, copy_data=.true.)
519 ! psib = psib + deltat * phi1(-i*deltat*H) inh_psib
520 call exponential_lanczos_function_batch(te, namespace, mesh, hm, tmpb, deltat, phi1, vmagnus)
521 call batch_axpy(mesh%np, real(deltat, real64), tmpb, psib)
522 call tmpb%end()
523 safe_deallocate_a(tmpb)
524 else
525 call exponential_lanczos_function_batch(te, namespace, mesh, hm, psib, deltat, exponential, vmagnus)
526 end if
527
529 end subroutine exponential_lanczos_batch
530
531 ! ---------------------------------------------------------
543 subroutine exponential_lanczos_function_batch(te, namespace, mesh, hm, psib, deltat, fun, vmagnus)
544 type(exponential_t), intent(inout) :: te
545 type(namespace_t), intent(in) :: namespace
546 class(mesh_t), intent(in) :: mesh
547 class(hamiltonian_abst_t), intent(inout) :: hm
548 class(batch_t), intent(inout) :: psib
549 complex(real64), intent(in) :: deltat
550 interface
551 complex(real64) function fun(z)
552 import real64
553 complex(real64), intent(in) :: z
554 end
555 end interface
556 real(real64), optional, intent(in) :: vmagnus(:,:,:) !(mesh%np, hm%d%nspin, 2)
557
558 integer :: iter, l, ii, ist, max_initialized
559 complex(real64), allocatable :: hamilt(:,:,:), expo(:,:,:)
560 real(real64), allocatable :: beta(:), res(:), norm(:)
561 integer, parameter :: max_order = 200
562 type(batch_p_t), allocatable :: vb(:) ! Krylov subspace vectors
563
565 call profiling_in("EXP_LANCZOS_FUN_BATCH")
566
567 if (te%exp_method /= exp_lanczos) then
568 message(1) = "The exponential method needs to be set to Lanzcos (TDExponentialMethod=lanczos)."
569 call messages_fatal(1)
570 end if
571
572 safe_allocate(beta(1:psib%nst))
573 safe_allocate(res(1:psib%nst))
574 safe_allocate(norm(1:psib%nst))
575 safe_allocate(vb(1:max_order))
576 call psib%clone_to(vb(1)%p)
577 max_initialized = 1
578
579 call psib%copy_data_to(mesh%np, vb(1)%p, async=.true.)
580 call mesh_batch_nrm2(mesh, vb(1)%p, beta)
581
582 if (te%full_batch) beta = norm2(beta)
583
584 ! If we have a null vector, no need to compute the action of the exponential.
585 if (all(abs(beta) <= 1.0e-12_real64)) then
586 safe_deallocate_a(beta)
587 safe_deallocate_a(res)
588 safe_deallocate_a(norm)
589 call vb(1)%p%end()
590 safe_deallocate_a(vb)
591 call profiling_out("EXP_LANCZOS_FUN_BATCH")
593 return
594 end if
595
596 call batch_scal(mesh%np, m_one/beta, vb(1)%p, a_full = .false.)
597
598 safe_allocate(hamilt(1:max_order+1, 1:max_order+1, 1:psib%nst))
599 safe_allocate( expo(1:max_order+1, 1:max_order+1, 1:psib%nst))
600
601 ! This is the Lanczos loop...
602 do iter = 1, max_order-1
603 call psib%clone_to(vb(iter + 1)%p)
604 max_initialized = iter + 1
605
606 ! to apply the Hamiltonian
607 call operate_batch(hm, namespace, mesh, vb(iter)%p, vb(iter+1)%p, vmagnus)
608
609 ! We use either the Lanczos method (Hermitian case) or the Arnoldi method
610 if (hm%is_hermitian()) then
611 l = max(1, iter - 1)
612 hamilt(1:max(l-1, 1), iter, 1:psib%nst) = m_zero
613 else
614 l = 1
615 if (iter > 2) then
616 hamilt(iter, 1:iter-2, 1:psib%nst) = m_zero
617 end if
618 end if
619
620 ! Orthogonalize against previous vectors
621 call zmesh_batch_orthogonalization(mesh, iter - l + 1, vb(l:iter), vb(iter+1)%p, &
622 normalize = .false., overlap = hamilt(l:iter, iter, 1:psib%nst), norm = hamilt(iter + 1, iter, 1:psib%nst), &
623 gs_scheme = te%arnoldi_gs, full_batch = te%full_batch)
624
625 ! We now need to compute exp(Hm), where Hm is the projection of the linear transformation
626 ! of the Hamiltonian onto the Krylov subspace Km
627 ! See Eq. 4
628 do ii = 1, psib%nst
629 call zlalg_matrix_function(iter, -m_zi*deltat, hamilt(:,:,ii), expo(:,:,ii), fun, hm%is_hermitian())
630 res(ii) = abs(hamilt(iter + 1, iter, ii) * abs(expo(iter, 1, ii)))
631 end do !ii
632
633 ! We now estimate the error we made. This is given by the formula denoted Er2 in Sec. 5.2
634 if (all(abs(hamilt(iter + 1, iter, :)) < 1.0e4_real64*m_epsilon)) exit ! "Happy breakdown"
635 ! We normalize only if the norm is non-zero
636 ! see http://www.netlib.org/utk/people/JackDongarra/etemplates/node216.html#alg:arn0
637 norm = m_one
638 do ist = 1, psib%nst
639 if (abs(hamilt(iter + 1, iter, ist)) >= 1.0e4_real64 * m_epsilon) then
640 norm(ist) = m_one / abs(hamilt(iter + 1, iter, ist))
641 end if
642 end do
643
644 call batch_scal(mesh%np, norm, vb(iter+1)%p, a_full = .false.)
645
646 if (iter > 3 .and. all(res < te%lanczos_tol)) exit
647
648 end do !iter
649
650 if (iter == max_order) then ! Here one should consider the possibility of the happy breakdown.
651 write(message(1),'(a,i5,a,es9.2)') 'Lanczos exponential expansion did not converge after ', iter, &
652 ' iterations. Residual: ', maxval(res)
653 call messages_warning(1, namespace=namespace)
654 else
655 write(message(1),'(a,i5)') 'Debug: Lanczos exponential iterations: ', iter
656 call messages_info(1, namespace=namespace, debug_only=.true.)
657 end if
658
659 ! See Eq. 4 for the expression here
660 ! zpsi = nrm * V * expo(1:iter, 1) = nrm * V * expo * V^(T) * zpsi
661 call batch_scal(mesh%np, expo(1,1,1:psib%nst), psib, a_full = .false.)
662 ! TODO: We should have a routine batch_gemv for improved performance (see #1070 on gitlab)
663 do ii = 2, iter
664 call batch_axpy(mesh%np, beta(1:psib%nst)*expo(ii,1,1:psib%nst), vb(ii)%p, psib, a_full = .false.)
665 ! In order to apply the two exponentials, we must store the eigenvalues and eigenvectors given by zlalg_exp
666 ! And to recontruct here the exp(i*dt*H) for deltat2
667 end do
668
669 do ii = 1, max_initialized
670 call vb(ii)%p%end()
671 end do
672
673 safe_deallocate_a(vb)
674 safe_deallocate_a(hamilt)
675 safe_deallocate_a(expo)
676 safe_deallocate_a(beta)
677 safe_deallocate_a(res)
678 safe_deallocate_a(norm)
679
680 call accel_finish()
681
682 call profiling_out("EXP_LANCZOS_FUN_BATCH")
683
686
687
688 ! ---------------------------------------------------------
700 subroutine exponential_cheby_batch(te, namespace, mesh, hm, psib, deltat, chebyshev_function, vmagnus)
701 type(exponential_t), intent(inout) :: te
702 type(namespace_t), intent(in) :: namespace
703 class(mesh_t), intent(in) :: mesh
704 class(hamiltonian_abst_t), intent(inout) :: hm
705 class(batch_t), intent(inout) :: psib
706 real(real64), intent(in) :: deltat
707 class(chebyshev_function_t), intent(in) :: chebyshev_function
708 real(real64), optional, intent(in) :: vmagnus(:,:,:) !(mesh%np, hm%d%nspin, 2)
709
710 integer :: j, order_needed
711 complex(real64) :: coefficient
712 complex(real64), allocatable :: coefficients(:)
713 real(real64) :: error
714 class(batch_t), allocatable, target :: psi0, psi1, psi2
715 class(batch_t), pointer :: psi_n, psi_n1, psi_n2
716 integer, parameter :: max_order = 200
717
719 call profiling_in("EXP_CHEBY_BATCH")
720
721 call psib%clone_to(psi0)
722 call psib%clone_to(psi1)
723 call psib%clone_to(psi2)
724 call psib%copy_data_to(mesh%np, psi0)
725
726 order_needed = max_order
727 do j = 1, max_order
728 error = chebyshev_function%get_error(j)
729 if (error > m_zero .and. error < te%chebyshev_tol) then
730 order_needed = j
731 exit
732 end if
733 end do
734
735 call chebyshev_function%get_coefficients(j, coefficients)
736
737 ! zero-order term
738 call batch_scal(mesh%np, coefficients(0), psib)
739 ! first-order term
740 ! shifted Hamiltonian
741 call operate_batch(hm, namespace, mesh, psi0, psi1, vmagnus)
742 call batch_axpy(mesh%np, -hm%spectral_middle_point, psi0, psi1)
743 call batch_scal(mesh%np, m_one/hm%spectral_half_span, psi1)
744 ! accumulate result
745 call batch_axpy(mesh%np, coefficients(1), psi1, psib)
746
747 ! use pointers to avoid copies
748 psi_n => psi2
749 psi_n1 => psi1
750 psi_n2 => psi0
751
752 do j = 2, order_needed
753 ! compute shifted Hamiltonian and Chebyshev recurrence formula
754 call operate_batch(hm, namespace, mesh, psi_n1, psi_n, vmagnus)
755 call batch_axpy(mesh%np, -hm%spectral_middle_point, psi_n1, psi_n)
756 call batch_xpay(mesh%np, psi_n2, -m_two/hm%spectral_half_span, psi_n)
757 call batch_scal(mesh%np, -m_one, psi_n)
758
759 ! accumulate result
760 call batch_axpy(mesh%np, coefficients(j), psi_n, psib)
761
762 ! shift pointers for the three-term recurrence, this avoids copies
763 if (mod(j, 3) == 2) then
764 psi_n => psi0
765 psi_n1 => psi2
766 psi_n2 => psi1
767 else if (mod(j, 3) == 1) then
768 psi_n => psi2
769 psi_n1 => psi1
770 psi_n2 => psi0
771 else
772 psi_n => psi1
773 psi_n1 => psi0
774 psi_n2 => psi2
775 end if
776 end do
777
778 if (order_needed == max_order) then
779 write(message(1),'(a,i5,a,es9.2)') 'Chebyshev exponential expansion did not converge after ', j, &
780 ' iterations. Coefficient: ', coefficient
781 call messages_warning(1, namespace=namespace)
782 else
783 write(message(1),'(a,i5)') 'Debug: Chebyshev exponential iterations: ', j
784 call messages_info(1, namespace=namespace, debug_only=.true.)
785 end if
786
787 call psi0%end()
788 call psi1%end()
789 call psi2%end()
790 safe_deallocate_a(psi0)
791 safe_deallocate_a(psi1)
792 safe_deallocate_a(psi2)
794 safe_deallocate_a(coefficients)
795
796 call profiling_out("EXP_CHEBY_BATCH")
797
799 end subroutine exponential_cheby_batch
800
801
802 subroutine operate_batch(hm, namespace, mesh, psib, hpsib, vmagnus)
803 class(hamiltonian_abst_t), intent(inout) :: hm
804 type(namespace_t), intent(in) :: namespace
805 class(mesh_t), intent(in) :: mesh
806 class(batch_t), intent(inout) :: psib
807 class(batch_t), intent(inout) :: hpsib
808 real(real64), optional, intent(in) :: vmagnus(:, :, :)
809
810 push_sub(operate_batch)
811
812 if (present(vmagnus)) then
813 call hm%zmagnus_apply(namespace, mesh, psib, hpsib, vmagnus)
814 else
815 call hm%zapply(namespace, mesh, psib, hpsib)
816 end if
817
818 pop_sub(operate_batch)
819 end subroutine operate_batch
820
821 ! ---------------------------------------------------------
825 subroutine exponential_apply_all(te, namespace, mesh, hm, st, deltat, order)
826 type(exponential_t), intent(inout) :: te
827 type(namespace_t), intent(in) :: namespace
828 class(mesh_t), intent(inout) :: mesh
829 type(hamiltonian_elec_t), intent(inout) :: hm
830 type(states_elec_t), intent(inout) :: st
831 real(real64), intent(in) :: deltat
832 integer, optional, intent(inout) :: order
833
834 integer :: ik, ib, i
835 real(real64) :: zfact
836
837 type(states_elec_t) :: st1, hst1
838
839 push_sub(exponential_apply_all)
840
841 assert(te%exp_method == exp_taylor)
842
843 call states_elec_copy(st1, st)
844 call states_elec_copy(hst1, st)
845
846 zfact = m_one
847 do i = 1, te%exp_order
848 zfact = zfact * deltat / i
849
850 if (i == 1) then
851 call zhamiltonian_elec_apply_all(hm, namespace, mesh, st, hst1)
852 else
853 call zhamiltonian_elec_apply_all(hm, namespace, mesh, st1, hst1)
854 end if
855
856 do ik = st%d%kpt%start, st%d%kpt%end
857 do ib = st%group%block_start, st%group%block_end
858 call batch_set_zero(st1%group%psib(ib, ik))
859 call batch_axpy(mesh%np, -m_zi, hst1%group%psib(ib, ik), st1%group%psib(ib, ik))
860 call batch_axpy(mesh%np, zfact, st1%group%psib(ib, ik), st%group%psib(ib, ik))
861 end do
862 end do
863
864 end do
865 ! End of Taylor expansion loop.
866
867 call states_elec_end(st1)
868 call states_elec_end(hst1)
869
870 ! We now add the inhomogeneous part, if present.
871 if (hamiltonian_elec_inh_term(hm)) then
872 !write(*, *) 'Now we apply the inhomogeneous term...'
873
874 call states_elec_copy(st1, hm%inh_st)
875 call states_elec_copy(hst1, hm%inh_st)
876
877
878 do ik = st%d%kpt%start, st%d%kpt%end
879 do ib = st%group%block_start, st%group%block_end
880 call batch_axpy(mesh%np, deltat, st1%group%psib(ib, ik), st%group%psib(ib, ik))
881 end do
882 end do
883
884 zfact = m_one
885 do i = 1, te%exp_order
886 zfact = zfact * deltat / (i+1)
887
888 if (i == 1) then
889 call zhamiltonian_elec_apply_all(hm, namespace, mesh, hm%inh_st, hst1)
890 else
891 call zhamiltonian_elec_apply_all(hm, namespace, mesh, st1, hst1)
892 end if
893
894 do ik = st%d%kpt%start, st%d%kpt%end
895 do ib = st%group%block_start, st%group%block_end
896 call batch_set_zero(st1%group%psib(ib, ik))
897 call batch_axpy(mesh%np, -m_zi, hst1%group%psib(ib, ik), st1%group%psib(ib, ik))
898 call batch_axpy(mesh%np, deltat * zfact, st1%group%psib(ib, ik), st%group%psib(ib, ik))
899 end do
900 end do
901
902 end do
903
904 call states_elec_end(st1)
905 call states_elec_end(hst1)
906
907 end if
908
909 if (present(order)) order = te%exp_order*st%nik*st%nst ! This should be the correct number
910
911 pop_sub(exponential_apply_all)
912 end subroutine exponential_apply_all
913
914 subroutine exponential_apply_phi_batch(te, namespace, mesh, hm, psib, deltat, k, vmagnus)
915 class(exponential_t), intent(inout) :: te
916 type(namespace_t), intent(in) :: namespace
917 class(mesh_t), intent(in) :: mesh
918 class(hamiltonian_abst_t), intent(inout) :: hm
919 class(batch_t), intent(inout) :: psib
920 real(real64), intent(in) :: deltat
921 integer, intent(in) :: k
922 real(real64), optional, intent(in) :: vmagnus(:,:,:) !(mesh%np, hm%d%nspin, 2)
923
924 class(chebyshev_function_t), pointer :: chebyshev_function
925 complex(real64) :: deltat_
926
927 push_sub_with_profile(exponential_apply_phi_batch)
928
929 assert(psib%type() == type_cmplx)
930 if (present(vmagnus)) then
931 select type(hm)
932 class is (hamiltonian_elec_t)
933 assert(size(vmagnus, 1) >= mesh%np)
934 assert(size(vmagnus, 2) == hm%d%nspin)
935 assert(size(vmagnus, 3) == 2)
936 class default
937 write(message(1), '(a)') 'Magnus operators only implemented for electrons at the moment'
938 call messages_fatal(1, namespace=namespace)
939 end select
940 end if
941
942 if (.not. hm%is_hermitian() .and. te%exp_method == exp_chebyshev) then
943 write(message(1), '(a)') 'The Chebyshev expansion for the exponential will only converge if the imaginary'
944 write(message(2), '(a)') 'eigenvalues are small enough compared to the span of the real eigenvalues,'
945 write(message(3), '(a)') 'i.e., for ratios smaller than about 1e-3.'
946 write(message(4), '(a)') 'The Lanczos method ("TDExponentialMethod = lanczos") is guaranteed to'
947 write(message(5), '(a)') 'always converge in this case.'
948 call messages_warning(5, namespace=namespace)
949 end if
950
951 deltat_ = cmplx(deltat, m_zero, real64)
952
953 select case (te%exp_method)
954 case (exp_taylor)
955 call exponential_taylor_series_batch(te, namespace, mesh, hm, psib, deltat_, vmagnus=vmagnus, phik_shift=k)
956
957 case (exp_lanczos)
958 if (k == 1) then
959 call exponential_lanczos_function_batch(te, namespace, mesh, hm, psib, deltat_, phi1, vmagnus)
960 else if (k == 2) then
961 call exponential_lanczos_function_batch(te, namespace, mesh, hm, psib, deltat_, phi2, vmagnus)
962 else
963 write(message(1), '(a)') 'Lanczos expansion not implemented for phi_k, k > 2'
964 call messages_fatal(1, namespace=namespace)
965 end if
966
967 case (exp_chebyshev)
968 if (k == 1) then
969 chebyshev_function => chebyshev_numerical_t(hm%spectral_half_span, hm%spectral_middle_point, deltat, phi1)
970 else if (k == 2) then
971 chebyshev_function => chebyshev_numerical_t(hm%spectral_half_span, hm%spectral_middle_point, deltat, phi2)
972 else
973 write(message(1), '(a)') 'Chebyshev expansion not implemented for phi_k, k > 2'
974 call messages_fatal(1, namespace=namespace)
975 end if
976 call exponential_cheby_batch(te, namespace, mesh, hm, psib, deltat, chebyshev_function, vmagnus)
977 deallocate(chebyshev_function)
978
979 end select
980 pop_sub_with_profile(exponential_apply_phi_batch)
981 end subroutine exponential_apply_phi_batch
982
983end module exponential_oct_m
984
985!! Local Variables:
986!! mode: f90
987!! coding: utf-8
988!! End:
batchified version of the BLAS axpy routine:
Definition: batch_ops.F90:154
scale a batch by a constant or vector
Definition: batch_ops.F90:162
batchified version of
Definition: batch_ops.F90:170
subroutine, public accel_finish()
Definition: accel.F90:1298
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
subroutine, public batch_set_zero(this, np, async)
fill all mesh functions of the batch with zero
Definition: batch_ops.F90:242
This module contains interfaces for BLAS routines You should not use these routines directly....
Definition: blas.F90:118
subroutine, public exponential_copy(teo, tei)
subroutine, public exponential_apply_all(te, namespace, mesh, hm, st, deltat, order)
Note that this routine not only computes the exponential, but also an extra term if there is a inhomo...
subroutine operate_batch(hm, namespace, mesh, psib, hpsib, vmagnus)
subroutine exponential_apply_phi_batch(te, namespace, mesh, hm, psib, deltat, k, vmagnus)
subroutine, public exponential_init(te, namespace, full_batch)
integer, parameter, public exp_taylor
subroutine exponential_lanczos_batch(te, namespace, mesh, hm, psib, deltat, vmagnus, inh_psib)
subroutine exponential_apply_single(te, namespace, mesh, hm, zpsi, ist, ik, deltat, vmagnus, imag_time)
Wrapper to batchified routine for applying exponential to an array.
subroutine exponential_taylor_series_batch(te, namespace, mesh, hm, psib, deltat, psib2, deltat2, vmagnus, inh_psib, phik_shift)
subroutine exponential_apply_batch(te, namespace, mesh, hm, psib, deltat, psib2, deltat2, vmagnus, imag_time, inh_psib)
This routine performs the operation:
subroutine, public exponential_lanczos_function_batch(te, namespace, mesh, hm, psib, deltat, fun, vmagnus)
Compute fun(H) psib, i.e. the application of a function of the Hamiltonian to a batch.
integer, parameter, public exp_chebyshev
subroutine exponential_cheby_batch(te, namespace, mesh, hm, psib, deltat, chebyshev_function, vmagnus)
Calculates the exponential of the Hamiltonian through an expansion in Chebyshev polynomials.
real(real64), parameter, public m_two
Definition: global.F90:190
real(real64), parameter, public m_zero
Definition: global.F90:188
complex(real64), parameter, public m_zi
Definition: global.F90:202
real(real64), parameter, public m_epsilon
Definition: global.F90:204
complex(real64), parameter, public m_z1
Definition: global.F90:199
real(real64), parameter, public m_one
Definition: global.F90:189
This module defines an abstract class for Hamiltonians.
subroutine, public zhamiltonian_elec_apply_all(hm, namespace, mesh, st, hst)
pure logical function, public hamiltonian_elec_inh_term(hm)
subroutine, public zlalg_matrix_function(n, factor, a, fun_a, fun, hermitian)
This routine calculates a function of a matrix by using an eigenvalue decomposition.
Definition: lalg_adv.F90:2182
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:115
complex(real64) pure function, public phi2(z)
Compute phi2(z) = (phi1(z)-1)/z = (exp(z) - z - 1)/z^2.
Definition: math.F90:974
complex(real64) pure function, public exponential(z)
Wrapper for exponential.
Definition: math.F90:943
complex(real64) pure function, public phi1(z)
Compute phi1(z) = (exp(z)-1)/z.
Definition: math.F90:955
This module defines functions over batches of mesh functions.
Definition: mesh_batch.F90:116
subroutine, public mesh_batch_nrm2(mesh, aa, nrm2, reduce)
Calculate the norms (norm2, not the square!) of a batch of mesh functions.
Definition: mesh_batch.F90:177
subroutine, public zmesh_batch_orthogonalization(mesh, nst, psib, phib, normalize, overlap, norm, gs_scheme, full_batch)
Orthonormalizes states of phib to the orbitals of nst batches of psi.
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_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:538
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:161
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:415
subroutine, public messages_input_error(namespace, var, details, row, column)
Definition: messages.F90:714
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:617
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 states_elec_end(st)
finalize the 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
type(type_t), parameter, public type_cmplx
Definition: types.F90:134
Definition: xc.F90:114
Class defining batches of mesh functions.
Definition: batch.F90:159
The abstract Hamiltonian class defines a skeleton for specific implementations.
Describes mesh distribution to nodes.
Definition: mesh.F90:186
The states_elec_t class contains all electronic wave functions.
batches of electronic states
Definition: wfs_elec.F90:138
int true(void)