92    logical                      :: pack_states
 
   93    logical                      :: parallel_in_states = .false. 
 
   94    integer, 
public              :: nst
 
   95    logical, 
public              :: packed
 
   97    complex(real64), 
allocatable           :: rs_state_plane_waves(:,:)
 
   98    complex(real64), 
allocatable           :: rs_state(:,:)
 
   99    complex(real64), 
allocatable           :: rs_state_prev(:,:)
 
  100    complex(real64), 
allocatable           :: rs_state_trans(:,:)
 
  101    complex(real64), 
allocatable           :: rs_state_long(:,:)
 
  102    complex(real64), 
allocatable           :: rs_current_density_t1(:,:)
 
  103    complex(real64), 
allocatable           :: rs_current_density_t2(:,:)
 
  105    logical                      :: rs_current_density_restart = .false.
 
  106    complex(real64), 
allocatable           :: rs_current_density_restart_t1(:,:)
 
  107    complex(real64), 
allocatable           :: rs_current_density_restart_t2(:,:)
 
  109    type(batch_t)                :: rs_stateb
 
  110    type(batch_t)                :: rs_state_prevb
 
  111    type(batch_t)                :: inhomogeneousb
 
  112    type(batch_t)                :: rs_state_plane_wavesb
 
  114    real(real64), 
allocatable           :: ep(:)
 
  115    real(real64), 
allocatable           :: mu(:)
 
  117    integer, 
allocatable         :: rs_state_fft_map(:,:,:)
 
  118    integer, 
allocatable         :: rs_state_fft_map_inv(:,:)
 
  120    real(real64)                 :: energy_rate
 
  121    real(real64)                 :: delta_energy
 
  122    real(real64)                 :: energy_via_flux_calc
 
  124    real(real64)                 :: trans_energy_rate
 
  125    real(real64)                 :: trans_delta_energy
 
  126    real(real64)                 :: trans_energy_via_flux_calc
 
  128    real(real64)                 :: plane_waves_energy_rate
 
  129    real(real64)                 :: plane_waves_delta_energy
 
  130    real(real64)                 :: plane_waves_energy_via_flux_calc
 
  132    real(real64)                 :: poynting_vector_box_surface(1:2,1:3,1:3) = 
m_zero 
  133    real(real64)                 :: poynting_vector_box_surface_plane_waves(1:2,1:3,1:3) = 
m_zero 
  134    real(real64)                 :: electric_field_box_surface(1:2,1:3,1:3) = 
m_zero 
  135    real(real64)                 :: electric_field_box_surface_plane_waves(1:2,1:3,1:3) = 
m_zero 
  136    real(real64)                 :: magnetic_field_box_surface(1:2,1:3,1:3) = 
m_zero 
  137    real(real64)                 :: magnetic_field_box_surface_plane_waves(1:2,1:3,1:3) = 
m_zero 
  139    logical                      :: rs_state_const_external = .false.
 
  140    complex(real64), 
allocatable           :: rs_state_const(:)
 
  141    complex(real64), 
allocatable           :: rs_state_const_amp(:,:)
 
  142    type(tdf_t), 
allocatable     :: rs_state_const_td_function(:)
 
  144    integer                      :: inner_points_number
 
  145    integer, 
allocatable         :: inner_points_map(:)
 
  146    logical, 
allocatable         :: inner_points_mask(:)
 
  147    integer                      :: boundary_points_number
 
  148    integer, 
allocatable         :: boundary_points_map(:)
 
  149    logical, 
allocatable         :: boundary_points_mask(:)
 
  150    type(accel_mem_t)            :: buff_inner_points_map, buff_boundary_points_map
 
  152    integer                      :: surface_points_number(3)
 
  153    integer, 
allocatable         :: surface_points_map(:,:,:)
 
  154    real(real64)                 :: surface_element(3)
 
  156    integer                      :: surface_grid_rows_number(3)
 
  157    integer, 
allocatable         :: surface_grid_points_number(:,:,:)
 
  158    integer(int64), 
allocatable     :: surface_grid_points_map(:,:,:,:,:)
 
  159    integer, 
