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 if (debug%info) then
653 write(message(1),'(a,i5)') 'Lanczos exponential iterations: ', iter
654 call messages_info(1, namespace=namespace)
655 end if
656 end if
657
658 ! See Eq. 4 for the expression here
659 ! zpsi = nrm * V * expo(1:iter, 1) = nrm * V * expo * V^(T) * zpsi
660 call batch_scal(mesh%np, expo(1,1,1:psib%nst), psib, a_full = .false.)
661 ! TODO: We should have a routine batch_gemv for improved performance (see #1070 on gitlab)
662 do ii = 2, iter
663 call batch_axpy(mesh%np, beta(1:psib%nst)*expo(ii,1,1:psib%nst), vb(ii)%p, psib, a_full = .false.)
664 ! In order to apply the two exponentials, we must store the eigenvalues and eigenvectors given by zlalg_exp
665 ! And to recontruct here the exp(i*dt*H) for deltat2
666 end do
667
668 do ii = 1, max_initialized
669 call vb(ii)%p%end()
670 end do
671
672 safe_deallocate_a(vb)
673 safe_deallocate_a(hamilt)
674 safe_deallocate_a(expo)
675 safe_deallocate_a(beta)
676 safe_deallocate_a(res)
677 safe_deallocate_a(norm)
678
679 call accel_finish()
680
681 call profiling_out("EXP_LANCZOS_FUN_BATCH")
682
685
686
687 ! ---------------------------------------------------------
699 subroutine exponential_cheby_batch(te, namespace, mesh, hm, psib, deltat, chebyshev_function, vmagnus)
700 type(exponential_t), intent(inout) :: te
701 type(namespace_t), intent(in) :: namespace
702 class(mesh_t), intent(in) :: mesh
703 class(hamiltonian_abst_t), intent(inout) :: hm
704 class(batch_t), intent(inout) :: psib
705 real(real64), intent(in) :: deltat
706 class(chebyshev_function_t), intent(in) :: chebyshev_function
707 real(real64), optional, intent(in) :: vmagnus(:,:,:) !(mesh%np, hm%d%nspin, 2)
708
709 integer :: j, order_needed
710 complex(real64) :: coefficient
711 complex(real64), allocatable :: coefficients(:)
712 real(real64) :: error
713 class(batch_t), allocatable, target :: psi0, psi1, psi2
714 class(batch_t), pointer :: psi_n, psi_n1, psi_n2
715 integer, parameter :: max_order = 200
716
718 call profiling_in("EXP_CHEBY_BATCH")
719
720 call psib%clone_to(psi0)
721 call psib%clone_to(psi1)
722 call psib%clone_to(psi2)
723 call psib%copy_data_to(mesh%np, psi0)
724
725 order_needed = max_order
726 do j = 1, max_order
727 error = chebyshev_function%get_error(j)
728 if (error > m_zero .and. error < te%chebyshev_tol) then
729 order_needed = j
730 exit
731 end if
732 end do
733
734 call chebyshev_function%get_coefficients(j, coefficients)
735
736 ! zero-order term
737 call batch_scal(mesh%np, coefficients(0), psib)
738 ! first-order term
739 ! shifted Hamiltonian
740 call operate_batch(hm, namespace, mesh, psi0, psi1, vmagnus)
741 call batch_axpy(mesh%np, -hm%spectral_middle_point, psi0, psi1)
742 call batch_scal(mesh%np, m_one/hm%spectral_half_span, psi1)
743 ! accumulate result
744 call batch_axpy(mesh%np, coefficients(1), psi1, psib)
745
746 ! use pointers to avoid copies
747 psi_n => psi2
748 psi_n1 => psi1
749 psi_n2 => psi0
750
751 do j = 2, order_needed
752 ! compute shifted Hamiltonian and Chebyshev recurrence formula
753 call operate_batch(hm, namespace, mesh, psi_n1, psi_n, vmagnus)
754 call batch_axpy(mesh%np, -hm%spectral_middle_point, psi_n1, psi_n)
755 call batch_xpay(mesh%np, psi_n2, -m_two/hm%spectral_half_span, psi_n)
756 call batch_scal(mesh%np, -m_one, psi_n)
757
758 ! accumulate result
759 call batch_axpy(mesh%np, coefficients(j), psi_n, psib)
760
761 ! shift pointers for the three-term recurrence, this avoids copies
762 if (mod(j, 3) == 2) then
763 psi_n => psi0
764 psi_n1 => psi2
765 psi_n2 => psi1
766 else if (mod(j, 3) == 1) then
767 psi_n => psi2
768 psi_n1 => psi1
769 psi_n2 => psi0
770 else
771 psi_n => psi1
772 psi_n1 => psi0
773 psi_n2 => psi2
774 end if
775 end do
776
777 if (order_needed == max_order) then
778 write(message(1),'(a,i5,a,es9.2)') 'Chebyshev exponential expansion did not converge after ', j, &
779 ' iterations. Coefficient: ', coefficient
780 call messages_warning(1, namespace=namespace)
781 else
782 if (debug%info) then
783 write(message(1),'(a,i5)') 'Chebyshev exponential iterations: ', j
784 call messages_info(1, namespace=namespace)
785 end if
786 end if
787
788 call psi0%end()
789 call psi1%end()
790 call psi2%end()
791 safe_deallocate_a(psi0)
792 safe_deallocate_a(psi1)
793 safe_deallocate_a(psi2)
794
795 safe_deallocate_a(coefficients)
796
797 call profiling_out("EXP_CHEBY_BATCH")
798
800 end subroutine exponential_cheby_batch
801
802
803 subroutine operate_batch(hm, namespace, mesh, psib, hpsib, vmagnus)
804 class(hamiltonian_abst_t), intent(inout) :: hm
805 type(namespace_t), intent(in) :: namespace
806 class(mesh_t), intent(in) :: mesh
807 class(batch_t), intent(inout) :: psib
808 class(batch_t), intent(inout) :: hpsib
809 real(real64), optional, intent(in) :: vmagnus(:, :, :)
810
811 push_sub(operate_batch)
812
813 if (present(vmagnus)) then
814 call hm%zmagnus_apply(namespace, mesh, psib, hpsib, vmagnus)
815 else
816 call hm%zapply(namespace, mesh, psib, hpsib)
817 end if
818
819 pop_sub(operate_batch)
820 end subroutine operate_batch
821
822 ! ---------------------------------------------------------
826 subroutine exponential_apply_all(te, namespace, mesh, hm, st, deltat, order)
827 type(exponential_t), intent(inout) :: te
828 type(namespace_t), intent(in) :: namespace
829 class(mesh_t), intent(inout) :: mesh
830 type(hamiltonian_elec_t), intent(inout) :: hm
831 type(states_elec_t), intent(inout) :: st
832 real(real64), intent(in) :: deltat
833 integer, optional, intent(inout) :: order
834
835 integer :: ik, ib, i
836 real(real64) :: zfact
837
838 type(states_elec_t) :: st1, hst1
839
840 push_sub(exponential_apply_all)
841
842 assert(te%exp_method == exp_taylor)
843
844 call states_elec_copy(st1, st)
845 call states_elec_copy(hst1, st)
846
847 zfact = m_one
848 do i = 1, te%exp_order
849 zfact = zfact * deltat / i
850
851 if (i == 1) then
852 call zhamiltonian_elec_apply_all(hm, namespace, mesh, st, hst1)
853 else
854 call zhamiltonian_elec_apply_all(hm, namespace, mesh, st1, hst1)
855 end if
856
857 do ik = st%d%kpt%start, st%d%kpt%end
858 do ib = st%group%block_start, st%group%block_end
859 call batch_set_zero(st1%group%psib(ib, ik))
860 call batch_axpy(mesh%np, -m_zi, hst1%group%psib(ib, ik), st1%group%psib(ib, ik))
861 call batch_axpy(mesh%np, zfact, st1%group%psib(ib, ik), st%group%psib(ib, ik))
862 end do
863 end do
864
865 end do
866 ! End of Taylor expansion loop.
867
868 call states_elec_end(st1)
869 call states_elec_end(hst1)
870
871 ! We now add the inhomogeneous part, if present.
872 if (hamiltonian_elec_inh_term(hm)) then
873 !write(*, *) 'Now we apply the inhomogeneous term...'
874
875 call states_elec_copy(st1, hm%inh_st)
876 call states_elec_copy(hst1, hm%inh_st)
877
878
879 do ik = st%d%kpt%start, st%d%kpt%end
880 do ib = st%group%block_start, st%group%block_end
881 call batch_axpy(mesh%np, deltat, st1%group%psib(ib, ik), st%group%psib(ib, ik))
882 end do
883 end do
884
885 zfact = m_one
886 do i = 1, te%exp_order
887 zfact = zfact * deltat / (i+1)
888
889 if (i == 1) then
890 call zhamiltonian_elec_apply_all(hm, namespace, mesh, hm%inh_st, hst1)
891 else
892 call zhamiltonian_elec_apply_all(hm, namespace, mesh, st1, hst1)
893 end if
894
895 do ik = st%d%kpt%start, st%d%kpt%end
896 do ib = st%group%block_start, st%group%block_end
897 call batch_set_zero(st1%group%psib(ib, ik))
898 call batch_axpy(mesh%np, -m_zi, hst1%group%psib(ib, ik), st1%group%psib(ib, ik))
899 call batch_axpy(mesh%np, deltat * zfact, st1%group%psib(ib, ik), st%group%psib(ib, ik))
900 end do
901 end do
902
903 end do
904
905 call states_elec_end(st1)
906 call states_elec_end(hst1)
907
908 end if
909
910 if (present(order)) order = te%exp_order*st%nik*st%nst ! This should be the correct number
911
912 pop_sub(exponential_apply_all)
913 end subroutine exponential_apply_all
914
915 subroutine exponential_apply_phi_batch(te, namespace, mesh, hm, psib, deltat, k, vmagnus)
916 class(exponential_t), intent(inout) :: te
917 type(namespace_t), intent(in) :: namespace
918 class(mesh_t), intent(in) :: mesh
919 class(hamiltonian_abst_t), intent(inout) :: hm
920 class(batch_t), intent(inout) :: psib
921 real(real64), intent(in) :: deltat
922 integer, intent(in) :: k
923 real(real64), optional, intent(in) :: vmagnus(:,:,:) !(mesh%np, hm%d%nspin, 2)
924
925 class(chebyshev_function_t), pointer :: chebyshev_function
926 complex(real64) :: deltat_
927
928 push_sub_with_profile(exponential_apply_phi_batch)
929
930 assert(psib%type() == type_cmplx)
931 if (present(vmagnus)) then
932 select type(hm)
933 class is (hamiltonian_elec_t)
934 assert(size(vmagnus, 1) >= mesh%np)
935 assert(size(vmagnus, 2) == hm%d%nspin)
936 assert(size(vmagnus, 3) == 2)
937 class default
938 write(message(1), '(a)') 'Magnus operators only implemented for electrons at the moment'
939 call messages_fatal(1, namespace=namespace)
940 end select
941 end if
942
943 if (.not. hm%is_hermitian() .and. te%exp_method == exp_chebyshev) then
944 write(message(1), '(a)') 'The Chebyshev expansion for the exponential will only converge if the imaginary'
945 write(message(2), '(a)') 'eigenvalues are small enough compared to the span of the real eigenvalues,'
946 write(message(3), '(a)') 'i.e., for ratios smaller than about 1e-3.'
947 write(message(4), '(a)') 'The Lanczos method ("TDExponentialMethod = lanczos") is guaranteed to'
948 write(message(5), '(a)') 'always converge in this case.'
949 call messages_warning(5, namespace=namespace)
950 end if
951
952 deltat_ = cmplx(deltat, m_zero, real64)
953
954 select case (te%exp_method)
955 case (exp_taylor)
956 call exponential_taylor_series_batch(te, namespace, mesh, hm, psib, deltat_, vmagnus=vmagnus, phik_shift=k)
957
958 case (exp_lanczos)
959 if (k == 1) then
960 call exponential_lanczos_function_batch(te, namespace, mesh, hm, psib, deltat_, phi1, vmagnus)
961 else if (k == 2) then
962 call exponential_lanczos_function_batch(te, namespace, mesh, hm, psib, deltat_, phi2, vmagnus)
963 else
964 write(message(1), '(a)') 'Lanczos expansion not implemented for phi_k, k > 2'
965 call messages_fatal(1, namespace=namespace)
966 end if
967
968 case (exp_chebyshev)
969 if (k == 1) then
970 chebyshev_function => chebyshev_numerical_t(hm%spectral_half_span, hm%spectral_middle_point, deltat, phi1)
971 else if (k == 2) then
972 chebyshev_function => chebyshev_numerical_t(hm%spectral_half_span, hm%spectral_middle_point, deltat, phi2)
973 else
974 write(message(1), '(a)') 'Chebyshev expansion not implemented for phi_k, k > 2'
975 call messages_fatal(1, namespace=namespace)
976 end if
977 call exponential_cheby_batch(te, namespace, mesh, hm, psib, deltat, chebyshev_function, vmagnus)
978 deallocate(chebyshev_function)
979
980 end select
981 pop_sub_with_profile(exponential_apply_phi_batch)
982 end subroutine exponential_apply_phi_batch
983
984end module exponential_oct_m
985
986!! Local Variables:
987!! mode: f90
988!! coding: utf-8
989!! End:
batchified version of the BLAS axpy routine:
Definition: batch_ops.F90:152
scale a batch by a constant or vector
Definition: batch_ops.F90:160
batchified version of
Definition: batch_ops.F90:168
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:240
This module contains interfaces for BLAS routines You should not use these routines directly....
Definition: blas.F90:118
type(debug_t), save, public debug
Definition: debug.F90:154
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
subroutine, public messages_info(no_lines, iunit, verbose_limit, stress, all_nodes, namespace)
Definition: messages.F90:624
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 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)