Octopus
propagation_oct_m Module Reference

Data Types

type  oct_prop_t
 

Functions/Subroutines

subroutine, public propagation_mod_init (niter, eta, delta, number_checkpoints, zbr98, gradients)
 This subroutine must be called before any QOCT propagations are done. It simply stores in the module some data that is needed for the propagations, and which should stay invariant during the whole run. There is no need for any propagation_mod_close. More...
 
subroutine, public propagate_forward (sys, td, par, tg, qcpsi, prop, write_iter)
 
subroutine, public propagate_backward (sys, td, qcpsi, prop)
 
subroutine, public fwd_step (sys, td, tg, par, par_chi, qcpsi, prop_chi, prop_psi)
 
subroutine, public bwd_step (sys, td, tg, par, par_chi, qcchi, prop_chi, prop_psi)
 
subroutine, public bwd_step_2 (sys, td, tg, par, par_chi, qcchi, prop_chi, prop_psi)
 
subroutine update_hamiltonian_elec_chi (iter, namespace, space, gr, ks, hm, ext_partners, td, tg, par_chi, ions, st, qtildehalf)
 
subroutine update_hamiltonian_elec_psi (iter, namespace, space, gr, ks, hm, ext_partners, td, tg, par, st, ions)
 
subroutine calculate_g (space, gr, hm, lasers, psi, chi, dl, dq)
 
subroutine update_field (iter, cp, space, gr, hm, ext_partners, ions, qcpsi, qcchi, cpp, dir)
 Calculates the value of the control functions at iteration iter, from the state psi and the Lagrange-multiplier chi. More...
 
subroutine, public oct_prop_init (prop, namespace, dirname, mesh, mc)
 
subroutine, public oct_prop_end (prop)
 
subroutine, public oct_prop_check (prop, namespace, space, psi, mesh, kpoints, iter)
 
subroutine oct_prop_dump_states (prop, space, iter, psi, mesh, kpoints, ierr)
 
subroutine oct_prop_load_states (prop, namespace, space, psi, mesh, kpoints, iter, ierr)
 
subroutine vlaser_operator_quadratic (laser, mesh, space, psi, hpsi)
 
subroutine vlaser_operator_linear (laser, der, std, psi, hpsi, ik, gyromagnetic_ratio, a_static)
 

Variables

integer niter_
 Module variables. More...
 
integer number_checkpoints_
 
real(real64) eta_
 
real(real64) delta_
 
logical zbr98_
 
logical gradients_
 

Function/Subroutine Documentation

◆ propagation_mod_init()

subroutine, public propagation_oct_m::propagation_mod_init ( integer, intent(in)  niter,
real(real64), intent(in)  eta,
real(real64), intent(in)  delta,
integer, intent(in)  number_checkpoints,
logical, intent(in)  zbr98,
logical, intent(in)  gradients 
)

This subroutine must be called before any QOCT propagations are done. It simply stores in the module some data that is needed for the propagations, and which should stay invariant during the whole run. There is no need for any propagation_mod_close.

Definition at line 207 of file propagation.F90.

◆ propagate_forward()

subroutine, public propagation_oct_m::propagate_forward ( type(electrons_t), intent(inout)  sys,
type(td_t), intent(inout)  td,
type(controlfunction_t), intent(in)  par,
type(target_t), intent(inout)  tg,
type(opt_control_state_t), intent(inout)  qcpsi,
type(oct_prop_t), intent(inout), optional  prop,
logical, intent(in), optional  write_iter 
)

Performs a full propagation of state psi, with the laser field specified in par. If write_iter is present and is

set to .true., writes down through the td_write module.

Definition at line 236 of file propagation.F90.

◆ propagate_backward()

subroutine, public propagation_oct_m::propagate_backward ( type(electrons_t), intent(inout)  sys,
type(td_t), intent(inout)  td,
type(opt_control_state_t), intent(inout)  qcpsi,
type(oct_prop_t), intent(inout)  prop 
)

Performs a full backward propagation of state psi, with the

external fields specified in Hamiltonian h.

Definition at line 378 of file propagation.F90.

◆ fwd_step()

subroutine, public propagation_oct_m::fwd_step ( type(electrons_t), intent(inout)  sys,
type(td_t), intent(inout)  td,
type(target_t), intent(inout)  tg,
type(controlfunction_t), intent(inout)  par,
type(controlfunction_t), intent(in)  par_chi,
type(opt_control_state_t), intent(inout)  qcpsi,
type(oct_prop_t), intent(inout)  prop_chi,
type(oct_prop_t), intent(inout)  prop_psi 
)