allocatable         :: surface_grid_center(:,:,:,:)
 
  160    real(real64)                 :: surface_grid_element(3)
 
  162    type(mesh_plane_t)           :: surface(2,3)
 
  164    integer                      :: selected_points_number
 
  165    real(real64), 
allocatable           :: selected_points_coordinate(:,:)
 
  166    complex(real64), 
allocatable           :: selected_points_rs_state(:,:)
 
  167    complex(real64), 
allocatable           :: selected_points_rs_state_long(:,:)
 
  168    complex(real64), 
allocatable           :: selected_points_rs_state_trans(:,:)
 
  169    integer, 
allocatable         :: selected_points_map(:)
 
  170    type(accel_mem_t)            :: buff_selected_points_map
 
  171    real(real64)                 :: rs_state_trans_var
 
  173    real(real64), 
allocatable           :: grid_rho(:,:)
 
  174    complex(real64), 
allocatable           :: kappa_psi(:,:)
 
  176    character(len=1024), 
allocatable :: user_def_e_field(:)
 
  177    character(len=1024), 
allocatable :: user_def_b_field(:)
 
  179    integer                      :: energy_incident_waves_calc_iter
 
  180    logical                      :: energy_incident_waves_calc
 
  183    integer                      :: external_current_number
 
  184    integer,             
allocatable :: external_current_modus(:)
 
  185    character(len=1024), 
allocatable :: external_current_string(:,:)
 
  186    real(real64),        
allocatable :: external_current_amplitude(:,:,:)
 
  187    type(
tdf_t),         
allocatable :: external_current_td_function(:)
 
  188    type(
tdf_t),         
allocatable :: external_current_td_phase(:)
 
  189    real(real64),        
allocatable :: external_current_omega(:)
 
  190    real(real64),        
allocatable :: external_current_phase(:)
 
  193    character(len=1024), 
allocatable :: user_def_states(:,:,:)
 
  194    logical                     :: fromscratch = .
true.
 
  202    logical                     :: scalapack_compatible
 
  204    integer                     :: st_start, st_end
 
  205    integer, 
allocatable        :: node(:)
 
  208    integer                     :: transverse_field_mode
 
  212  integer, 
public, 
parameter ::      &
 
  225    integer :: idim, nlines, ncols, il
 
  226    real(real64), 
allocatable   :: pos(:)
 
  227    integer :: ix_max, iy_max, iz_max
 
  235    assert(space%dim == 3)
 
  239    safe_allocate(st%user_def_e_field(1:st%dim))
 
  240    safe_allocate(st%user_def_b_field(1:st%dim))
 
  246    safe_allocate(st%node(1:st%nst))
 
  247    st%node(1:st%nst) = 0
 
  250    st%parallel_in_states = .false.
 
  280    safe_allocate(pos(1:st%dim))
 
  281    st%selected_points_number = 1
 
  282    if (
parse_block(namespace, 
'MaxwellFieldsCoordinate', blk) == 0) 
then 
  284      st%selected_points_number = nlines
 
  285      safe_allocate(st%selected_points_coordinate(1:st%dim,1:nlines))
 
  286      safe_allocate(st%selected_points_rs_state(1:st%dim,1:nlines))
 
  287      safe_allocate(st%selected_points_rs_state_long(1:st%dim,1:nlines))
 
  288      safe_allocate(st%selected_points_rs_state_trans(1:st%dim,1:nlines))
 
  289      safe_allocate(st%selected_points_map(1:nlines))
 
  292        if (ncols < 3 .or. ncols > 3) 
then 
  293          message(1) = 
'MaxwellFieldCoordinate must have 3 columns.' 
  299        st%selected_points_coordinate(:,il) = pos
 
  300        st%selected_points_rs_state(:,il)  = 
m_z0 
  301        st%selected_points_rs_state_long(:,il) = 
m_z0 
  302        st%selected_points_rs_state_trans(:,il) = 
