37  use, 
intrinsic :: iso_fortran_env
 
   64  type(namespace_t),        
pointer, 
private :: namespace_p
 
   65  class(mesh_t),            
pointer, 
private :: mesh_p
 
   66  type(hamiltonian_elec_t), 
pointer, 
private :: hm_p
 
   67  type(propagator_base_t),  
pointer, 
private :: tr_p
 
   68  type(states_elec_t),      
pointer, 
private :: st_p
 
   69  integer,                           
private :: ik_op, ist_op, dim_op
 
   70  real(real64),                      
private :: t_op, dt_op
 
   76  subroutine td_crank_nicolson(hm, namespace, space, gr, st, tr, time, dt, ions_dyn, ions, ext_partners, use_sparskit)
 
   77    type(hamiltonian_elec_t), 
target, 
intent(inout) :: hm
 
   78    type(namespace_t),        
target, 
intent(in)    :: namespace
 
   79    type(electron_space_t),           
intent(in)    :: space
 
   80    type(grid_t),             
target, 
intent(inout) :: gr
 
   81    type(states_elec_t),      
target, 
intent(inout) :: st
 
   82    type(propagator_base_t),  
target, 
intent(inout) :: tr
 
   83    real(real64),                     
intent(in)    :: time
 
   84    real(real64),                     
intent(in)    :: dt
 
   85    type(ion_dynamics_t),             
intent(inout) :: ions_dyn
 
   86    type(ions_t),                     
intent(inout) :: ions
 
   87    type(partner_list_t),             
intent(in)    :: ext_partners
 
   88    logical,                          
intent(in)    :: use_sparskit
 
   90    complex(real64), 
allocatable :: zpsi(:), rhs(:)
 
   91    integer :: ik, ist, idim, np, max_iter
 
   92    real(real64) :: cgtol = 1.0e-12_real64
 
   93    integer :: ib, minst, maxst
 
   94    integer, 
allocatable :: iter_used(:)
 
   95    real(real64), 
allocatable :: residue(:)
 
   96    type(density_calc_t) :: dens_calc
 
   97    type(wfs_elec_t) :: psib_rhs
 
   98    type(gauge_field_t), 
pointer :: gfield
 
  104    if(
associated(gfield)) 
then 
  110    if (use_sparskit) 
then 
  111      message(1) = 
"Cannot use SPARSKIT in Crank-Nicolson propagator: not compiled with SPARSKIT support." 
  119    namespace_p => namespace
 
  125    t_op  = time - dt/
m_two 
  133    safe_allocate(zpsi(1:np*st%d%dim))
 
  134    safe_allocate(rhs(1:np*st%d%dim))
 
  141      time, dt, time -dt/
m_two, hm%vhxc, vtau = hm%vtau)
 
  148    do ik = st%d%kpt%start, st%d%kpt%end
 
  150      do ib = st%group%block_start, st%group%block_end
 
  154        call st%group%psib(ib, ik)%copy_to(psib_rhs)
 
  160          call batch_axpy(np, dt, hm%inh_st%group%psib(ib, ik), psib_rhs)
 
  163        if(use_sparskit) 
