29  use, 
intrinsic :: iso_fortran_env
 
   60    class(mesh_t),     
intent(in)    :: mesh
 
   61    type(namespace_t), 
intent(in)    :: namespace
 
   62    class(space_t),    
intent(in)    :: space
 
   63    type(target_t),    
intent(inout) :: tg
 
   64    type(td_t),        
intent(in)    :: td
 
   65    type(restart_t),   
intent(in)    :: restart
 
   66    type(kpoints_t),   
intent(in)    :: kpoints
 
   72    message(1) =  
'Info: Using Ground State for TargetOperator' 
   77      message(1) = 
"Unable to read wavefunctions." 
   81    tg%move_ions = td%ions_dyn%ions_move()
 
   90    type(target_t),      
intent(in) :: tg
 
   91    type(namespace_t),   
intent(in) :: namespace
 
   92    class(space_t),      
intent(in)    :: space
 
   93    type(grid_t),        
intent(in) :: gr
 
   94    character(len=*),    
intent(in) :: dir
 
   95    type(ions_t),        
intent(in) :: ions
 
   96    type(hamiltonian_elec_t), 
intent(in) :: hm
 
   97    type(output_t),      
intent(in) :: outp
 
  102    call output_states(outp, namespace, space, trim(dir), tg%st, gr, ions, hm, -1)
 
  111  real(real64) function target_j1_groundstate(tg, gr, psi) result(j1)
 
  112    type(target_t),      
intent(in) :: tg
 
  113    type(grid_t),        
intent(in) :: gr
 
  114    type(states_elec_t), 
intent(in) :: psi
 
  117    complex(real64), 
allocatable :: zpsi(:, :), zst(:, :)
 
  119    push_sub(target_j1_groundstate)
 
  121    safe_allocate(zpsi(1:gr%np, 1:tg%st%d%dim))
 
  122    safe_allocate(zst(1:gr%np, 1:tg%st%d%dim))
 
  127      do ist = psi%st_start, psi%st_end
 
  131        j1 = j1 + psi%occ(ist, ik)*abs(
zmf_dotp(gr, psi%d%dim, zpsi, zst))**2
 
  136    safe_deallocate_a(zpsi)
 
  137    safe_deallocate_a(zst)
 
  139    pop_sub(target_j1_groundstate)
 
  146    type(target_t),       
intent(in)    :: tg
 
  147    type(grid_t),         
intent(in)    :: gr
 
  148    type(states_elec_t),  
intent(in)    :: psi_in
 
  149    type(states_elec_t),  
intent(inout) :: chi_out
 
  152    complex(real64) :: olap
 
  153    complex(real64), 
allocatable :: zpsi(:, :), zst(:, :), zchi(:, :)
 
  157    safe_allocate(zpsi(1:gr%np, 1:tg%st%d%dim))
 
  158    safe_allocate(zst(1:gr%np, 1:tg%st%d%dim))
 
  159    safe_allocate(zchi(1:gr%np, 1:tg%st%d%dim))
 
  161    do ik = 1, psi_in%nik
 
  162      do ist = psi_in%st_start, psi_in%st_end
 
  164        call states_elec_get_state(psi_in, gr, ist, ik, zpsi)
 
  165        call states_elec_get_state(tg%st, gr, ist, ik, zst)
 
  167        olap = zmf_dotp(gr, zst(:, 1), zpsi(:, 1))
 
  168        zchi(1:gr%np, 1:tg%st%d%dim) = olap*zst(1:gr%np, 1:tg%st%d%dim)
 
  170        call states_elec_set_state(chi_out, gr, ist, ik, zchi)
 
  175    safe_deallocate_a(zpsi)
 
  176    safe_deallocate_a(zst)
 
  177    safe_deallocate_a(zchi)
 
real(real64), parameter, public m_zero
 
This module implements the underlying real-space grid.
 
subroutine, public io_mkdir(fname, namespace, parents)
 
This module defines various routines, operating on mesh functions.
 
This module defines the meshes, which are used in Octopus.
 
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 messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
 
this module contains the low-level part of the output system
 
this module contains the output system
 
subroutine, public output_states(outp, namespace, space, dir, st, gr, ions, hm, iter)
 
This module handles reading and writing restart information for the states_elec_t.
 
subroutine, public states_elec_load(restart, namespace, space, st, mesh, kpoints, ierr, iter, lr, lowest_missing, label, verbose, skip)
returns in ierr: <0 => Fatal error, or nothing read =0 => read all wavefunctions >0 => could only rea...
 
subroutine, public target_output_groundstate(tg, namespace, space, gr, dir, ions, hm, outp)
 
subroutine, public target_init_groundstate(mesh, namespace, space, tg, td, restart, kpoints)
 
subroutine, public target_chi_groundstate(tg, gr, psi_in, chi_out)
 
real(real64) function, public target_j1_groundstate(tg, gr, psi)