m_z0 
  306      safe_allocate(st%selected_points_coordinate(1:st%dim, 1))
 
  307      safe_allocate(st%selected_points_rs_state(1:st%dim, 1))
 
  308      safe_allocate(st%selected_points_rs_state_long(1:st%dim, 1))
 
  309      safe_allocate(st%selected_points_rs_state_trans(1:st%dim, 1))
 
  310      safe_allocate(st%selected_points_map(1))
 
  311      st%selected_points_coordinate(:,:) = 
m_zero 
  312      st%selected_points_rs_state(:,:) = 
m_z0 
  313      st%selected_points_rs_state_long(:,:) = 
m_z0 
  314      st%selected_points_rs_state_trans(:,:) = 
m_z0 
  315      st%selected_points_map(:) = -1
 
  318    safe_deallocate_a(pos)
 
  320    st%surface_grid_rows_number(1) = 3
 
  321    ix_max  = st%surface_grid_rows_number(1)
 
  322    st%surface_grid_rows_number(2) = 3
 
  323    iy_max  = st%surface_grid_rows_number(2)
 
  324    st%surface_grid_rows_number(3) = 3
 
  325    iz_max  = st%surface_grid_rows_number(3)
 
  327    safe_allocate(st%surface_grid_center(1:2, 1:st%dim, 1:ix_max, 1:iy_max))
 
  328    safe_allocate(st%surface_grid_points_number(1:st%dim, 1:ix_max, 1:iy_max))
 
  342      st%transverse_field_mode)
 
  354    class(
mesh_t),          
intent(in)      :: mesh
 
  361    safe_allocate(st%rs_state(1:mesh%np_part, 1:st%dim))
 
  362    st%rs_state(:,:) = 
m_z0 
  364    safe_allocate(st%rs_state_prev(1:mesh%np_part, 1:st%dim))
 
  365    st%rs_state_prev(:,:) = 
m_z0 
  367    safe_allocate(st%rs_state_trans(1:mesh%np_part, 1:st%dim))
 
  368    st%rs_state_trans(:,:) = 
m_z0 
  370    safe_allocate(st%rs_state_long(1:mesh%np_part, 1:st%dim))
 
  371    st%rs_state_long(:,:) = 
m_z0 
  373    safe_allocate(st%rs_state_plane_waves(1:mesh%np_part, 1:st%dim))
 
  374    st%rs_state_plane_waves(:,:) = 
m_z0 
  376    safe_allocate(st%rs_current_density_t1(1:mesh%np, 1:st%dim))
 
  377    st%rs_current_density_t1 = 
m_z0 
  379    safe_allocate(st%rs_current_density_t2(1:mesh%np, 1:st%dim))
 
  380    st%rs_current_density_t2 = 
m_z0 
  382    safe_allocate(st%rs_current_density_restart_t1(1:mesh%np_part, 1:st%dim))
 
  383    st%rs_current_density_restart_t1 = 
m_z0 
  385    safe_allocate(st%rs_current_density_restart_t2(1:mesh%np_part, 1:st%dim))
 
  386    st%rs_current_density_restart_t2 = 
