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-5
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-5_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)) < 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), inh_psib, psib)
451
452 if (present(psib2)) then
453 zfact2 = zfact2*deltat2
454 call batch_axpy(mesh%np, real(zfact2), 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), hpsi1b, psib)
484 if (present(psib2)) call batch_axpy(mesh%np, real(zfact2), 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 else
524 call exponential_lanczos_function_batch(te, namespace, mesh, hm, psib, deltat, exponential, vmagnus)
525 end if
526
528 end subroutine exponential_lanczos_batch
529
530 ! ---------------------------------------------------------
542 subroutine exponential_lanczos_function_batch(te, namespace, mesh, hm, psib, deltat, fun, vmagnus)
543 type(exponential_t), intent(inout) :: te
544 type(namespace_t), intent(in) :: namespace
545 class(mesh_t), intent(in) :: mesh
546 class(hamiltonian_abst_t), intent(inout) :: hm
547 class(batch_t), intent(inout) :: psib
548 complex(real64), intent(in) :: deltat
549 interface
550 complex(real64) function fun(z)
551 import real64
552 complex(real64), intent(in) :: z
553 end
554 end interface
555 real(real64), optional, intent(in) :: vmagnus(:,:,:) !(mesh%np, hm%d%nspin, 2)
556
557 integer :: iter, l, ii, ist, max_initialized
558 complex(real64), allocatable :: hamilt(:,:,:), expo(:,:,:)
559 real(real64), allocatable :: beta(:), res(:), norm(:)
560 integer, parameter :: max_order = 200
561 type(batch_p_t), allocatable :: vb(:) ! Krylov subspace vectors
562
564 call profiling_in("EXP_LANCZOS_FUN_BATCH")
565
566 if (te%exp_method /= exp_lanczos) then
567 message(1) = "The exponential method needs to be set to Lanzcos (TDExponentialMethod=lanczos)."
568 call messages_fatal(1)
569 end if
570
571 safe_allocate(beta(1:psib%nst))
572 safe_allocate(res(1:psib%nst))
573 safe_allocate(norm(1:psib%nst))
574 safe_allocate(vb(1:max_order))
575 call psib%clone_to(vb(1)%p)
576 max_initialized = 1
577
578 call psib%copy_data_to(mesh%np, vb(1)%p, async=.true.)
579 call mesh_batch_nrm2(mesh, vb(1)%p, beta)
580
581 if (te%full_batch) beta = norm2(beta)
582
583 ! If we have a null vector, no need to compute the action of the exponential.
584 if (all(abs(beta) <= 1.0e-12_real64)) then
585 safe_deallocate_a(beta)
586 safe_deallocate_a(res)
587 safe_deallocate_a(norm)
588 call vb(1)%p%end()
589 safe_deallocate_a(vb)
590 call profiling_out("EXP_LANCZOS_FUN_BATCH")
592 return
593 end if
594
595 call batch_scal(mesh%np, m_one/beta, vb(1)%p, a_full = .false.)
597 safe_allocate(hamilt(1:max_order+1, 1:max_order+1, 1:psib%nst))
598 safe_allocate( expo(1:max_order+1, 1:max_order+1, 1:psib%nst))
599 hamilt = m_z0
600 expo = m_z0
601
602 ! This is the Lanczos loop...
603 do iter = 1, max_order
604 call psib%clone_to(vb(iter + 1)%p)
605 max_initialized = iter + 1
606
607 ! to apply the Hamiltonian
608 call operate_batch(hm, namespace, mesh, vb(iter)%p, vb(iter+1)%p, vmagnus)
609
610 ! We use either the Lanczos method (Hermitian case) or the Arnoldi method
611 if (hm%is_hermitian()) then
612 l = max(1, iter - 1)
613 else
614 l = 1
615 end if
616
617 ! Orthogonalize against previous vectors
618 call zmesh_batch_orthogonalization(mesh, iter - l + 1, vb(l:iter), vb(iter+1)%p, &
619 normalize = .false., overlap = hamilt(l:iter, iter, 1:psib%nst), norm = hamilt(iter + 1, iter, 1:psib%nst), &
620 gs_scheme = te%arnoldi_gs, full_batch = te%full_batch)
621
622 ! We now need to compute exp(Hm), where Hm is the projection of the linear transformation
623 ! of the Hamiltonian onto the Krylov subspace Km
624 ! See Eq. 4
625 do ii = 1, psib%nst
626 call zlalg_matrix_function(iter, -m_zi*deltat, hamilt(:,:,ii), expo(:,:,ii), fun, hm%is_hermitian())
627 res(ii) = abs(hamilt(iter + 1, iter, ii) * abs(expo(iter, 1, ii)))
628 end do !ii
629
630 ! We now estimate the error we made. This is given by the formula denoted Er2 in Sec. 5.2
631 if (all(abs(hamilt(iter + 1, iter, :)) < 1.0e4_real64*m_epsilon)) exit ! "Happy breakdown"
632 ! We normalize only if the norm is non-zero
633 ! see http://www.netlib.org/utk/people/JackDongarra/etemplates/node216.html#alg:arn0
634 norm = m_one
635 do ist = 1, psib%nst
636 if (abs(hamilt(iter + 1, iter, ist)) >= 1.0e4_real64 * m_epsilon) then
637 norm(ist) = m_one / abs(hamilt(iter + 1, iter, ist))
638 end if
639 end do
640
641 call batch_scal(mesh%np, norm, vb(iter+1)%p, a_full = .false.)
642
643 if (iter > 3 .and. all(res < te%lanczos_tol)) exit
644
645 end do !iter
646
647 if (iter == max_order) then ! Here one should consider the possibility of the happy breakdown.
648 write(message(1),'(a,i5,a,es9.2)') 'Lanczos exponential expansion did not converge after ', iter, &
649 ' iterations. Residual: ', maxval(res)
650 call messages_warning(1, namespace=namespace)
651 else
652 write(message(1),'(a,i5)') 'Debug: Lanczos exponential iterations: ', iter
653 call messages_info(1, namespace=namespace, debug_only=.true.)
654 end if
655
656 ! See Eq. 4 for the expression here
657 ! zpsi = nrm * V * expo(1:iter, 1) = nrm * V * expo * V^(T) * zpsi
658 call batch_scal(mesh%np, expo(1,1,1:psib%nst), psib, a_full = .false.)
659 ! TODO: We should have a routine batch_gemv for improved performance (see #1070 on gitlab)
660 do ii = 2, iter
661 call batch_axpy(mesh%np, beta(1:psib%nst)*expo(ii,1,1:psib%nst), vb(ii)%p, psib, a_full = .false.)
662 ! In order to apply the two exponentials, we must store the eigenvalues and eigenvectors given by zlalg_exp
663 ! And to recontruct here the exp(i*dt*H) for deltat2
664 end do
665
666 do ii = 1, max_initialized
667 call vb(ii)%p%end()
668 end do
669
670 safe_deallocate_a(vb)
671 safe_deallocate_a(hamilt)
672 safe_deallocate_a(expo)
673 safe_deallocate_a(beta)
674 safe_deallocate_a(res)
675 safe_deallocate_a(norm)
676
677 call accel_finish()
678
679 call profiling_out("EXP_LANCZOS_FUN_BATCH")
680
683
684
685 ! ---------------------------------------------------------
697 subroutine exponential_cheby_batch(te, namespace, mesh, hm, psib, deltat, chebyshev_function, vmagnus)
698 type(exponential_t), intent(inout) :: te
699 type(namespace_t), intent(in) :: namespace
700 class(mesh_t), intent(in) :: mesh
701 class(hamiltonian_abst_t), intent(inout) :: hm
702 class(batch_t), intent(inout) :: psib
703 real(real64), intent(in) :: deltat
704 class(chebyshev_function_t), intent(in) :: chebyshev_function
705 real(real64), optional, intent(in) :: vmagnus(:,:,:) !(mesh%np, hm%d%nspin, 2)
706
707 integer :: j, order_needed
708 complex(real64) :: coefficient
709 complex(real64), allocatable :: coefficients(:)
710 real(real64) :: error
711 class(batch_t), allocatable, target :: psi0, psi1, psi2
712 class(batch_t), pointer :: psi_n, psi_n1, psi_n2
713 integer, parameter :: max_order = 200
714
716 call profiling_in("EXP_CHEBY_BATCH")
717
718 call psib%clone_to(psi0)
719 call psib%clone_to(psi1)
720 call psib%clone_to(psi2)
721 call psib%copy_data_to(mesh%np, psi0)
722
723 order_needed = max_order
724 do j = 1, max_order
725 error = chebyshev_function%get_error(j)
726 if (error > m_zero .and. error < te%chebyshev_tol) then
727 order_needed = j
728 exit
729 end if
730 end do
731
732 call chebyshev_function%get_coefficients(j, coefficients)
733
734 ! zero-order term
735 call batch_scal(mesh%np, coefficients(0), psib)
736 ! first-order term
737 ! shifted Hamiltonian
738 call operate_batch(hm, namespace, mesh, psi0, psi1, vmagnus)
739 call batch_axpy(mesh%np, -hm%spectral_middle_point, psi0, psi1)
740 call batch_scal(mesh%np, m_one/hm%spectral_half_span, psi1)
741 ! accumulate result
742 call batch_axpy(mesh%np, coefficients(1), psi1, psib)
743
744 ! use pointers to avoid copies
745 psi_n => psi2
746 psi_n1 => psi1
747 psi_n2 => psi0
748
749 do j = 2, order_needed
750 ! compute shifted Hamiltonian and Chebyshev recurrence formula
751 call operate_batch(hm, namespace, mesh, psi_n1, psi_n, vmagnus)
752 call batch_axpy(mesh%np, -hm%spectral_middle_point, psi_n1, psi_n)
753 call batch_xpay(mesh%np, psi_n2, -m_two/hm%spectral_half_span, psi_n)
754 call batch_scal(mesh%np, -m_one, psi_n)
755
756 ! accumulate result
757 call batch_axpy(mesh%np, coefficients(j), psi_n, psib)
758
759 ! shift pointers for the three-term recurrence, this avoids copies
760 if (mod(j, 3) == 2) then
761 psi_n => psi0
762 psi_n1 => psi2
763 psi_n2 => psi1
764 else if (mod(j, 3) == 1) then
765 psi_n => psi2
766 psi_n1 => psi1
767 psi_n2 => psi0
768 else
769 psi_n => psi1
770 psi_n1 => psi0
771 psi_n2 => psi2
772 end if
773 end do
774
775 if (order_needed == max_order) then
776 write(message(1),'(a,i5,a,es9.2)') 'Chebyshev exponential expansion did not converge after ', j, &
777 ' iterations. Coefficient: ', coefficient
778 call messages_warning(1, namespace=namespace)
779 else
780 write(message(1),'(a,i5)') 'Debug: Chebyshev exponential iterations: ', j
781 call messages_info(1, namespace=namespace, debug_only=.true.)
782 end if
783
784 call psi0%end()
785 call psi1%end()
786 call psi2%end()
787 safe_deallocate_a(psi0)
788 safe_deallocate_a(psi1)
789 safe_deallocate_a(psi2)
791 safe_deallocate_a(coefficients)
792
793 call profiling_out("EXP_CHEBY_BATCH")
794
796 end subroutine exponential_cheby_batch
797
798
799 subroutine operate_batch(hm, namespace, mesh, psib, hpsib, vmagnus)
800 class(hamiltonian_abst_t), intent(inout) :: hm
801 type(namespace_t), intent(in) :: namespace
802 class(mesh_t), intent(in) :: mesh
803 class(batch_t), intent(inout) :: psib
804 class(batch_t), intent(inout) :: hpsib
805 real(real64), optional, intent(in) :: vmagnus(:, :, :)
806
807 push_sub(operate_batch)
808
809 if (present(vmagnus)) then
810 call hm%zmagnus_apply(namespace, mesh, psib, hpsib, vmagnus)
811 else
812 call hm%zapply(namespace, mesh, psib, hpsib)
813 end if
814
815 pop_sub(operate_batch)
816 end subroutine operate_batch
817
818 ! ---------------------------------------------------------
822 subroutine exponential_apply_all(te, namespace, mesh, hm, st, deltat, order)
823 type(exponential_t), intent(inout) :: te
824 type(namespace_t), intent(in) :: namespace
825 class(mesh_t), intent(inout) :: mesh
826 type(hamiltonian_elec_t), intent(inout) :: hm
827 type(states_elec_t), intent(inout) :: st
828 real(real64), intent(in) :: deltat
829 integer, optional, intent(inout) :: order
830
831 integer :: ik, ib, i
832 real(real64) :: zfact
833
834 type(states_elec_t) :: st1, hst1
835
836 push_sub(exponential_apply_all)
837
838 assert(te%exp_method == exp_taylor)
839
840 call states_elec_copy(st1, st)
841 call states_elec_copy(hst1, st)
842
843 zfact = m_one
844 do i = 1, te%exp_order
845 zfact = zfact * deltat / i
846
847 if (i == 1) then
848 call zhamiltonian_elec_apply_all(hm, namespace, mesh, st, hst1)
849 else
850 call zhamiltonian_elec_apply_all(hm, namespace, mesh, st1, hst1)
851 end if
852
853 do ik = st%d%kpt%start, st%d%kpt%end
854 do ib = st%group%block_start, st%group%block_end
855 call batch_set_zero(st1%group%psib(ib, ik))
856 call batch_axpy(mesh%np, -m_zi, hst1%group%psib(ib, ik), st1%group%psib(ib, ik))
857 call batch_axpy(mesh%np, zfact, st1%group%psib(ib, ik), st%group%psib(ib, ik))
858 end do
859 end do
860
861 end do
862 ! End of Taylor expansion loop.
863
864 call states_elec_end(st1)
865 call states_elec_end(hst1)
866
867 ! We now add the inhomogeneous part, if present.
868 if (hamiltonian_elec_inh_term(hm)) then
869 !write(*, *) 'Now we apply the inhomogeneous term...'
870
871 call states_elec_copy(st1, hm%inh_st)
872 call states_elec_copy(hst1, hm%inh_st)
873
874
875 do ik = st%d%kpt%start, st%d%kpt%end
876 do ib = st%group%block_start, st%group%block_end
877 call batch_axpy(mesh%np, deltat, st1%group%psib(ib, ik), st%group%psib(ib, ik))
878 end do
879 end do
880
881 zfact = m_one
882 do i = 1, te%exp_order
883 zfact = zfact * deltat / (i+1)
884
885 if (i == 1) then
886 call zhamiltonian_elec_apply_all(hm, namespace, mesh, hm%inh_st, hst1)
887 else
888 call zhamiltonian_elec_apply_all(hm, namespace, mesh, st1, hst1)
889 end if
890
891 do ik = st%d%kpt%start, st%d%kpt%end
892 do ib = st%group%block_start, st%group%block_end
893 call batch_set_zero(st1%group%psib(ib, ik))
894 call batch_axpy(mesh%np, -m_zi, hst1%group%psib(ib, ik), st1%group%psib(ib, ik))
895 call batch_axpy(mesh%np, deltat * zfact, st1%group%psib(ib, ik), st%group%psib(ib, ik))
896 end do
897 end do
898
899 end do
900
901 call states_elec_end(st1)
902 call states_elec_end(hst1)
903
904 end if
905
906 if (present(order)) order = te%exp_order*st%nik*st%nst ! This should be the correct number
907
908 pop_sub(exponential_apply_all)
909 end subroutine exponential_apply_all
910
911 subroutine exponential_apply_phi_batch(te, namespace, mesh, hm, psib, deltat, k, vmagnus)
912 class(exponential_t), intent(inout) :: te
913 type(namespace_t), intent(in) :: namespace
914 class(mesh_t), intent(in) :: mesh
915 class(hamiltonian_abst_t), intent(inout) :: hm
916 class(batch_t), intent(inout) :: psib
917 real(real64), intent(in) :: deltat
918 integer, intent(in) :: k
919 real(real64), optional, intent(in) :: vmagnus(:,:,:) !(mesh%np, hm%d%nspin, 2)
920
921 class(chebyshev_function_t), pointer :: chebyshev_function
922 complex(real64) :: deltat_
923
924 push_sub_with_profile(exponential_apply_phi_batch)
925
926 assert(psib%type() == type_cmplx)
927 if (present(vmagnus)) then
928 select type(hm)
929 class is (hamiltonian_elec_t)
930 assert(size(vmagnus, 1) >= mesh%np)
931 assert(size(vmagnus, 2) == hm%d%nspin)
932 assert(size(vmagnus, 3) == 2)
933 class default
934 write(message(1), '(a)') 'Magnus operators only implemented for electrons at the moment'
935 call messages_fatal(1, namespace=namespace)
936 end select
937 end if
938
939 if (.not. hm%is_hermitian() .and. te%exp_method == exp_chebyshev) then
940 write(message(1), '(a)') 'The Chebyshev expansion for the exponential will only converge if the imaginary'
941 write(message(2), '(a)') 'eigenvalues are small enough compared to the span of the real eigenvalues,'
942 write(message(3), '(a)') 'i.e., for ratios smaller than about 1e-3.'
943 write(message(4), '(a)') 'The Lanczos method ("TDExponentialMethod = lanczos") is guaranteed to'
944 write(message(5), '(a)') 'always converge in this case.'
945 call messages_warning(5, namespace=namespace)
946 end if
947
948 deltat_ = cmplx(deltat, m_zero, real64)
949
950 select case (te%exp_method)
951 case (exp_taylor)
952 call exponential_taylor_series_batch(te, namespace, mesh, hm, psib, deltat_, vmagnus=vmagnus, phik_shift=k)
953
954 case (exp_lanczos)
955 if (k == 1) then
956 call exponential_lanczos_function_batch(te, namespace, mesh, hm, psib, deltat_, phi1, vmagnus)
957 else if (k == 2) then
958 call exponential_lanczos_function_batch(te, namespace, mesh, hm, psib, deltat_, phi2, vmagnus)
959 else
960 write(message(1), '(a)') 'Lanczos expansion not implemented for phi_k, k > 2'
961 call messages_fatal(1, namespace=namespace)
962 end if
963
964 case (exp_chebyshev)
965 if (k == 1) then
966 chebyshev_function => chebyshev_numerical_t(hm%spectral_half_span, hm%spectral_middle_point, deltat, phi1)
967 else if (k == 2) then
968 chebyshev_function => chebyshev_numerical_t(hm%spectral_half_span, hm%spectral_middle_point, deltat, phi2)
969 else
970 write(message(1), '(a)') 'Chebyshev expansion not implemented for phi_k, k > 2'
971 call messages_fatal(1, namespace=namespace)
972 end if
973 call exponential_cheby_batch(te, namespace, mesh, hm, psib, deltat, chebyshev_function, vmagnus)
974 deallocate(chebyshev_function)
975
976 end select
977 pop_sub_with_profile(exponential_apply_phi_batch)
978 end subroutine exponential_apply_phi_batch
979
980end module exponential_oct_m
981
982!! Local Variables:
983!! mode: f90
984!! coding: utf-8
985!! 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:1296
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:189
real(real64), parameter, public m_zero
Definition: global.F90:187
complex(real64), parameter, public m_z0
Definition: global.F90:197
complex(real64), parameter, public m_zi
Definition: global.F90:201
real(real64), parameter, public m_epsilon
Definition: global.F90:203
complex(real64), parameter, public m_z1
Definition: global.F90:198
real(real64), parameter, public m_one
Definition: global.F90:188
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:2089
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:970
complex(real64) pure function, public exponential(z)
Wrapper for exponential.
Definition: math.F90:939
complex(real64) pure function, public phi1(z)
Compute phi1(z) = (exp(z)-1)/z.
Definition: math.F90:951
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:543
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_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:624
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), 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)