then  
  166          do ist = minst, maxst
 
  168            do idim = 1, st%d%dim
 
  169              call batch_get_state(st%group%psib(ib, ik), (/ist, idim/), np, zpsi((idim - 1)*np+1:idim*np))
 
  170              call batch_get_state(psib_rhs, (/ist, idim/), np, rhs((idim - 1)*np+1:idim*np))
 
  178            do idim = 1, st%d%dim
 
  179              call batch_set_state(st%group%psib(ib, ik), (/ist, idim/), gr%np, zpsi((idim - 1)*np+1:idim*np))
 
  186          safe_allocate(iter_used(psib_rhs%nst))
 
  187          safe_allocate(residue(psib_rhs%nst))
 
  189          call zbatch_qmr_dotu(namespace, gr, st, st%group%psib(ib, ik), psib_rhs, &
 
  192          if (any(iter_used == max_iter)) 
then 
  193            write(
message(1),
'(a)') 
'The linear solver used for the Crank-Nicolson' 
  194            write(
message(2),
'(a)') 
'propagator did not converge: ' 
  195            do ist=1, psib_rhs%nst
 
  196              write(
message(ist+2),
'(a,i2,a,es14.4)') 
'State = ', ist, 
' Residual = ', residue(ist)
 
  201          safe_deallocate_a(iter_used)
 
  202          safe_deallocate_a(residue)
 
  221    safe_deallocate_a(zpsi)
 
  222    safe_deallocate_a(rhs)
 
  230  subroutine td_zop(xre, xim, yre, yim)
 
  231    real(real64), 
intent(in)  :: xre(:)
 
  232    real(real64), 
intent(in)  :: xim(:)
 
  233    real(real64), 
intent(out) :: yre(:)
 
  234    real(real64), 
intent(out) :: yim(:)
 
  237    complex(real64), 
allocatable :: zpsi(:, :)
 
  241    safe_allocate(zpsi(1:mesh_p%np_part, 1:dim_op))
 
  244      zpsi(1:mesh_p%np, idim) = xre((idim-1)*mesh_p%np+1:idim*mesh_p%np) + &
 
  245        m_zi * xim((idim-1)*mesh_p%np+1:idim*mesh_p%np)
 
  248    call tr_p%te%apply_single(namespace_p, mesh_p, hm_p, zpsi, ist_op, ik_op, -dt_op/
m_two)
 
  251      yre((idim-1)*mesh_p%np+1:idim*mesh_p%np) = real(zpsi(1:mesh_p%np, idim), real64)
 
  252      yim((idim-1)*mesh_p%np+1:idim*mesh_p%np) = aimag(zpsi(1:mesh_p%np, idim))
 
  255    safe_deallocate_a(zpsi)
 
  264  subroutine td_zopt(xre, xim, yre, yim)
 
  265    real(real64), 
intent(in)  :: xre(:)
 
  266    real(real64), 
intent(in)  :: xim(:)
 
  267    real(real64), 
intent(out) :: yre(:)
 
  268    real(real64), 
intent(out) :: yim(:)
 
  271    complex(real64), 
allocatable :: zpsi(:, :)
 
  278    safe_allocate(zpsi(1:mesh_p%np_part, 1:dim_op))
 
  281      zpsi(1:mesh_p%np, idim) = xre((idim-1)*mesh_p%np+1:idim*mesh_p%np) - &
 
  282        m_zi * xim((idim-1)*mesh_p%np+1:idim*mesh_p%np)
 
  285    call tr_p%te%apply_single(namespace_p, mesh_p, hm_p, zpsi, ist_op, ik_op, -dt_op/
m_two)
 
  288      yre((idim-1)*mesh_p%np+1:idim*mesh_p%np) =   real(zpsi(1:mesh_p%np, idim), real64)
 
  289      yim((idim-1)*mesh_p%np+1:idim*mesh_p%np) = - aimag(zpsi(1:mesh_p%np, idim))
 
  292    safe_deallocate_a(zpsi)
 
  306    call xxb%copy_data_to(mesh_p%np, yyb)
 
  307    call tr_p%te%apply_batch(namespace_p, mesh_p, hm_p, yyb, -dt_op/
m_two)
 
  320  subroutine zbatch_qmr_dotu(namespace, mesh, st, xb, bb, op, max_iter, iter_used, &
 
  321    residue, threshold, use_initial_guess)
 
  323    class(
mesh_t),            
intent(in)    :: mesh
 
  335    integer,                  
intent(in)    :: max_iter
 
  336    integer,                  
intent(out)   :: iter_used(:)
 
  337    real(real64), 
contiguous,        
intent(out)   :: residue(:)
 
  338    real(real64),             
intent(in)    :: threshold
 
  339    logical, 
optional,        
intent(in)    :: use_initial_guess
 
  342    type(
wfs_elec_t) :: vvb, res, qqb, ppb, deltax, deltar
 
  344    real(real64), 
allocatable  :: rho(:), oldrho(:), norm_b(:), xsi(:), gamma(:), alpha(:), theta(:), oldtheta(:), saved_res(:)
 
  345    real(real64), 
allocatable  :: oldgamma(:)
 
  346    complex(real64), 
allocatable :: eta(:), beta(:), delta(:), eps(:), exception_saved(:, :, :)
 
  347    integer, 
allocatable :: status(:), saved_iter(:)
 
  349    integer, 
parameter ::        &
 
  350      QMR_NOT_CONVERGED    = 0,  &
 
  354      qmr_breakdown_pb     = 4,  &
 
  355      qmr_breakdown_vz     = 5,  &
 
  356      qmr_breakdown_qp     = 6,  &
 
  357      qmr_breakdown_gamma  = 7
 
  361    safe_allocate(rho(1:xb%nst))
 
  362    safe_allocate(oldrho(1:xb%nst))
 
  363    safe_allocate(norm_b(1:xb%nst))
 
  364    safe_allocate(xsi(1:xb%nst))
 
  365    safe_allocate(gamma(1:xb%nst))
 
  366    safe_allocate(oldgamma(1:xb%nst))
 
  367    safe_allocate(alpha(1:xb%nst))
 
  368    safe_allocate(eta(1:xb%nst))
 
  369    safe_allocate(theta(1:xb%nst))
 
  370    safe_allocate(oldtheta(1:xb%nst))
 
  371    safe_allocate(beta(1:xb%nst))
 
  372    safe_allocate(delta(1:xb%nst))
 
  373    safe_allocate(eps(1:xb%nst))
 
  374    safe_allocate(saved_res(1:xb%nst))
 
  376    safe_allocate(status(1:xb%nst))
 
  377    safe_allocate(saved_iter(1:xb%nst))
 
  379    safe_allocate(exception_saved(1:mesh%np, 1:st%d%dim, 1:xb%nst))
 
  385    call xb%copy_to(deltax)
 
  386    call xb%copy_to(deltar)
 
  397      call batch_xpay(mesh%np, bb, -1.0_real64, vvb)
 
  398      call vvb%copy_data_to(mesh%np, res)
 
  406        if (abs(rho(ii)) <= threshold) 
then 
  407          status(ii) = qmr_res_zero
 
  408          residue(ii) = rho(ii)
 
  411          saved_res(ii) = residue(ii)
 
  418      call bb%copy_data_to(mesh%np, vvb)
 
  419      call vvb%copy_data_to(mesh%np, res)
 
  425    status = qmr_not_converged
 
  431      residue(ii) = rho(ii)
 
  433      if (status(ii) == qmr_not_converged .and. abs(norm_b(ii)) <= 
m_tiny) 
then 
  435        status(ii) = qmr_b_zero
 
  436        residue(ii) = norm_b(ii)
 
  437        saved_iter(ii) = iter
 
  438        saved_res(ii) = residue(ii)
 
  450    do while(iter < max_iter)
 
  454      if (all(status /= qmr_not_converged)) 
exit 
  458        if (status(ii) == qmr_not_converged .and. (abs(rho(ii)) < 
m_tiny .or. abs(xsi(ii)) < 
m_tiny)) 
then 
  460          status(ii) = qmr_breakdown_pb
 
  461          saved_iter(ii) = iter
 
  462          saved_res(ii) = residue(ii)
 
  465        alpha(ii) = alpha(ii)*xsi(ii)/rho(ii)
 
  473      delta = delta / alpha
 
  478        if (status(ii) == qmr_not_converged .and. abs(delta(ii)) < 
m_tiny) 
then 
  480          status(ii) = qmr_breakdown_vz
 
  481          saved_iter(ii) = iter
 
  482          saved_res(ii) = residue(ii)
 
  488        call vvb%copy_data_to(mesh%np, qqb)
 
  491        call batch_xpay(mesh%np, vvb, -rho*delta/eps, qqb, a_full = .false.)
 
  497      call batch_scal(mesh%np, alpha, ppb, a_full = .false.)
 
  504        if (status(ii) == qmr_not_converged .and. abs(eps(ii)) < 
m_tiny) 
then 
  506          status(ii) = qmr_breakdown_qp
 
  507          saved_iter(ii) = iter
 
  508          saved_res(ii) = residue(ii)
 
  511        beta(ii) = eps(ii)/delta(ii)
 
  515      call batch_xpay(mesh%np, ppb, -beta, vvb, a_full = .false.)
 
  525      xsi = rho / (alpha**2)
 
  529        oldtheta(ii) = theta(ii)
 
  531        theta(ii) = rho(ii)/(gamma(ii)*abs(beta(ii)))
 
  533        oldgamma(ii) = gamma(ii)
 
  538        if (status(ii) == qmr_not_converged .and. abs(gamma(ii)) < 
m_tiny) 
then 
  540          status(ii) = qmr_breakdown_gamma
 
  541          saved_iter(ii) = iter
 
  542          saved_res(ii) = residue(ii)
 
  546        eta(ii) = -eta(ii)*oldrho(ii)*gamma(ii)**2/(beta(ii)*oldgamma(ii)**2)
 
  552        call qqb%copy_data_to(mesh%np, deltax)
 
  553        call batch_scal(mesh%np, eta*alpha, deltax, a_full = .false.)
 
  556        call ppb%copy_data_to(mesh%np, deltar)
 
  557        call batch_scal(mesh%np, eta, deltar, a_full = .false.)
 
  562        call batch_scal(mesh%np, (oldtheta*gamma)**2, deltax, a_full = .false.)
 
  563        call batch_axpy(mesh%np, eta*alpha, qqb, deltax, a_full = .false.)
 
  566        call batch_scal(mesh%np, (oldtheta*gamma)**2, deltar, a_full = .false.)
 
  567        call batch_axpy(mesh%np, eta, ppb, deltar, a_full = .false.)
 
  576      call batch_axpy(mesh%np, -1.0_real64, deltar, res)
 
  581        residue(ii) = residue(ii)/norm_b(ii)
 
  586        if (status(ii) == qmr_not_converged .and. residue(ii) < threshold) 
then 
  587          status(ii) = qmr_converged
 
  589            write(
message(1),*) 
'Debug: State ', xb%ist(ii), 
' converged in ', iter, 
' iterations.' 
  598      if (status(ii) == qmr_not_converged .or. status(ii) == qmr_converged) 
then 
  600        iter_used(ii) = iter -1
 
  603        iter_used(ii) = saved_iter(ii)
 
  604        residue(ii) = saved_res(ii)
 
  607      select case (status(ii))
 
  608      case (qmr_not_converged)
 
  609        write(
message(1), 
'(a)') 
"QMR solver not converged!" 
  610        write(
message(2), 
'(a)') 
"Try increasing the maximum number of iterations or the tolerance." 
  612      case (qmr_breakdown_pb)
 
  613        write(
message(1), 
'(a)') 
"QMR breakdown, cannot continue: b or P*b is the zero vector!" 
  615      case (qmr_breakdown_vz)
 
  616        write(
message(1), 
'(a)') 
"QMR breakdown, cannot continue: v^T*z is zero!" 
  618      case (qmr_breakdown_qp)
 
  619        write(
message(1), 
'(a)') 
"QMR breakdown, cannot continue: q^T*p is zero!" 
  621      case (qmr_breakdown_gamma)
 
  622        write(
message(1), 
'(a)') 
"QMR breakdown, cannot continue: gamma is zero!" 
  635    safe_deallocate_a(exception_saved)
 
  636    safe_deallocate_a(rho)
 
  637    safe_deallocate_a(oldrho)
 
  638    safe_deallocate_a(norm_b)
 
  639    safe_deallocate_a(xsi)
 
  640    safe_deallocate_a(gamma)
 
  641    safe_deallocate_a(oldgamma)
 
  642    safe_deallocate_a(alpha)
 
  643    safe_deallocate_a(eta)
 
  644    safe_deallocate_a(theta)
 
  645    safe_deallocate_a(oldtheta)
 
  646    safe_deallocate_a(beta)
 
  647    safe_deallocate_a(delta)
 
  648    safe_deallocate_a(eps)
 
  650    safe_deallocate_a(saved_res)
 
  651    safe_deallocate_a(status)
 
  652    safe_deallocate_a(saved_iter)
 
batchified version of the BLAS axpy routine:
 
scale a batch by a constant or vector
 
There are several ways how to call batch_set_state and batch_get_state:
 
double sqrt(double __x) __attribute__((__nothrow__
 
subroutine, public accel_set_stream(stream_number)
 
This module implements common operations on batches of mesh functions.
 
subroutine, public batch_set_zero(this, np, async)
fill all mesh functions of the batch with zero
 
type(debug_t), save, public debug
 
This module implements a calculator for the density and defines related functions.
 
subroutine, public density_calc_end(this, allreduce, symmetrize)
Finalize the density calculation.
 
subroutine, public density_calc_accumulate(this, psib, async)
Accumulate weighted orbital densities for a batch psib.
 
subroutine, public density_calc_init(this, st, gr, density)
initialize the density calculator
 
integer, parameter, public exp_taylor
 
type(gauge_field_t) function, pointer, public list_get_gauge_field(partners)
 
logical pure function, public gauge_field_is_propagated(this)
 
real(real64), parameter, public m_two
 
real(real64), parameter, public m_zero
 
real(real64), parameter, public m_tiny
 
complex(real64), parameter, public m_z0
 
complex(real64), parameter, public m_zi
 
real(real64), parameter, public m_half
 
real(real64), parameter, public m_one
 
This module implements the underlying real-space grid.
 
pure logical function, public hamiltonian_elec_inh_term(hm)
 
This module defines classes and functions for interaction partners.
 
This module defines functions over batches of mesh functions.
 
subroutine, public mesh_batch_nrm2(mesh, aa, nrm2, reduce)
Calculate the norms (norm2, not the square!) of a batch of mesh functions.
 
subroutine, public zmesh_batch_dotp_vector(mesh, aa, bb, dot, reduce, cproduct)
calculate the vector of dot-products of mesh functions between two batches
 
This module defines various routines, operating on mesh functions.
 
This module defines the meshes, which are used in Octopus.
 
subroutine, public messages_warning(no_lines, all_nodes, namespace)
 
subroutine, public messages_info(no_lines, iunit, verbose_limit, stress, all_nodes, namespace)
 
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
 
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
 
subroutine, public potential_interpolation_interpolate(potential_interpolation, order, time, dt, t, vhxc, vtau)
 
subroutine, public propagation_ops_elec_restore_ions(wo, ions_dyn, ions)
 
subroutine, public propagation_ops_elec_move_ions(wo, gr, hm, st, namespace, space, ions_dyn, ions, ext_partners, time, dt, save_pos)
 
subroutine, public propagation_ops_do_unpack(st, hm, ib, ik)
 
subroutine, public propagation_ops_elec_update_hamiltonian(namespace, space, st, mesh, hm, ext_partners, time)
 
subroutine, public propagation_ops_do_pack(st, hm, ib, ik)
 
subroutine, public propagation_ops_finish_unpack(st, hm, ib, ik)
 
subroutine td_zop(xre, xim, yre, yim)
operators for Crank-Nicolson scheme
 
subroutine zbatch_qmr_dotu(namespace, mesh, st, xb, bb, op, max_iter, iter_used, residue, threshold, use_initial_guess)
for complex symmetric matrices W Chen and B Poirier, J Comput Phys 219, 198-209 (2006)
 
subroutine td_zopt(xre, xim, yre, yim)
Transpose of H (called e.g. by bi-conjugate gradient solver)
 
subroutine, public td_crank_nicolson(hm, namespace, space, gr, st, tr, time, dt, ions_dyn, ions, ext_partners, use_sparskit)
Crank-Nicolson propagator.
 
subroutine propagator_qmr_op_batch(xxb, yyb)
operators for Crank-Nicolson scheme
 
This module is intended to contain "only mathematical" functions and procedures.
 
subroutine, public zsparskit_solver_run(namespace, sk, op, opt, sol, rhs)
 
integer pure function, public states_elec_block_max(st, ib)
return index of last state in block ib
 
integer pure function, public states_elec_block_min(st, ib)
return index of first state in block ib
 
Describes mesh distribution to nodes.
 
The states_elec_t class contains all electronic wave functions.
 
batches of electronic states
 
subroutine preconditioner(x, hx)
Simple preconditioner Here we need to approximate P^-1 We use the Jacobi approximation and that \nabl...