m_z0 
  388    safe_allocate(st%ep(1:mesh%np_part))
 
  389    safe_allocate(st%mu(1:mesh%np_part))
 
  407    safe_deallocate_a(st%rs_state)
 
  408    safe_deallocate_a(st%rs_state_prev)
 
  409    safe_deallocate_a(st%rs_state_trans)
 
  410    safe_deallocate_a(st%selected_points_coordinate)
 
  411    safe_deallocate_a(st%selected_points_rs_state)
 
  412    safe_deallocate_a(st%selected_points_rs_state_long)
 
  413    safe_deallocate_a(st%selected_points_rs_state_trans)
 
  414    safe_deallocate_a(st%rs_current_density_t1)
 
  415    safe_deallocate_a(st%rs_current_density_t2)
 
  416    safe_deallocate_a(st%rs_state_long)
 
  417    safe_deallocate_a(st%rs_current_density_restart_t1)
 
  418    safe_deallocate_a(st%rs_current_density_restart_t2)
 
  419    safe_deallocate_a(st%user_def_e_field)
 
  420    safe_deallocate_a(st%user_def_b_field)
 
  422    safe_deallocate_a(st%rs_state_const)
 
  423    safe_deallocate_a(st%rs_state_const_td_function)
 
  424    safe_deallocate_a(st%rs_state_const_amp)
 
  425    safe_deallocate_a(st%rs_state_plane_waves)
 
  427    safe_deallocate_a(st%surface_grid_center)
 
  428    safe_deallocate_a(st%surface_grid_points_number)
 
  429    safe_deallocate_a(st%surface_grid_points_map)
 
  430    safe_deallocate_a(st%inner_points_map)
 
  431    safe_deallocate_a(st%boundary_points_map)
 
  432    safe_deallocate_a(st%inner_points_mask)
 
  433    safe_deallocate_a(st%boundary_points_mask)
 
  434    safe_deallocate_a(st%ep)
 
  435    safe_deallocate_a(st%mu)
 
  444    safe_deallocate_a(st%external_current_modus)
 
  445    safe_deallocate_a(st%external_current_string)
 
  446    safe_deallocate_a(st%external_current_amplitude)
 
  447    safe_deallocate_a(st%external_current_td_function)
 
  448    safe_deallocate_a(st%external_current_omega)
 
  449    safe_deallocate_a(st%external_current_td_phase)
 
  452    safe_deallocate_a(st%node)
 
  461  subroutine build_rs_element(e_element, b_element, rs_sign, rs_element, ep_element, mu_element)
 
  462    real(real64),      
intent(in)    :: e_element, b_element
 
  463    complex(real64),   
intent(inout) :: rs_element
 
  464    integer,           
intent(in)    :: rs_sign
 
  465    real(real64),   
optional, 
intent(in)    :: ep_element
 
  466    real(real64),   
optional, 
intent(in)    :: mu_element
 
  470    if (
present(ep_element) .and. 
present(mu_element)) 
then 
  480  subroutine build_rs_vector(e_vector, b_vector, rs_sign, rs_vector, ep_element, mu_element)
 
  481    real(real64),      
intent(in)    :: e_vector(:), b_vector(:)
 
  482    complex(real64),   
intent(inout) :: rs_vector(:)
 
  483    integer,           
intent(in)    :: rs_sign
 
  484    real(real64),   
optional, 
intent(in)    :: ep_element
 
  485    real(real64),   
optional, 
intent(in)    :: mu_element
 
  489    if (
present(ep_element) .and. 
present(mu_element)) 
then 
  499  subroutine build_rs_state(e_field, b_field, rs_sign, rs_state, mesh, ep_field, mu_field, np)
 
  500    real(real64),      
intent(in)    :: e_field(:,:), b_field(:,:)
 
  501    complex(real64),   
intent(inout) :: rs_state(:,:)
 
  502    integer,           
intent(in)    :: rs_sign
 
  503    class(
mesh_t),     
intent(in)    :: mesh
 
  504    real(real64),   
optional, 
intent(in)    :: ep_field(:)
 
  505    real(real64),   
optional, 
intent(in)    :: mu_field(:)
 
  506    integer, 
optional, 
intent(in)    :: np
 
  517      if (
present(ep_field) .and. 
present(mu_field)) 
then 
  518        rs_state(ip, :) = 
sqrt(ep_field(ip)/
m_two) * e_field(ip, :) &
 
  535    real(real64),    
intent(in)    :: current_element
 
  536    complex(real64), 
intent(inout) :: rs_current_element
 
  537    real(real64), 
optional, 
intent(in)    :: ep_element
 
  541    if (
present(ep_element)) 
then 
  542      rs_current_element = 
m_one/
sqrt(
m_two*ep_element) * current_element
 
  552    real(real64),    
intent(in)    :: current_vector(:)
 
  553    complex(real64), 
intent(inout) :: rs_current_vector(:)
 
  554    real(real64), 
optional, 
intent(in)    :: ep_element
 
  557    if (
present(ep_element)) 
then 
  568    real(real64),      
intent(in)    :: current_state(:,:)
 
  569    class(
mesh_t),     
intent(in)    :: mesh
 
  570    complex(real64),   
intent(inout) :: rs_current_state(:,:)
 
  571    real(real64),   
optional, 
intent(in)    :: ep_field(:)
 
  572    integer, 