Performs a forward propagation on the state psi and on the Lagrange-multiplier state chi. It also updates the control function par, according to the following scheme:

|chi> --> U[par_chi](T, 0)|chi> par = par[|psi>, |chi>] |psi> --> U[par](T, 0)|psi>

Note that the control functions "par" are updated on the fly, so that the propagation of psi is performed with the

"new" control functions.

Definition at line 444 of file propagation.F90.

◆ bwd_step()

subroutine, public propagation_oct_m::bwd_step ( type(electrons_t), intent(inout)  sys,
type(td_t), intent(inout)  td,
type(target_t), intent(inout)  tg,
type(controlfunction_t), intent(in)  par,
type(controlfunction_t), intent(inout)  par_chi,
type(opt_control_state_t), intent(inout)  qcchi,
type(oct_prop_t), intent(inout)  prop_chi,
type(oct_prop_t), intent(inout)  prop_psi 
)

Performs a backward propagation on the state psi and on the Lagrange-multiplier state chi, according to the following scheme:

|psi> --> U[par](0, T)|psi> par_chi = par_chi[|psi>, |chi>]

|chi> --> U[par_chi](0, T)|chi>

Definition at line 573 of file propagation.F90.

◆ bwd_step_2()

subroutine, public propagation_oct_m::bwd_step_2 ( type(electrons_t), intent(inout)  sys,
type(td_t), intent(inout)  td,
type(target_t), intent(inout)  tg,
type(controlfunction_t), intent(in)  par,
type(controlfunction_t), intent(inout)  par_chi,
type(opt_control_state_t), intent(inout)  qcchi,
type(oct_prop_t), intent(inout)  prop_chi,
type(oct_prop_t), intent(inout)  prop_psi 
)

Performs a backward propagation on the state psi and on the Lagrange-multiplier state chi, according to the following scheme:

|psi> --> U[par](0, T)|psi> |chi> --> U[par](0, T)|chi>

It also calculates during the propagation, a new "output" field:

par_chi = par_chi[|psi>, |chi>]

Definition at line 675 of file propagation.F90.

◆ update_hamiltonian_elec_chi()

subroutine propagation_oct_m::update_hamiltonian_elec_chi ( integer, intent(in)  iter,
type(namespace_t), intent(in)  namespace,
class(space_t), intent(in)  space,
type(grid_t), intent(inout)  gr,
type(v_ks_t), intent(inout)  ks,
type(hamiltonian_elec_t), intent(inout)  hm,
type(partner_list_t), intent(in)  ext_partners,
type(td_t), intent(inout)  td,
type(target_t), intent(inout)  tg,
type(controlfunction_t), intent(in)  par_chi,
type(ions_t), intent(in)  ions,
type(states_elec_t), intent(inout)  st,
real(real64), dimension(:, :), intent(in), optional  qtildehalf 
)
private

Definition at line 886 of file propagation.F90.

◆ update_hamiltonian_elec_psi()

subroutine propagation_oct_m::update_hamiltonian_elec_psi ( integer, intent(in)  iter,
type(namespace_t), intent(in)  namespace,
type(electron_space_t), intent(in)  space,
type(grid_t), intent(inout)  gr,
type(v_ks_t), intent(inout)  ks,
type(hamiltonian_elec_t), intent(inout)  hm,
type(partner_list_t), intent(in)  ext_partners,
type(td_t), intent(inout)  td,
type(target_t), intent(inout)  tg,
type(controlfunction_t), intent(in)  par,
type(states_elec_t), intent(inout)  st,
type(ions_t), intent(in)  ions 
)
private

Definition at line 979 of file propagation.F90.

◆ calculate_g()

subroutine propagation_oct_m::calculate_g ( class(space_t), intent(in)  space,
type(grid_t), intent(inout)  gr,
type(hamiltonian_elec_t), intent(in)  hm,
type(lasers_t), intent(in)  lasers,
type(states_elec_t), intent(inout)  psi,
type(states_elec_t), intent(inout)  chi,
complex(real64), dimension(:), intent(inout)  dl,
complex(real64), dimension(:), intent(inout)  dq 
)
private

Definition at line 1028 of file propagation.F90.

◆ update_field()

subroutine propagation_oct_m::update_field ( integer, intent(in)  iter,
type(controlfunction_t), intent(inout)  cp,
class(space_t), intent(in)  space,
type(grid_t), intent(inout)  gr,
type(hamiltonian_elec_t), intent(in)  hm,
type(partner_list_t), intent(in)  ext_partners,
type(ions_t), intent(in)  ions,
type(opt_control_state_t), intent(inout)  qcpsi,
type(opt_control_state_t), intent(inout)  qcchi,
type(controlfunction_t), intent(in)  cpp,
character(len=1), intent(in)  dir 
)
private

