101 logical :: bc_add_ab_region = .false.
102 logical :: bc_zero = .false.
103 logical :: bc_constant = .false.
104 logical :: bc_mirror_pec = .false.
105 logical :: bc_mirror_pmc = .false.
106 logical :: bc_plane_waves = .false.
107 logical :: bc_medium = .false.
108 type(exponential_t) :: te
109 logical :: plane_waves_in_box
110 integer :: tr_etrs_approx
113 integer,
public,
parameter :: &
114 RS_TRANS_FORWARD = 1, &
117 integer,
parameter :: &
118 MXWLL_ETRS_FULL = 0, &
125 type(grid_t),
intent(in) :: gr
126 type(namespace_t),
intent(in) :: namespace
127 type(states_mxll_t),
intent(inout) :: st
128 type(hamiltonian_mxll_t),
intent(inout) :: hm
129 type(propagator_mxll_t),
intent(inout) :: tr
138 select case (hm%bc%bc_type(idim))
143 tr%bc_constant = .
true.
144 tr%bc_add_ab_region = .
true.
145 hm%bc_constant = .
true.
146 hm%bc_add_ab_region = .
true.
148 tr%bc_mirror_pec = .
true.
149 hm%bc_mirror_pec = .
true.
151 tr%bc_mirror_pmc = .
true.
152 hm%bc_mirror_pmc = .
true.
154 tr%bc_plane_waves = .
true.
155 tr%bc_add_ab_region = .
true.
156 hm%plane_waves = .
true.
157 hm%bc_plane_waves = .
true.
158 hm%bc_add_ab_region = .
true.
164 safe_allocate(st%rs_state_const(1:st%dim))
165 st%rs_state_const =
m_z0
179 call parse_variable(namespace,
'MaxwellTDETRSApprox', mxwll_etrs_full, tr%tr_etrs_approx)
189 call parse_variable(namespace,
'MaxwellPlaneWavesInBox', .false., tr%plane_waves_in_box)
190 if (tr%plane_waves_in_box .and. .not. hm%bc%do_plane_waves)
then
205 subroutine mxll_propagation_step(hm, namespace, gr, space, st, tr, rs_stateb, ff_rs_inhom_t1, ff_rs_inhom_t2, time, dt)
207 type(namespace_t),
intent(in) :: namespace
209 class(
space_t),
intent(in) :: space
213 complex(real64),
contiguous,
intent(in) :: ff_rs_inhom_t1(:,:)
214 complex(real64),
contiguous,
intent(in) :: ff_rs_inhom_t2(:,:)
215 real(real64),
intent(in) :: time
216 real(real64),
intent(in) :: dt
218 integer :: ii, ff_dim, idim, istate, inter_steps
219 real(real64) :: inter_dt, inter_time
223 type(
batch_t) :: ff_rs_stateb, ff_rs_state_pmlb
224 type(
batch_t) :: ff_rs_inhom_1b, ff_rs_inhom_2b, ff_rs_inhom_meanb
225 complex(real64),
allocatable :: rs_state(:, :)
232 if (hm%ma_mx_coupling_apply)
then
233 message(1) =
"Maxwell-matter coupling not implemented yet"
236 safe_allocate(rs_state(gr%np, st%dim))
238 if (tr%plane_waves_in_box)
then
242 safe_deallocate_a(rs_state)
248 if (hm%bc%bc_ab_type(idim) == option__maxwellabsorbingboundaries__cpml)
then
254 if (pml_check .and. .not. hm%bc%pml%parameters_initialized) &
261 inter_dt =
m_one / inter_steps * dt
263 call zbatch_init(ff_rs_stateb, 1, 1, hm%dim, gr%np_part)
264 if (st%pack_states)
call ff_rs_stateb%do_pack(copy=.false.)
267 call ff_rs_stateb%copy_to(ff_rs_state_pmlb)
271 if ((hm%ma_mx_coupling_apply .or. hm%current_density_ext_flag .or. hm%current_density_from_medium) .and. &
273 call ff_rs_stateb%copy_to(ff_rs_inhom_1b)
274 call ff_rs_stateb%copy_to(ff_rs_inhom_2b)
275 call ff_rs_stateb%copy_to(ff_rs_inhom_meanb)
277 do istate = 1, hm%dim
278 call batch_set_state(ff_rs_inhom_meanb, istate, gr%np, ff_rs_inhom_t1(:, istate))
279 call batch_set_state(ff_rs_inhom_2b, istate, gr%np, ff_rs_inhom_t2(:, istate))
285 call ff_rs_inhom_meanb%copy_data_to(gr%np, ff_rs_inhom_1b)
286 call ff_rs_inhom_meanb%copy_data_to(gr%np, ff_rs_inhom_2b)
289 hm%cpml_hamiltonian = .false.
290 call tr%te%apply_batch(namespace, gr, hm, ff_rs_inhom_2b, inter_dt)
294 call ff_rs_inhom_meanb%copy_data_to(gr%np, ff_rs_inhom_2b)
296 call tr%te%apply_batch(namespace, gr, hm, ff_rs_inhom_2b, inter_dt*
m_half)
299 call ff_rs_inhom_meanb%copy_data_to(gr%np, ff_rs_inhom_2b)
301 call tr%te%apply_batch(namespace, gr, hm, ff_rs_inhom_2b, -inter_dt*
m_half)
304 call ff_rs_inhom_2b%end()
305 call ff_rs_inhom_meanb%end()
308 do ii = 1, inter_steps
311 inter_time = time + inter_dt * (ii-1)
322 hm%cpml_hamiltonian = pml_check
323 call tr%te%apply_batch(namespace, gr, hm, ff_rs_stateb, dt)
324 hm%cpml_hamiltonian = .false.
331 if ((hm%ma_mx_coupling_apply) .or. hm%current_density_ext_flag .or. hm%current_density_from_medium)
then
332 if (tr%tr_etrs_approx == mxwll_etrs_full)
then
333 call ff_rs_stateb%copy_to(ff_rs_inhom_1b)
334 call ff_rs_stateb%copy_to(ff_rs_inhom_2b)
335 call ff_rs_stateb%copy_to(ff_rs_inhom_meanb)
338 do istate = 1, hm%dim
339 call batch_set_state(ff_rs_inhom_meanb, istate, gr%np, ff_rs_inhom_t2(:, istate))
340 call batch_set_state(ff_rs_inhom_1b, istate, gr%np, ff_rs_inhom_t1(:, istate))
344 call ff_rs_inhom_1b%copy_data_to(gr%np, ff_rs_inhom_2b)
345 call batch_axpy(gr%np, ii / real(inter_steps, real64) , ff_rs_inhom_meanb, ff_rs_inhom_2b)
346 call batch_axpy(gr%np, (ii-1) / real(inter_steps, real64) , ff_rs_inhom_meanb, ff_rs_inhom_1b)
348 hm%cpml_hamiltonian = .false.
349 call tr%te%apply_batch(namespace, gr, hm, ff_rs_inhom_1b, inter_dt)
354 do istate = 1, hm%dim
355 call batch_set_state(ff_rs_inhom_1b, istate, gr%np, ff_rs_inhom_t1(:, istate))
356 call batch_set_state(ff_rs_inhom_2b, istate, gr%np, ff_rs_inhom_t2(:, istate))
360 call ff_rs_inhom_1b%copy_data_to(gr%np, ff_rs_inhom_2b)
362 call tr%te%apply_batch(namespace, gr, hm, ff_rs_inhom_1b, inter_dt/
m_two)
363 call tr%te%apply_batch(namespace, gr, hm, ff_rs_inhom_2b, -inter_dt/
m_two)
369 call ff_rs_inhom_1b%end()
370 call ff_rs_inhom_2b%end()
371 call ff_rs_inhom_meanb%end()
385 if (tr%bc_constant)
then
388 if (st%rs_state_const_external)
then
410 if (any(hm%bc%bc_ab_type == option__maxwellabsorbingboundaries__mask))
then
416 if (tr%bc_plane_waves)
then
425 if (tr%tr_etrs_approx == option__maxwelltdetrsapprox__const_steps)
then
426 call ff_rs_inhom_1b%end()
429 call ff_rs_stateb%end()
432 call ff_rs_state_pmlb%end()
435 safe_deallocate_a(rs_state)
445 type(namespace_t),
intent(in) :: namespace
446 type(
grid_t),
intent(inout) :: gr
449 real(real64),
intent(in) :: time
450 real(real64),
intent(in) :: dt
451 integer,
intent(in) :: counter
457 call st%rs_stateb%copy_to(rs_state_tmpb)
465 if (any(hm%bc%bc_ab_type == option__maxwellabsorbingboundaries__cpml))
then
476 if (counter == 0)
then
478 call batch_xpay(gr%np, st%rs_stateb, dt, rs_state_tmpb)
485 call st%rs_stateb%copy_data_to(gr%np, st%rs_state_prevb)
487 call rs_state_tmpb%copy_data_to(gr%np, st%rs_stateb)
490 if (any(hm%bc%bc_ab_type == option__maxwellabsorbingboundaries__cpml))
then
494 call rs_state_tmpb%end()
515 type(namespace_t),
intent(in) :: namespace
516 type(
grid_t),
intent(inout) :: gr
519 real(real64),
intent(in) :: time
520 real(real64),
intent(in) :: dt
526 call st%rs_stateb%copy_to(rs_state_tmpb)
535 if (any(hm%bc%bc_ab_type == option__maxwellabsorbingboundaries__cpml))
then
541 call hm%zapply(namespace, gr, st%rs_stateb, rs_state_tmpb)
544 if (hm%current_density_ext_flag .or. hm%current_density_from_medium)
then
546 call mxll_set_batch(st%inhomogeneousb, st%rs_current_density_t1, gr%np, st%dim)
551 call tr%te%apply_phi_batch(namespace, gr, hm, rs_state_tmpb, dt, 1)
553 call batch_axpy(gr%np, dt, rs_state_tmpb, st%rs_stateb)
555 call rs_state_tmpb%end()
558 if (any(hm%bc%bc_ab_type == option__maxwellabsorbingboundaries__cpml))
then
587 type(namespace_t),
intent(in) :: namespace
588 type(
grid_t),
intent(inout) :: gr
591 real(real64),
intent(in) :: time
592 real(real64),
intent(in) :: dt
598 call st%rs_stateb%copy_to(rs_state_tmpb)
607 if (any(hm%bc%bc_ab_type == option__maxwellabsorbingboundaries__cpml))
then
613 call hm%zapply(namespace, gr, st%rs_stateb, rs_state_tmpb)
616 if (hm%current_density_ext_flag .or. hm%current_density_from_medium)
then
618 call mxll_set_batch(st%inhomogeneousb, st%rs_current_density_t1, gr%np, st%dim)
622 call mxll_set_batch(st%inhomogeneousb, st%rs_current_density_t2, gr%np, st%dim)
627 call tr%te%apply_phi_batch(namespace, gr, hm, rs_state_tmpb, dt, 1)
629 call batch_axpy(gr%np, dt, rs_state_tmpb, st%rs_stateb)
630 if (hm%current_density_ext_flag .or. hm%current_density_from_medium)
then
633 call mxll_set_batch(st%inhomogeneousb, st%rs_current_density_t1, gr%np, st%dim)
637 call mxll_set_batch(st%inhomogeneousb, st%rs_current_density_t2, gr%np, st%dim)
642 call tr%te%apply_phi_batch(namespace, gr, hm, rs_state_tmpb, dt, 2)
644 call batch_axpy(gr%np, dt, rs_state_tmpb, st%rs_stateb)
648 if (any(hm%bc%bc_ab_type == option__maxwellabsorbingboundaries__cpml))
then
652 call rs_state_tmpb%end()
660 type(
grid_t),
intent(in) :: gr
663 integer :: ip, ip_in, il, idim
667 assert(
allocated(st%ep) .and.
allocated(st%mu))
671 if (hm%calc_medium_box)
then
672 do il = 1,
size(hm%medium_boxes)
673 assert(.not. hm%medium_boxes(il)%has_mapping)
674 do ip = 1, hm%medium_boxes(il)%points_number
675 if (abs(hm%medium_boxes(il)%c(ip)) <=
m_epsilon) cycle
676 st%ep(ip) = hm%medium_boxes(il)%ep(ip)
677 st%mu(ip) = hm%medium_boxes(il)%mu(ip)
684 do ip_in = 1, hm%bc%medium(idim)%points_number
685 ip = hm%bc%medium(idim)%points_map(ip_in)
686 st%ep(ip) = hm%bc%medium(idim)%ep(ip_in)
687 st%mu(ip) = hm%bc%medium(idim)%mu(ip_in)
700 type(
grid_t),
intent(in) :: gr
702 type(
batch_t),
intent(inout) :: rs_stateb
703 type(
batch_t),
intent(inout) :: ff_rs_stateb
704 integer,
intent(in) :: sign
706 complex(real64),
allocatable :: rs_state(:,:)
707 complex(real64),
allocatable :: rs_state_tmp(:,:)
717 safe_allocate(rs_state(1:gr%np, 1:st%dim))
720 if (sign == rs_trans_forward)
then
729 safe_allocate(rs_state_tmp(1:gr%np, 1:st%dim))
733 rs_state(1:np, ii) =
m_half * (rs_state(1:np, ii) + conjg(rs_state_tmp(1:np, ii)))
736 safe_deallocate_a(rs_state_tmp)
739 if (sign == rs_trans_forward)
then
740 call rs_stateb%copy_data_to(gr%np, ff_rs_stateb)
742 call ff_rs_stateb%copy_data_to(gr%np, rs_stateb)
745 safe_deallocate_a(rs_state)
756 class(
mesh_t),
intent(in) :: mesh
757 complex(real64),
intent(inout) :: rs_current_density(:,:)
758 complex(real64),
intent(inout) :: ff_density(:,:)
759 integer,
intent(in) :: sign
763 assert(
size(rs_current_density, dim=2) == 3)
770 if (sign == rs_trans_forward)
then
776 if (sign == rs_trans_forward)
then
777 ff_density(1:mesh%np, 1:3) = rs_current_density(1:mesh%np, 1:3)
779 rs_current_density(1:mesh%np, 1:3) = ff_density(1:mesh%np, 1:3)
791 class(
mesh_t),
intent(in) :: mesh
792 complex(real64),
intent(in) :: rs_current_density(:,:)
793 complex(real64),
intent(inout) :: rs_density_6x6(:,:)
797 assert(
size(rs_current_density, dim=2) == 3)
798 assert(
size(rs_density_6x6, dim=2) == 6)
802 rs_density_6x6(1:mesh%np, ii) = rs_current_density(1:mesh%np, ii)
803 rs_density_6x6(1:mesh%np, ii+3) = rs_current_density(1:mesh%np, ii)
810 class(
mesh_t),
intent(in) :: mesh
811 complex(real64),
intent(in) :: rs_density_6x6(:,:)
812 complex(real64),
intent(inout) :: rs_current_density(:,:)
816 assert(
size(rs_current_density, dim=2) == 3)
817 assert(
size(rs_density_6x6, dim=2) == 6)
821 rs_current_density(1:mesh%np, ii) =
m_half * &
822 real(rs_density_6x6(1:mesh%np, ii) + rs_density_6x6(1:mesh%np, ii+3), real64)
829 type(
grid_t),
intent(in) :: gr_mxll
832 type(
grid_t),
intent(in) :: gr_elec
834 complex(real64),
intent(inout) :: rs_state_matter(:,:)
836 complex(real64),
allocatable :: tmp_pot_mx_gr(:,:), tmp_grad_mx_gr(:,:)
838 safe_allocate(tmp_pot_mx_gr(1:gr_mxll%np_part,1))
839 safe_allocate(tmp_grad_mx_gr(1:gr_mxll%np,1:gr_mxll%box%dim))
844 tmp_pot_mx_gr(:,:) =
m_zero
845 tmp_grad_mx_gr(:,:) =
m_zero
846 call zderivatives_grad(gr_mxll%der, tmp_pot_mx_gr(:,1), tmp_grad_mx_gr(:,:), set_bc = .false.)
847 tmp_grad_mx_gr = - tmp_grad_mx_gr
850 call build_rs_state(real(tmp_grad_mx_gr(1:gr_mxll%np,:)), aimag(tmp_grad_mx_gr(1:gr_mxll%np,:)), st_mxll%rs_sign, &
851 rs_state_matter(1:gr_mxll%np,:), gr_mxll, st_mxll%ep(1:gr_mxll%np), st_mxll%mu(1:gr_mxll%np), &
854 safe_deallocate_a(tmp_pot_mx_gr)
855 safe_deallocate_a(tmp_grad_mx_gr)
862 poisson_solver, helmholtz, field, transverse_field, vector_potential)
863 type(namespace_t),
intent(in) :: namespace
864 type(
grid_t),
intent(in) :: gr_mxll
869 type(
poisson_t),
intent(in) :: poisson_solver
871 complex(real64),
intent(inout) :: field(:,:)
872 complex(real64),
intent(inout) :: transverse_field(:,:)
873 real(real64),
intent(inout) :: vector_potential(:,:)
876 complex(real64),
allocatable :: rs_state_plane_waves(:, :)
880 transverse_field =
m_z0
885 if (hm_mxll%ma_mx_coupling)
then
890 if (tr_mxll%bc_plane_waves .and. hm_mxll%plane_waves_apply)
then
891 safe_allocate(rs_state_plane_waves(1:gr_mxll%np, 1:st_mxll%dim))
892 call mxll_get_batch(st_mxll%rs_state_plane_wavesb, rs_state_plane_waves, gr_mxll%np, st_mxll%dim)
896 if (tr_mxll%bc_plane_waves .and. hm_mxll%plane_waves_apply)
then
897 transverse_field(1:np,:) = field(1:np,:) - rs_state_plane_waves(1:np,:)
899 transverse_field(1:np,:) = field(1:np,:)
902 call helmholtz%get_trans_field(namespace, transverse_field, total_field=field)
904 if (tr_mxll%bc_plane_waves .and. hm_mxll%plane_waves_apply)
then
905 transverse_field(1:np,:) = transverse_field(1:np,:) + rs_state_plane_waves(1:np,:)
906 safe_deallocate_a(rs_state_plane_waves)
911 transverse_field(1:np,:) = field
922 type(namespace_t),
intent(in) :: namespace
924 type(
grid_t),
intent(in) :: gr
926 complex(real64),
intent(in) :: field(:,:)
927 real(real64),
contiguous,
intent(inout) :: vector_potential(:,:)
930 real(real64),
allocatable :: dtmp(:,:)
932 safe_allocate(dtmp(1:gr%np_part,1:3))
937 dtmp = vector_potential
940 call dpoisson_solve(poisson_solver, namespace, dtmp(:,idim), vector_potential(:,idim), .
true.)
944 safe_deallocate_a(dtmp)
949 subroutine energy_mxll_calc(gr, st, hm, energy_mxll, rs_field, rs_field_plane_waves)
950 type(
grid_t),
intent(in) :: gr
954 complex(real64),
intent(in) :: rs_field(:,:)
955 complex(real64),
optional,
intent(in) :: rs_field_plane_waves(:,:)
957 real(real64),
allocatable :: energy_density(:), e_energy_density(:), b_energy_density(:), energy_density_plane_waves(:)
963 safe_allocate(energy_density(1:gr%np))
964 safe_allocate(e_energy_density(1:gr%np))
965 safe_allocate(b_energy_density(1:gr%np))
966 if (
present(rs_field_plane_waves) .and. hm%plane_waves)
then
967 safe_allocate(energy_density_plane_waves(1:gr%np))
971 b_energy_density, hm%plane_waves, rs_field_plane_waves, energy_density_plane_waves)
972 energy_mxll%energy =
dmf_integrate(gr, energy_density, mask=st%inner_points_mask)
973 energy_mxll%e_energy =
dmf_integrate(gr, e_energy_density, mask=st%inner_points_mask)
974 energy_mxll%b_energy =
dmf_integrate(gr, b_energy_density, mask=st%inner_points_mask)
975 if (
present(rs_field_plane_waves) .and. hm%plane_waves)
then
976 energy_mxll%energy_plane_waves =
dmf_integrate(gr, energy_density_plane_waves, mask=st%inner_points_mask)
978 energy_mxll%energy_plane_waves =
m_zero
981 energy_mxll%boundaries =
dmf_integrate(gr, energy_density, mask=st%boundary_points_mask)
983 safe_deallocate_a(energy_density)
984 safe_deallocate_a(e_energy_density)
985 safe_deallocate_a(b_energy_density)
986 if (
present(rs_field_plane_waves) .and. hm%plane_waves)
then
987 safe_deallocate_a(energy_density_plane_waves)
997 type(
grid_t),
intent(in) :: gr
1001 type(
batch_t),
intent(in) :: rs_fieldb
1002 type(
batch_t),
intent(in) :: rs_field_plane_wavesb
1004 type(
batch_t) :: e_fieldb, b_fieldb, e_field_innerb, b_field_innerb, rs_field_plane_waves_innerb
1005 real(real64) :: tmp(1:st%dim)
1006 complex(real64) :: ztmp(1:st%dim)
1013 if (st%pack_states)
then
1014 call e_fieldb%do_pack(copy=.false.)
1016 call e_fieldb%copy_to(b_fieldb)
1017 call e_fieldb%copy_to(e_field_innerb)
1018 call e_fieldb%copy_to(b_field_innerb)
1026 call batch_copy_with_map(st%inner_points_number, st%buff_inner_points_map, e_fieldb, e_field_innerb)
1027 call batch_copy_with_map(st%inner_points_number, st%buff_inner_points_map, b_fieldb, b_field_innerb)
1029 call batch_copy_with_map(st%inner_points_number, st%inner_points_map, e_fieldb, e_field_innerb)
1030 call batch_copy_with_map(st%inner_points_number, st%inner_points_map, b_fieldb, b_field_innerb)
1033 energy_mxll%e_energy = sum(tmp)
1035 energy_mxll%b_energy = sum(tmp)
1036 energy_mxll%energy = energy_mxll%e_energy + energy_mxll%b_energy
1039 energy_mxll%boundaries = sum(tmp)
1041 energy_mxll%boundaries = energy_mxll%boundaries + sum(tmp)
1042 energy_mxll%boundaries = energy_mxll%boundaries - energy_mxll%energy
1044 if (hm%plane_waves)
then
1045 call rs_field_plane_wavesb%copy_to(rs_field_plane_waves_innerb)
1049 rs_field_plane_wavesb, rs_field_plane_waves_innerb)
1052 rs_field_plane_wavesb, rs_field_plane_waves_innerb)
1055 energy_mxll%energy_plane_waves = sum(real(tmp, real64) )
1056 call rs_field_plane_waves_innerb%end()
1058 energy_mxll%energy_plane_waves =
m_zero
1063 call e_field_innerb%end()
1064 call b_field_innerb%end()
1073 type(namespace_t),
intent(in) :: namespace
1074 type(
grid_t),
intent(in) :: gr
1078 real(real64),
intent(in) :: time
1079 real(real64),
intent(in) :: dt
1080 real(real64),
intent(in) :: time_delay
1081 complex(real64),
intent(inout) :: rs_state(:,:)
1083 integer :: ip, ip_in, idim
1084 logical :: mask_check
1089 mask_check = .false.
1092 if (hm%bc%bc_ab_type(idim) == option__maxwellabsorbingboundaries__mask)
then
1097 if (mask_check)
then
1098 if (tr%bc_plane_waves .and. hm%plane_waves_apply)
then
1100 call mxll_get_batch(st%rs_state_plane_wavesb, st%rs_state_plane_waves, gr%np, st%dim)
1101 rs_state = rs_state - st%rs_state_plane_waves
1103 rs_state = rs_state + st%rs_state_plane_waves
1104 else if (tr%bc_constant .and. hm%spatial_constant_apply)
then
1107 do ip_in=1, hm%bc%constant_points_number
1108 ip = hm%bc%constant_points_map(ip_in)
1109 rs_state(ip,:) = rs_state(ip,:) - st%rs_state_const(:)
1112 do ip_in=1, hm%bc%constant_points_number
1113 ip = hm%bc%constant_points_map(ip_in)
1114 rs_state(ip,:) = rs_state(ip,:) + st%rs_state_const(:)
1129 complex(real64),
intent(inout) :: rs_state(:,:)
1131 integer :: ip, ip_in, idim
1138 if (hm%bc%bc_ab_type(idim) == option__maxwellabsorbingboundaries__mask)
then
1139 do ip_in = 1, hm%bc%mask_points_number(idim)
1140 ip = hm%bc%mask_points_map(ip_in,idim)
1141 rs_state(ip,:) = rs_state(ip,:) * hm%bc%mask(ip_in,idim)
1154 type(
grid_t),
intent(in) :: gr
1157 type(
batch_t),
intent(in) :: ff_rs_stateb
1158 type(
batch_t),
intent(inout) :: ff_rs_state_pmlb
1161 complex(real64),
allocatable :: rs_state_constant(:,:)
1162 type(
batch_t) :: rs_state_constantb
1168 if (tr%bc_plane_waves .and. hm%plane_waves_apply)
then
1170 ff_rs_state_pmlb, rs_trans_forward)
1172 else if (tr%bc_constant .and. hm%spatial_constant_apply)
then
1177 safe_allocate(rs_state_constant(1:gr%np,1:3))
1179 rs_state_constant(1:gr%np, ii) = st%rs_state_const(ii)
1181 call ff_rs_stateb%copy_to(rs_state_constantb)
1182 call mxll_set_batch(rs_state_constantb, rs_state_constant, gr%np, 3)
1185 ff_rs_state_pmlb, rs_trans_forward)
1188 call rs_state_constantb%end()
1190 safe_deallocate_a(rs_state_constant)
1193 call ff_rs_stateb%copy_data_to(gr%np, ff_rs_state_pmlb)
1204 type(namespace_t),
intent(in) :: namespace
1205 type(
grid_t),
intent(in) :: gr
1208 real(real64),
intent(in) :: time
1209 real(real64),
intent(in) :: dt
1210 real(real64),
intent(in) :: time_delay
1211 type(
batch_t),
intent(inout) :: ff_rs_state_pmlb
1212 type(
batch_t),
intent(inout) :: ff_rs_stateb
1214 integer :: ii, ff_dim
1215 complex(real64),
allocatable :: rs_state_constant(:,:), ff_rs_state_constant(:,:)
1216 type(
batch_t) :: ff_rs_state_plane_wavesb, ff_rs_constantb, rs_state_constantb
1222 if (tr%bc_plane_waves .and. hm%plane_waves_apply)
then
1223 hm%cpml_hamiltonian = .
true.
1224 call tr%te%apply_batch(namespace, gr, hm, ff_rs_state_pmlb, dt)
1225 hm%cpml_hamiltonian = .false.
1228 call ff_rs_stateb%copy_to(ff_rs_state_plane_wavesb)
1234 ff_rs_state_pmlb, ff_rs_state_plane_wavesb, ff_rs_stateb)
1236 call batch_add_with_map(hm%bc%plane_wave%points_number, hm%bc%plane_wave%points_map, &
1237 ff_rs_state_pmlb, ff_rs_state_plane_wavesb, ff_rs_stateb)
1240 call ff_rs_state_plane_wavesb%end()
1242 else if (tr%bc_constant .and. hm%spatial_constant_apply)
then
1243 hm%cpml_hamiltonian = .
true.
1244 call tr%te%apply_batch(namespace, gr, hm, ff_rs_state_pmlb, dt)
1245 hm%cpml_hamiltonian = .false.
1247 call ff_rs_stateb%copy_to(ff_rs_constantb)
1248 ff_dim = ff_rs_stateb%nst_linear
1249 safe_allocate(rs_state_constant(1:gr%np, 1:st%dim))
1253 rs_state_constant(1:gr%np, ii) = st%rs_state_const(ii)
1255 call ff_rs_stateb%copy_to(rs_state_constantb)
1256 call mxll_set_batch(rs_state_constantb, rs_state_constant, gr%np, st%dim)
1261 call batch_add_with_map(hm%bc%constant_points_number, hm%bc%buff_constant_points_map, &
1262 ff_rs_state_pmlb, ff_rs_constantb, ff_rs_stateb)
1265 ff_rs_state_pmlb, ff_rs_constantb, ff_rs_stateb)
1268 call ff_rs_constantb%end()
1269 call rs_state_constantb%end()
1271 safe_deallocate_a(rs_state_constant)
1272 safe_deallocate_a(ff_rs_state_constant)
1283 type(
grid_t),
intent(in) :: gr
1284 type(
batch_t),
intent(inout) :: ff_rs_state_pmlb
1301 type(
grid_t),
intent(in) :: gr
1302 type(
batch_t),
intent(inout) :: ff_rs_state_pmlb
1304 integer :: ip, ip_in, np_part, rs_sign
1305 complex(real64) :: pml_a, pml_b, pml_g, grad
1306 integer :: pml_dir, field_dir, ifield, idir
1307 integer,
parameter :: field_dirs(3, 2) = reshape([2, 3, 1, 3, 1, 2], [3, 2])
1308 logical :: with_medium
1309 type(
batch_t) :: gradb(gr%der%dim)
1311 integer :: bsize, gsize
1317 assert(hm%dim == 3 .or. hm%dim == 6)
1319 np_part = gr%np_part
1320 rs_sign = hm%rs_sign
1324 with_medium = hm%dim == 6
1326 do pml_dir = 1, hm%st%dim
1327 select case (gradb(pml_dir)%status())
1329 do ip_in=1, hm%bc%pml%points_number
1330 ip = hm%bc%pml%points_map(ip_in)
1331 pml_a = hm%bc%pml%a(ip_in,pml_dir)
1332 pml_b = hm%bc%pml%b(ip_in,pml_dir)
1334 field_dir = field_dirs(pml_dir, ifield)
1335 grad = gradb(pml_dir)%zff_linear(ip, field_dir)
1336 pml_g = hm%bc%pml%conv_plus(ip_in, pml_dir, field_dir)
1337 hm%bc%pml%conv_plus(ip_in, pml_dir, field_dir) = pml_a * grad + pml_b * pml_g
1338 if (with_medium)
then
1339 grad = gradb(pml_dir)%zff_linear(ip, field_dir+3)
1340 pml_g = hm%bc%pml%conv_minus(ip_in, pml_dir, field_dir)
1341 hm%bc%pml%conv_minus(ip_in, pml_dir, field_dir) = pml_a * grad + pml_b * pml_g
1346 do ip_in=1, hm%bc%pml%points_number
1347 ip = hm%bc%pml%points_map(ip_in)
1348 pml_a = hm%bc%pml%a(ip_in,pml_dir)
1349 pml_b = hm%bc%pml%b(ip_in,pml_dir)
1351 field_dir = field_dirs(pml_dir, ifield)
1352 grad = gradb(pml_dir)%zff_pack(field_dir, ip)
1353 pml_g = hm%bc%pml%conv_plus(ip_in, pml_dir, field_dir)
1354 hm%bc%pml%conv_plus(ip_in, pml_dir, field_dir) = pml_a * grad + pml_b * pml_g
1355 if (with_medium)
then
1356 grad = gradb(pml_dir)%zff_pack(field_dir+3, ip)
1357 pml_g = hm%bc%pml%conv_minus(ip_in, pml_dir, field_dir)
1358 hm%bc%pml%conv_minus(ip_in, pml_dir, field_dir) = pml_a * grad + pml_b * pml_g
1365 if (with_medium)
then
1388 do idir = 1, gr%der%dim
1389 call gradb(idir)%end()
1404 type(namespace_t),
intent(in) :: namespace
1408 integer :: il, nlines, idim, ncols, ierr
1409 real(real64) :: e_field(st%dim), b_field(st%dim)
1410 character(len=1024) :: mxf_expression
1433 if (
parse_block(namespace,
'UserDefinedConstantSpatialMaxwellField', blk) == 0)
then
1434 st%rs_state_const_external = .
true.
1436 safe_allocate(st%rs_state_const_td_function(1:nlines))
1437 safe_allocate(st%rs_state_const_amp(1:st%dim, 1:nlines))
1444 if (ncols /= 7)
then
1445 message(1) =
'Each line in the UserDefinedConstantSpatialMaxwellField block must have'
1456 call build_rs_vector(e_field, b_field, st%rs_sign, st%rs_state_const_amp(:,il))
1457 call tdf_read(st%rs_state_const_td_function(il), namespace, trim(mxf_expression), ierr)
1471 call parse_variable(namespace,
'PropagateSpatialMaxwellField', .
true., hm%spatial_constant_propagate)
1480 logical,
intent(in) :: constant_calc
1482 type(
grid_t),
intent(in) :: gr
1484 real(real64),
intent(in) :: time
1485 real(real64),
intent(in) :: dt
1486 real(real64),
intent(in) :: delay
1487 complex(real64),
intent(inout) :: rs_state(:,:)
1488 logical,
optional,
intent(in) :: set_initial_state
1490 integer :: ip, ic, icn
1491 real(real64) :: tf_old, tf_new
1492 logical :: set_initial_state_
1498 set_initial_state_ = .false.
1499 if (
present(set_initial_state)) set_initial_state_ = set_initial_state
1501 if (hm%spatial_constant_apply)
then
1502 if (constant_calc)
then
1503 icn =
size(st%rs_state_const_td_function(:))
1504 st%rs_state_const(:) =
m_z0
1506 tf_old =
tdf(st%rs_state_const_td_function(ic), time-delay-dt)
1507 tf_new =
tdf(st%rs_state_const_td_function(ic), time-delay)
1509 if (set_initial_state_ .or. (.not. hm%spatial_constant_propagate))
then
1510 rs_state(ip,:) = st%rs_state_const_amp(:,ic) * tf_new
1512 rs_state(ip,:) = rs_state(ip,:) + st%rs_state_const_amp(:,ic) * (tf_new - tf_old)
1515 st%rs_state_const(:) = st%rs_state_const(:) + st%rs_state_const_amp(:, ic)
1517 st%rs_state_const(:) = st%rs_state_const(:) * tf_new
1528 logical,
intent(in) :: constant_calc
1532 complex(real64),
intent(inout) :: rs_state(:,:)
1534 integer :: ip_in, ip
1539 if (hm%spatial_constant_apply)
then
1540 if (constant_calc)
then
1541 do ip_in = 1, bc%constant_points_number
1542 ip = bc%constant_points_map(ip_in)
1543 rs_state(ip,:) = st%rs_state_const(:)
1544 bc%constant_rs_state(ip_in,:) = st%rs_state_const(:)
1558 complex(real64),
intent(inout) :: rs_state(:,:)
1560 integer :: ip, ip_in, idim
1561 real(real64) :: e_field(st%dim), b_field(st%dim)
1567 do ip_in = 1, bc%mirror_points_number(idim)
1568 ip = bc%mirror_points_map(ip_in, idim)
1571 call build_rs_vector(e_field(:), b_field(:), st%rs_sign, rs_state(ip,:), st%ep(ip), st%mu(ip))
1583 complex(real64),
intent(inout) :: rs_state(:,:)
1585 integer :: ip, ip_in, idim
1586 real(real64) :: e_field(st%dim), b_field(st%dim)
1592 do ip_in = 1, bc%mirror_points_number(idim)
1593 ip = bc%mirror_points_map(ip_in,idim)
1596 call build_rs_vector(e_field(:), b_field(:), st%rs_sign, rs_state(ip,:), st%ep(ip), st%mu(ip))
1608 class(
mesh_t),
intent(in) :: mesh
1609 real(real64),
intent(in) :: time
1610 real(real64),
intent(in) :: time_delay
1611 complex(real64),
intent(inout) :: rs_state(:,:)
1613 integer :: ip, ip_in, wn
1614 real(real64) :: x_prop(mesh%box%dim), rr, vv(mesh%box%dim), k_vector(mesh%box%dim)
1615 real(real64) :: k_vector_abs, nn
1616 complex(real64) :: e0(mesh%box%dim)
1617 real(real64) :: e_field(mesh%box%dim), b_field(mesh%box%dim)
1618 complex(real64) :: rs_state_add(mesh%box%dim)
1619 complex(real64) :: mx_func
1625 if (hm%plane_waves_apply)
then
1626 do wn = 1, hm%bc%plane_wave%number
1627 k_vector(:) = hm%bc%plane_wave%k_vector(1:mesh%box%dim, wn)
1628 k_vector_abs = norm2(k_vector(1:mesh%box%dim))
1629 vv(:) = hm%bc%plane_wave%v_vector(1:mesh%box%dim, wn)
1630 e0(:) = hm%bc%plane_wave%e_field(1:mesh%box%dim, wn)
1631 do ip_in = 1, hm%bc%plane_wave%points_number
1632 ip = hm%bc%plane_wave%points_map(ip_in)
1633 if (wn == 1) rs_state(ip,:) =
m_z0
1635 x_prop(1:mesh%box%dim) = mesh%x(1:mesh%box%dim, ip) - vv(1:mesh%box%dim) * (time - time_delay)
1636 rr = norm2(x_prop(1:mesh%box%dim))
1637 if (hm%bc%plane_wave%modus(wn) == option__maxwellincidentwaves__plane_wave_mx_function)
then
1639 mx_func =
mxf(hm%bc%plane_wave%mx_function(wn), x_prop(1:mesh%box%dim))
1640 e_field(1:mesh%box%dim) = real(e0(1:mesh%box%dim) * mx_func, real64)
1643 call build_rs_vector(e_field, b_field, st%rs_sign, rs_state_add, st%ep(ip), st%mu(ip))
1644 rs_state(ip, :) = rs_state(ip, :) + rs_state_add(:)
1648 do ip_in = 1, hm%bc%plane_wave%points_number
1649 ip = hm%bc%plane_wave%points_map(ip_in)
1663 type(namespace_t),
intent(in) :: namespace
1665 type(
grid_t),
intent(in) :: gr
1666 real(real64),
intent(in) :: time
1667 real(real64),
intent(in) :: dt
1668 real(real64),
intent(in) :: time_delay
1678 call zbatch_init(ff_rs_stateb, 1, 1, hm%dim, gr%np_part)
1679 if (st%pack_states)
call ff_rs_stateb%do_pack(copy=.false.)
1685 hm%cpml_hamiltonian = .false.
1686 call tr%te%apply_batch(namespace, gr, hm, ff_rs_stateb, dt)
1689 call mxll_get_batch(st%rs_state_plane_wavesb, st%rs_state_plane_waves, gr%np, st%dim)
1691 call mxll_set_batch(st%rs_state_plane_wavesb, st%rs_state_plane_waves, gr%np, st%dim)
1692 call ff_rs_stateb%end()
1701 real(real64),
intent(in) :: time
1702 class(
mesh_t),
intent(in) :: mesh
1705 complex(real64),
intent(inout) :: rs_state(:,:)
1707 real(real64) :: e_field_total(mesh%np,st%dim), b_field_total(mesh%np,st%dim)
1708 complex(real64) :: rs_state_add(mesh%np,st%dim)
1718 rs_state_add(1:mesh%np,:), mesh, st%ep, st%mu)
1719 rs_state(1:mesh%np,:) = rs_state(1:mesh%np,:) + rs_state_add(1:mesh%np,:)
1731 type(
grid_t),
intent(in) :: gr
1732 type(namespace_t),
intent(in) :: namespace
1733 real(real64),
intent(in) :: time
1734 real(real64),
intent(in) :: dt
1735 type(
batch_t),
intent(inout) :: rs_stateb
1737 complex(real64),
allocatable :: rs_state(:, :)
1741 safe_allocate(rs_state(gr%np, st%dim))
1743 if (tr%bc_constant)
then
1746 if (st%rs_state_const_external)
then
1767 if (any(hm%bc%bc_ab_type == option__maxwellabsorbingboundaries__mask))
then
1774 if (tr%bc_plane_waves)
then
1781 safe_deallocate_a(rs_state)
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_kernel_start_call(this, file_name, kernel_name, flags)
subroutine, public accel_finish()
pure logical function, public accel_is_enabled()
integer pure function, public accel_max_block_size()
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
subroutine, public zbatch_init(this, dim, st_start, st_end, np, special, packed)
initialize a TYPE_CMPLX valued batch to given size without providing external memory
subroutine, public dbatch_init(this, dim, st_start, st_end, np, special, packed)
initialize a TYPE_FLOAT valued batch to given size without providing external memory
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.
subroutine, public batch_split_complex(np, xx, yy, zz)
extract the real and imaginary parts of a complex batch
subroutine, public batch_set_zero(this, np, async)
fill all mesh functions of the batch with zero
This module implements a calculator for the density and defines related functions.
This module calculates the derivatives (gradients, Laplacians, etc.) of a function.
subroutine, public dderivatives_curl(der, ff, op_ff, ghost_update, set_bc)
apply the curl operator to a vector of mesh functions
subroutine, public zderivatives_batch_grad(der, ffb, opffb, ghost_update, set_bc, to_cartesian, factor)
apply the gradient to a batch of mesh functions
subroutine, public zderivatives_grad(der, ff, op_ff, ghost_update, set_bc, to_cartesian)
apply the gradient to a mesh function
subroutine, public energy_density_calc(mesh, st, rs_field, energy_dens, e_energy_dens, b_energy_dens, plane_waves_check, rs_field_plane_waves, energy_dens_plane_waves)
subroutine, public exponential_init(te, namespace, full_batch)
subroutine, public external_waves_eval(external_waves, time, mesh, type_of_field, out_field_total, der)
Calculation of external waves from parsed formula.
subroutine, public external_waves_init(external_waves, namespace)
Here, plane wave is evaluated from analytical formulae on grid.
Fast Fourier Transform module. This module provides a single interface that works with different FFT ...
real(real64), parameter, public m_two
real(real64), parameter, public m_zero
real(real64), parameter, public m_four
real(real64), parameter, public m_pi
some mathematical constants
real(real64), parameter, public m_fourth
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_epsilon
real(real64), parameter, public m_half
real(real64), parameter, public p_c
Electron gyromagnetic ratio, see Phys. Rev. Lett. 130, 071801 (2023)
real(real64), parameter, public m_one
real(real64), parameter, public m_three
This module implements the underlying real-space grid.
subroutine, public hamiltonian_mxll_apply_simple(hm, namespace, mesh, psib, hpsib, terms, set_bc)
subroutine, public mxll_update_pml_simple(hm, rs_stateb)
integer, parameter, public faraday_ampere_medium
subroutine, public hamiltonian_mxll_update(this, time)
Maxwell Hamiltonian update (here only the time is updated, can maybe be added to another routine)
subroutine, public mxll_copy_pml_simple(hm, rs_stateb)
The Helmholtz decomposition is intended to contain "only mathematical" functions and procedures to co...
This module implements the index, used for the mesh points.
This module is intended to contain "only mathematical" functions and procedures.
pure real(real64) function, dimension(1:3), public dcross_product(a, b)
integer, parameter, public mxll_bc_mirror_pmc
integer, parameter, public mxll_bc_mirror_pec
integer, parameter, public mxll_bc_plane_waves
integer, parameter, public mxll_bc_medium
subroutine, public bc_mxll_generate_pml_parameters(bc, space, gr, c_factor, dt)
integer, parameter, public mxll_bc_constant
integer, parameter, public mxll_bc_zero
This module defines functions over batches 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
subroutine, public dmesh_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_print_with_emphasis(msg, iunit, 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)
this module contains the output system
Some general things and nomenclature:
subroutine, public parse_block_string(blk, l, c, res, convert_to_c)
integer function, public parse_block(namespace, name, blk, check_varinfo_)
subroutine, public dpoisson_solve(this, namespace, pot, rho, all_nodes, kernel, reset)
Calculates the Poisson equation. Given the density returns the corresponding potential.
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.
subroutine, public transform_rs_densities(hm, mesh, rs_current_density, ff_density, sign)
subroutine, public mirror_pec_boundaries_calculation(bc, st, rs_state)
subroutine td_function_mxll_init(st, namespace, hm)
subroutine, public get_vector_pot_and_transverse_field(namespace, gr_mxll, hm_mxll, st_mxll, tr_mxll, hm, poisson_solver, helmholtz, field, transverse_field, vector_potential)
subroutine, public plane_waves_in_box_calculation(bc, time, mesh, der, st, rs_state)
subroutine, public energy_mxll_calc_batch(gr, st, hm, energy_mxll, rs_fieldb, rs_field_plane_wavesb)
subroutine plane_waves_propagation(hm, tr, namespace, st, gr, time, dt, time_delay)
subroutine cpml_conv_function_update(hm, gr, ff_rs_state_pmlb)
subroutine transform_rs_state_batch(hm, gr, st, rs_stateb, ff_rs_stateb, sign)
subroutine, public mxll_propagate_leapfrog(hm, namespace, gr, st, tr, time, dt, counter)
subroutine, public constant_boundaries_calculation(constant_calc, bc, hm, st, rs_state)
subroutine, public mask_absorbing_boundaries(namespace, gr, hm, st, tr, time, dt, time_delay, rs_state)
subroutine, public calculate_vector_potential(namespace, poisson_solver, gr, st, field, vector_potential)
subroutine pml_propagation_stage_2_batch(hm, namespace, gr, st, tr, time, dt, time_delay, ff_rs_state_pmlb, ff_rs_stateb)
subroutine, public mxll_propagation_step(hm, namespace, gr, space, st, tr, rs_stateb, ff_rs_inhom_t1, ff_rs_inhom_t2, time, dt)
subroutine, public plane_waves_boundaries_calculation(hm, st, mesh, time, time_delay, rs_state)
subroutine, public spatial_constant_calculation(constant_calc, st, gr, hm, time, dt, delay, rs_state, set_initial_state)
subroutine, public set_medium_rs_state(st, gr, hm)
subroutine cpml_conv_function_update_via_riemann_silberstein(hm, gr, ff_rs_state_pmlb)
subroutine maxwell_mask(hm, rs_state)
integer, parameter mxwll_etrs_const
subroutine, public calculate_matter_longitudinal_field(gr_mxll, st_mxll, hm_mxll, gr_elec, hm_elec, rs_state_matter)
subroutine, public mxll_propagate_expgauss1(hm, namespace, gr, st, tr, time, dt)
Exponential propagation scheme with Gauss collocation points, s=1.
subroutine transform_rs_densities_to_6x6_rs_densities_forward(mesh, rs_current_density, rs_density_6x6)
subroutine, public energy_mxll_calc(gr, st, hm, energy_mxll, rs_field, rs_field_plane_waves)
subroutine, public mxll_propagate_expgauss2(hm, namespace, gr, st, tr, time, dt)
Exponential propagation scheme with Gauss collocation points, s=2.
integer, parameter, public rs_trans_backward
subroutine transform_rs_densities_to_6x6_rs_densities_backward(mesh, rs_density_6x6, rs_current_density)
subroutine, public mirror_pmc_boundaries_calculation(bc, st, rs_state)
subroutine, public mxll_apply_boundaries(tr, st, hm, gr, namespace, time, dt, rs_stateb)
subroutine pml_propagation_stage_1_batch(hm, gr, st, tr, ff_rs_stateb, ff_rs_state_pmlb)
subroutine, public propagator_mxll_init(gr, namespace, st, hm, tr)
This module defines the quantity_t class and the IDs for quantities, which can be exposed by a system...
subroutine, public mxll_set_batch(rs_stateb, rs_state, np, dim, offset)
subroutine, public get_electric_field_vector(rs_state_vector, electric_field_vector, ep_element)
subroutine, public build_rs_vector(e_vector, b_vector, rs_sign, rs_vector, ep_element, mu_element)
subroutine, public mxll_get_batch(rs_stateb, rs_state, np, dim, offset)
subroutine, public build_rs_state(e_field, b_field, rs_sign, rs_state, mesh, ep_field, mu_field, np)
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 tdf_read(f, namespace, function_name, ierr)
This function initializes "f" from the TDFunctions block.
Class defining batches of mesh functions.
class representing derivatives
Description of the grid, containing information on derivatives, stencil, and symmetries.
Describes mesh distribution to nodes.