optional, 
intent(in)    :: np
 
  574    integer :: ip, idim, np_, ff_dim
 
  581    ff_dim = 
size(current_state, dim=2)
 
  583    if (
present(ep_field)) 
then 
  586          rs_current_state(ip, idim) = 
m_one/
sqrt(
m_two*ep_field(ip)) * current_state(ip, idim)
 
  604    complex(real64),   
intent(in)    :: rs_state_vector(:)
 
  605    real(real64),      
intent(out)   :: electric_field_vector(:)
 
  606    real(real64),   
optional, 
intent(in)    :: ep_element
 
  610    if (
present(ep_element)) 
then 
  611      electric_field_vector(:) = 
sqrt(
m_two/ep_element) * real(rs_state_vector(:), real64)
 
  613      electric_field_vector(:) = 
sqrt(
m_two/
p_ep) * real(rs_state_vector(:), real64)
 
  621    complex(real64),   
intent(in)    :: rs_state_vector(:)
 
  622    integer,           
intent(in)    :: rs_sign
 
  623    real(real64),      
intent(out)   :: magnetic_field_vector(:)
 
  624    real(real64),   
optional, 
intent(in)    :: mu_element
 
  628    if (
present(mu_element)) 
then 
  629      magnetic_field_vector(:) = 
sqrt(
m_two*mu_element) * rs_sign * aimag(rs_state_vector(:))
 
  631      magnetic_field_vector(:) = 
sqrt(
m_two*
p_mu) * rs_sign * aimag(rs_state_vector(:))
 
  639    complex(real64),   
intent(in)    :: rs_state(:,:)
 
  640    class(
mesh_t),     
intent(in)    :: mesh
 
  641    real(real64),      
intent(out)   :: electric_field(:,:)
 
  642    real(real64),   
optional, 
intent(in)    :: ep_field(:)
 
  643    integer, 
optional, 
intent(in)    :: np
 
  654      if (
present(ep_field)) 
then 
  655        electric_field(ip, :) = 
sqrt(
m_two/ep_field(ip)) * real(rs_state(ip, :), real64)
 
  657        electric_field(ip,:) = 
sqrt(
m_two/
p_ep) * real(rs_state(ip, :), real64)
 
  670    complex(real64),   
intent(in)    :: rs_state(:,:)
 
  671    class(
mesh_t),     
intent(in)    :: mesh
 
  672    integer,           
intent(in)    :: rs_sign
 
  673    real(real64),      
intent(out)   :: magnetic_field(:,:)
 
  674    real(real64),   
optional, 
intent(in)    :: mu_field(:)
 
  675    integer, 
optional, 
intent(in)    :: np
 
  685    if (
present(mu_field)) 
then 
  687        magnetic_field(ip, :) = 
sqrt(
m_two*mu_field(ip)) * rs_sign * aimag(rs_state(ip, :))
 
  691        magnetic_field(ip, :) = 
sqrt(
m_two*
p_mu) * rs_sign * aimag(rs_state(ip, :))
 
  703    complex(real64), 
intent(in)    :: rs_current_element
 
  704    real(real64),    
intent(inout) :: current_element
 
  705    real(real64), 
optional, 
intent(in)    :: ep_element
 
  709    if (
present(ep_element)) 
then 
  710      current_element = 
sqrt(
m_two*ep_element) * real(rs_current_element, real64)
 
  712      current_element = 
sqrt(
m_two*
p_ep) * real(rs_current_element, real64)
 
  720    complex(real64), 
intent(in)    :: rs_current_vector(:)
 
  721    real(real64),    
intent(inout) :: current_vector(:)
 
  722    real(real64), 
optional, 
intent(in)    :: ep_element
 
  726    if (
present(ep_element)) 
then 
  727      current_vector(:) = 
sqrt(
m_two*ep_element) * real(rs_current_vector(:), real64)
 
  729      current_vector(:) = 
sqrt(
m_two*
p_ep) * real(rs_current_vector(:), real64)
 
  736  subroutine get_current_state(rs_current_field, current_field, mesh, ep_field, np)
 
  737    complex(real64),   