Calculates the value of the control functions at iteration iter, from the state psi and the Lagrange-multiplier chi.

If dir = 'f', the field must be updated for a forward propagation. In that case, the propagation step that is going to be done moves from (iter-1)*|dt| to iter*|dt|.

If dir = 'b', the field must be updated for a backward propagation. In taht case, the propagation step that is going to be done moves from iter*|dt| to (iter-1)*|dt|.

cp = (1-eta)*cpp - (eta/alpha) * <chi|V|Psi>

Definition at line 1113 of file propagation.F90.

◆ oct_prop_init()

subroutine, public propagation_oct_m::oct_prop_init ( type(oct_prop_t), intent(inout)  prop,
type(namespace_t), intent(in)  namespace,
character(len=*), intent(in)  dirname,
class(mesh_t), intent(in)  mesh,
type(multicomm_t), intent(in)  mc 
)

Definition at line 1205 of file propagation.F90.

◆ oct_prop_end()

subroutine, public propagation_oct_m::oct_prop_end ( type(oct_prop_t), intent(inout)  prop)

Definition at line 1238 of file propagation.F90.

◆ oct_prop_check()

subroutine, public propagation_oct_m::oct_prop_check ( type(oct_prop_t), intent(inout)  prop,
type(namespace_t), intent(in)  namespace,
class(space_t), intent(in)  space,
type(states_elec_t), intent(inout)  psi,
class(mesh_t), intent(in)  mesh,
type(kpoints_t), intent(in)  kpoints,
integer, intent(in)  iter 
)

Definition at line 1255 of file propagation.F90.

◆ oct_prop_dump_states()

subroutine propagation_oct_m::oct_prop_dump_states ( type(oct_prop_t), intent(inout)  prop,
class(space_t), intent(in)  space,
integer, intent(in)  iter,
type(states_elec_t), intent(in)  psi,
class(mesh_t), intent(in)  mesh,
type(kpoints_t), intent(in)  kpoints,
integer, intent(out)  ierr 
)
private

Definition at line 1307 of file propagation.F90.

◆ oct_prop_load_states()

subroutine propagation_oct_m::oct_prop_load_states ( type(oct_prop_t), intent(inout)  prop,
type(namespace_t), intent(in)  namespace,
class(space_t), intent(in)  space,
type(states_elec_t), intent(inout)  psi,
class(mesh_t), intent(in)  mesh,
type(kpoints_t), intent(in)  kpoints,
integer, intent(in)  iter,
integer, intent(out)  ierr 
)
private

Definition at line 1360 of file propagation.F90.

◆ vlaser_operator_quadratic()

subroutine propagation_oct_m::vlaser_operator_quadratic ( type(laser_t), intent(in)  laser,
class(mesh_t), intent(in)  mesh,
class(space_t), intent(in)  space,
complex(real64), dimension(:,:), intent(inout)  psi,
complex(real64), dimension(:,:), intent(inout)  hpsi 
)
private
Parameters
[in,out]psipsi(dermeshnp_part, hddim)
[in,out]hpsihpsi(dermeshnp_part, hddim)

Definition at line 1414 of file propagation.F90.

◆ vlaser_operator_linear()

subroutine propagation_oct_m::vlaser_operator_linear ( type(laser_t), intent(in)  laser,
type(derivatives_t), intent(in)  der,
type(states_elec_dim_t), intent(in)  std,
complex(real64), dimension(:,:), intent(inout)  psi,
complex(real64), dimension(:,:), intent(inout)  hpsi,
integer, intent(in)  ik,
real(real64), intent(in)  gyromagnetic_ratio,
real(real64), dimension(:,:), intent(in), optional  a_static 
)
private

Definition at line 1475 of file propagation.F90.

Variable Documentation

◆ niter_

integer propagation_oct_m::niter_
private

Module variables.

Definition at line 193 of file propagation.F90.

◆ number_checkpoints_

integer propagation_oct_m::number_checkpoints_
private

Definition at line 194 of file propagation.F90.

◆ eta_

real(real64) propagation_oct_m::eta_
private

Definition at line 195 of file propagation.F90.

◆ delta_

real(real64) propagation_oct_m::delta_
private

Definition at line 196 of file propagation.F90.

◆ zbr98_

logical propagation_oct_m::zbr98_
private

Definition at line 197 of file propagation.F90.

◆ gradients_

logical propagation_oct_m::gradients_
private

Definition at line 198 of file propagation.F90.