intent(in)    :: rs_current_field(:,:)
 
  738    real(real64),      
intent(inout) :: current_field(:,:)
 
  739    real(real64),   
optional, 
intent(in)    :: ep_field(:)
 
  740    class(
mesh_t),     
intent(in)    :: mesh
 
  741    integer, 
optional, 
intent(in)    :: np
 
  750      if (
present(ep_field)) 
then 
  751        current_field(ip, :) = 
sqrt(
m_two*ep_field(ip)) * real(rs_current_field(ip, :), real64)
 
  753        current_field(ip, :) = 
sqrt(
m_two*
p_ep) * real(rs_current_field(ip, :), real64)
 
  765    complex(real64),     
intent(inout)   :: rs_state_point(:,:)
 
  766    complex(real64),     
intent(in)      :: rs_state(:,:)
 
  767    real(real64),        
intent(in)      :: pos(:,:)
 
  769    class(
mesh_t),       
intent(in)      :: mesh
 
  771    integer :: ip, pos_index, rankmin
 
  773    complex(real64), 
allocatable :: ztmp(:)
 
  777    safe_allocate(ztmp(1:
size(rs_state, dim=2)))
 
  779    do ip = 1, st%selected_points_number
 
  781      if (mesh%mpi_grp%rank == rankmin) 
then 
  782        ztmp(:) = rs_state(pos_index, :)
 
  784      if (mesh%parallel_in_domains) 
then 
  785        call mesh%mpi_grp%bcast(ztmp, st%dim, mpi_double_complex, rankmin)
 
  787      rs_state_point(:, ip) = ztmp(:)
 
  790    safe_deallocate_a(ztmp)
 
  799    complex(real64), 
contiguous,   
intent(inout)   :: rs_state_point(:,:)
 
  800    type(
batch_t),       
intent(in)      :: rs_stateb
 
  802    class(
mesh_t),       
intent(in)      :: mesh
 
  805    complex(real64) :: rs_state_tmp(1:st%dim, 1:st%selected_points_number)
 
  808    integer(int64) :: localsize, dim3, dim2
 
  814    select case (rs_stateb%status())
 
  816      do ip_in = 1, st%selected_points_number
 
  817        ip = st%selected_points_map(ip_in)
 
  819          rs_state_tmp(1:st%dim, ip_in) = rs_stateb%zff_linear(ip, 1:st%dim)
 
  823      do ip_in = 1, st%selected_points_number
 
  824        ip = st%selected_points_map(ip_in)
 
  826          rs_state_tmp(1:st%dim, ip_in) = rs_stateb%zff_pack(1:st%dim, ip)
 
  833        st%selected_points_number*st%dim)
 
  848      call accel_kernel_run(kernel, (/rs_stateb%pack_size_real(1), dim2, dim3/), &
 
  849        (/rs_stateb%pack_size_real(1), localsize, 1_int64/))
 
  850      call accel_read_buffer(buff_points, st%selected_points_number*st%dim, rs_state_tmp)
 
  854    call mesh%mpi_grp%allreduce(rs_state_tmp, rs_state_point, st%selected_points_number*st%dim, mpi_double_complex, mpi_sum)
 
  861    type(
grid_t),      
intent(in)    :: gr
 
  862    real(real64), 
contiguous, 
intent(inout) :: field(:,:)
 
  863    real(real64), 
contiguous, 
intent(inout) :: field_div(:)
 
  864    logical,           
intent(in)    :: charge_density
 
  871      field_div = 
p_ep * field_div
 
  879  subroutine get_poynting_vector(mesh, st, rs_state, rs_sign, poynting_vector, ep_field, mu_field)
 
  880    class(
mesh_t),            
intent(in)    :: mesh
 
  882    complex(real64),          
intent(in)    :: rs_state(:,:)
 
  883    integer,                  
intent(in)    :: rs_sign
 
  884    real(real64),             
intent(out)   :: poynting_vector(:,:)
 
  885    real(real64),   
optional, 
intent(in)    :: ep_field(:)
 
  886    real(real64),   
optional, 
intent(in)    :: mu_field(:)
 
  891    if (
present(ep_field) .and. 
present(mu_field)) 
then 
  893        poynting_vector(ip, 1:3) = 
m_one/mu_field(ip) * 
sqrt(
m_two/ep_field(ip)) &
 
  896          rs_sign*aimag(rs_state(ip,1:3)))
 
  900        poynting_vector(ip, 1:3) = 
m_one/st%mu(ip) * 
sqrt(
m_two/st%ep(ip)) &
 
  903          rs_sign*aimag(rs_state(ip, 1:3)))
 
  913    class(
mesh_t),            
intent(in)    :: mesh
 
  915    integer,                  
intent(in)    :: rs_sign
 
  916    real(real64),             
intent(out)   :: poynting_vector(:,:)
 
  925        rs_sign*aimag(st%rs_state_plane_waves(ip,:)))
 
  934    type(
mesh_t),             
intent(in)    :: mesh
 
  936    real(real64),             
intent(in)    :: poynting_vector(:,:)
 
  937    real(real64),             
intent(out)   :: orbital_angular_momentum(:,:)
 
  944      orbital_angular_momentum(ip,1:3) = 
dcross_product(real(mesh%x(ip, 1:3), real64) , &
 
  945        poynting_vector(ip, 1:3))
 
  954    complex(real64), 
contiguous, 
intent(in)    :: rs_state(:, :)
 
  955    integer,           
intent(in)    :: np
 
  956    integer,           
intent(in)    :: dim
 
  957    integer, 
optional, 
intent(in)    :: offset
 
  959    integer :: offset_, idir
 
  965    do idir = offset_, offset_ + dim - 1
 
  974    type(
batch_t),     
intent(in)    :: rs_stateb
 
  975    complex(real64), 
contiguous, 
intent(out)   :: rs_state(:, :)
 
  976    integer,           
intent(in)    :: np
 
  977    integer,           
intent(in)    :: dim
 
  978    integer, 
optional, 
intent(in)    :: offset
 
  980    integer :: offset_, idir
 
  986    do idir = offset_, offset_ + dim - 1
 
 1004    select case (st%transverse_field_mode)
 
 1006      call helmholtz%get_trans_field(namespace, st%rs_state_trans, total_field=st%rs_state)
 
 1008      call helmholtz%get_long_field(namespace, st%rs_state_long, total_field=st%rs_state)
 
 1009      st%rs_state_trans = st%rs_state - st%rs_state_long
 
 1011      message(1) = 
'Unknown transverse field calculation mode.' 
There are several ways how to call batch_set_state and batch_get_state:
 
Prints out to iunit a message in the form: ["InputVariable" = value] where "InputVariable" is given b...
 
subroutine, public accel_kernel_start_call(this, file_name, kernel_name, flags)
 
integer, parameter, public accel_mem_read_write
 
integer pure function, public accel_max_size_per_dim(dim)
 
subroutine, public accel_release_buffer(this)
 
pure logical function, public accel_is_enabled()
 
integer function, public accel_kernel_workgroup_size(kernel)
 
This module implements batches of mesh functions.
 
integer, parameter, public batch_not_packed
functions are stored in CPU memory, unpacked order
 
integer, parameter, public batch_device_packed
functions are stored in device memory in packed order
 
integer, parameter, public batch_packed
functions are stored in CPU memory, in transposed (packed) order
 
This module implements common operations on batches of mesh functions.
 
This module provides the BLACS processor grid.
 
subroutine, public blacs_proc_grid_end(this)
 
This module calculates the derivatives (gradients, Laplacians, etc.) of a function.
 
subroutine, public dderivatives_div(der, ff, op_ff, ghost_update, set_bc, to_cartesian)
apply the divergence operator to a vector of mesh functions
 
subroutine, public distributed_end(this)
 
real(real64), parameter, public m_two
 
real(real64), parameter, public m_zero
 
real(real64), parameter, public p_mu
 
real(real64), parameter, public p_ep
 
complex(real64), parameter, public m_z0
 
complex(real64), parameter, public m_zi
 
real(real64), parameter, public m_one
 
This module implements the underlying real-space grid.
 
The Helmholtz decomposition is intended to contain "only mathematical" functions and procedures to co...
 
This module is intended to contain "only mathematical" functions and procedures.
 
pure real(real64) function, dimension(1:3), public dcross_product(a, b)
 
This module defines various routines, operating on mesh functions.
 
This module defines the meshes, which are used in Octopus.
 
integer function, public mesh_nearest_point(mesh, pos, dmin, rankmin)
Returns the index of the point which is nearest to a given vector position pos.
 
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
 
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
 
type(mpi_comm), parameter, public mpi_comm_undefined
used to indicate a communicator has not been initialized
 
subroutine mpi_grp_init(grp, comm)
Initialize MPI group instance.
 
This module handles the communicators for the various parallelization strategies.
 
integer function, public parse_block(namespace, name, blk, check_varinfo_)
 
subroutine, public profiling_out(label)
Increment out counter and sum up difference between entry and exit time.
 
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
 
This module handles spin dimensions of the states and the k-point distribution.
 
This module handles groups of electronic batches and their parallel distribution.
 
subroutine, public get_poynting_vector_plane_waves(mesh, st, rs_sign, poynting_vector)
 
subroutine, public mxll_set_batch(rs_stateb, rs_state, np, dim, offset)
 
integer, parameter, public transverse_from_helmholtz
 
subroutine, public get_rs_state_at_point(rs_state_point, rs_state, pos, st, mesh)
 
subroutine, public get_electric_field_vector(rs_state_vector, electric_field_vector, ep_element)
 
subroutine, public get_current_vector(rs_current_vector, current_vector, ep_element)
 
subroutine, public build_rs_vector(e_vector, b_vector, rs_sign, rs_vector, ep_element, mu_element)
 
subroutine, public get_transverse_rs_state(helmholtz, st, namespace)
 
subroutine, public build_rs_current_element(current_element, rs_current_element, ep_element)
 
subroutine, public mxll_get_batch(rs_stateb, rs_state, np, dim, offset)
 
subroutine, public get_rs_state_batch_selected_points(rs_state_point, rs_stateb, st, mesh)
 
subroutine, public states_mxll_end(st)
 
subroutine, public states_mxll_allocate(st, mesh)
Allocates the Maxwell states defined within a states_mxll_t structure.
 
subroutine, public build_rs_element(e_element, b_element, rs_sign, rs_element, ep_element, mu_element)
 
subroutine, public build_rs_current_state(current_state, mesh, rs_current_state, ep_field, np)
 
integer, parameter, public transverse_as_total_minus_long
 
subroutine, public build_rs_state(e_field, b_field, rs_sign, rs_state, mesh, ep_field, mu_field, np)
 
subroutine, public states_mxll_init(st, namespace, space)
 
subroutine, public get_current_element(rs_current_element, current_element, ep_element)
 
subroutine, public get_current_state(rs_current_field, current_field, mesh, ep_field, np)
 
subroutine, public get_electric_field_state(rs_state, mesh, electric_field, ep_field, np)
 
subroutine, public build_rs_current_vector(current_vector, rs_current_vector, ep_element)
 
subroutine, public get_poynting_vector(mesh, st, rs_state, rs_sign, poynting_vector, ep_field, mu_field)
 
subroutine, public get_divergence_field(gr, field, field_div, charge_density)
 
subroutine, public get_magnetic_field_state(rs_state, mesh, rs_sign, magnetic_field, mu_field, np)
 
subroutine, public get_magnetic_field_vector(rs_state_vector, rs_sign, magnetic_field_vector, mu_element)
 
subroutine, public get_orbital_angular_momentum(mesh, st, poynting_vector, orbital_angular_momentum)
 
type(type_t), public type_cmplx
 
type(type_t), public type_integer
 
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
 
This module defines the unit system, used for input and output.
 
type(unit_system_t), public units_inp
the units systems for reading and writing
 
Class defining batches of mesh functions.
 
Distribution of N instances over mpi_grpsize processes, for the local rank mpi_grprank....
 
Description of the grid, containing information on derivatives, stencil, and symmetries.
 
Describes mesh distribution to nodes.
 
This is defined even when running serial.