34 use,
intrinsic :: iso_fortran_env
74 logical :: adjoint = .false.
76 real(real64) :: current_time
77 logical :: apply_packed
81 type(nl_operator_t),
pointer :: operators(:)
84 type(derivatives_t),
pointer,
private :: der
85 type(states_mxll_t),
pointer :: st
89 logical :: propagation_apply = .false.
91 integer,
pointer :: rs_state_fft_map(:,:,:)
92 integer,
pointer :: rs_state_fft_map_inv(:,:)
94 logical :: mx_ma_coupling = .false.
95 logical :: mx_ma_coupling_apply = .false.
96 integer :: mx_ma_trans_field_calc_method
97 logical :: mx_ma_trans_field_calc_corr = .false.
98 integer :: mx_ma_coupling_points_number
99 real(real64),
allocatable :: mx_ma_coupling_points(:,:)
100 integer,
allocatable :: mx_ma_coupling_points_map(:)
101 integer :: mx_ma_coupling_order
102 logical :: ma_mx_coupling = .false.
103 logical :: ma_mx_coupling_apply = .false.
105 logical :: bc_add_ab_region = .false.
106 logical :: bc_zero = .false.
107 logical :: bc_constant = .false.
108 logical :: bc_mirror_pec = .false.
109 logical :: bc_mirror_pmc = .false.
110 logical :: bc_plane_waves = .false.
111 logical :: bc_medium = .false.
113 logical :: plane_waves = .false.
114 logical :: plane_waves_apply = .false.
115 logical :: spatial_constant = .false.
116 logical :: spatial_constant_apply = .false.
117 logical :: spatial_constant_propagate = .false.
119 logical :: calc_medium_box = .false.
120 type(single_medium_box_t),
allocatable :: medium_boxes(:)
121 logical :: medium_boxes_initialized = .false.
125 logical :: current_density_ext_flag = .false.
126 logical :: current_density_from_medium = .false.
128 type(energy_mxll_t) :: energy
130 logical :: cpml_hamiltonian = .false.
132 logical :: diamag_current = .false.
133 real(real64) :: c_factor
134 real(real64) :: current_factor
137 type(mesh_cube_parallel_map_t) :: mesh_cube_map
146 integer,
public,
parameter :: &
147 FARADAY_AMPERE = 1, &
156 type(hamiltonian_mxll_t),
intent(inout) :: hm
157 type(namespace_t),
intent(in) :: namespace
158 type(grid_t),
target,
intent(inout) :: gr
159 type(states_mxll_t),
target,
intent(inout) :: st
169 assert(
associated(gr%der%lapl))
171 hm%operators(1:3) => gr%der%grad(1:3)
173 hm%rs_sign = st%rs_sign
175 hm%mx_ma_coupling_apply = .false.
176 hm%mx_ma_coupling = .false.
177 hm%ma_mx_coupling_apply = .false.
178 hm%ma_mx_coupling = .false.
195 call parse_variable(namespace,
'MaxwellHamiltonianOperator', faraday_ampere, hm%operator)
208 call parse_variable(namespace,
'ExternalCurrent', .false., hm%current_density_ext_flag)
210 hm%plane_waves_apply = .false.
211 hm%spatial_constant_apply = .false.
212 hm%spatial_constant_propagate = .
true.
214 hm%propagation_apply = .false.
238 hm%rs_state_fft_map => st%rs_state_fft_map
239 hm%rs_state_fft_map_inv => st%rs_state_fft_map_inv
246 call hm%update_span(gr%spacing(1:gr%der%dim),
m_zero, namespace)
263 nullify(hm%operators)
267 if (
allocated(hm%medium_boxes))
then
268 do il = 1,
size(hm%medium_boxes)
285 if (any(hm%bc%bc_ab_type == option__maxwellabsorbingboundaries__cpml))
then
299 real(real64),
intent(in) :: delta(:)
300 real(real64),
intent(in) :: emin
304 real(real64) :: emax, fd_factor
305 real(real64),
parameter :: fd_factors(4) = [1.0_real64, 1.3723_real64, 1.5861_real64, 1.7307_real64]
312 if (hm%der%stencil_type ==
der_star .and. hm%der%order <= 4)
then
313 fd_factor = fd_factors(hm%der%order)
320 do i = 1,
size(delta)
321 emax = emax +
m_one / delta(i)**2
325 hm%spectral_middle_point =
m_zero
326 hm%spectral_half_span = emax
338 if (.not. hm%adjoint)
then
364 real(real64),
optional,
intent(in) :: time
368 this%current_time =
m_zero
369 if (
present(time)) this%current_time = time
379 time = this%current_time
387 class(mesh_t),
intent(in) :: mesh
389 apply = this%apply_packed
390 if (mesh%use_curvilinear) apply = .false.
397 type(namespace_t),
intent(in) :: namespace
398 type(derivatives_t),
intent(in) :: der
399 type(batch_t),
target,
intent(inout) :: psib
400 type(batch_t),
target,
intent(inout) :: hpsib
401 real(real64),
optional,
intent(in) :: time
402 integer,
optional,
intent(in) :: terms
403 logical,
optional,
intent(in) :: set_bc
405 type(batch_t) :: gradb(der%dim)
406 integer :: idir, ifield, field_dir, pml_dir, rs_sign
407 integer :: ip, ip_in, il
408 real(real64) :: pml_a, pml_b, pml_c
409 complex(real64) :: pml_g, grad
410 integer,
parameter :: field_dirs(3, 2) = reshape([2, 3, 1, 3, 1, 2], [3, 2])
411 real(real64) :: cc, aux_ep(3), aux_mu(3), sigma_e, sigma_m, ff_real(3), ff_imag(3), coeff_real, coeff_imag
412 complex(real64) :: ff_plus(3), ff_minus(3), hpsi(6)
413 integer,
parameter :: sign_medium(6) = [1, 1, 1, -1, -1, -1]
414 logical :: with_medium
418 call profiling_in(
"H_MXLL_APPLY_BATCH")
419 assert(psib%status() == hpsib%status())
421 assert(psib%nst == hpsib%nst)
422 assert(hm%st%dim == 3)
423 p_c_ = p_c * hm%c_factor
426 assert(.not.
present(terms))
429 if (
present(time))
then
430 if (abs(time - hm%current_time) > 1e-10_real64)
then
431 write(message(1),
'(a)')
'hamiltonian_apply_batch time assertion failed.'
432 write(message(2),
'(a,f12.6,a,f12.6)')
'time = ', time,
'; hm%current_time = ', hm%current_time
433 call messages_fatal(2, namespace=namespace)
439 call profiling_out(
"H_MXLL_APPLY_BATCH")
444 call zderivatives_batch_grad(der, psib, gradb, set_bc=set_bc)
446 if (hm%cpml_hamiltonian)
then
450 call zderivatives_batch_curl_from_gradient(der, hpsib, gradb)
454 if (hm%bc_constant .and. .not. with_medium)
then
458 if (with_medium)
then
460 if ((hm%bc%bc_type(idir) == mxll_bc_medium))
then
465 if (hm%calc_medium_box)
then
466 do il = 1,
size(hm%medium_boxes)
473 call gradb(idir)%end()
476 call profiling_out(
"H_MXLL_APPLY_BATCH")
481 type(batch_t) :: gradb(:)
482 type(accel_kernel_t),
save :: ker_pml
483 integer :: bsize, gsize
485 call profiling_in(
"APPLY_PML_BOUNDARY")
487 if (with_medium)
then
493 call batch_scal(der%mesh%np, rs_sign * p_c_, gradb(idir))
496 do pml_dir = 1, hm%st%dim
497 if (hm%bc%bc_ab_type(pml_dir) == mxll_ab_cpml)
then
498 select case (gradb(pml_dir)%status())
499 case (batch_not_packed)
501 field_dir = field_dirs(pml_dir, ifield)
503 do ip_in = 1, hm%bc%pml%points_number
504 ip = hm%bc%pml%points_map(ip_in)
505 pml_c = hm%bc%pml%c(ip_in, pml_dir)
506 pml_a = hm%bc%pml%a(ip_in, pml_dir)
507 pml_b = hm%bc%pml%b(ip_in, pml_dir)
508 pml_g = hm%bc%pml%conv_plus(ip_in, pml_dir, field_dir)
509 grad = gradb(pml_dir)%zff_linear(ip, field_dir)
510 gradb(pml_dir)%zff_linear(ip, field_dir) = pml_c * ((m_one + pml_a)*grad/p_c_ &
511 + rs_sign * pml_b*pml_g)
512 if (with_medium)
then
513 pml_g = hm%bc%pml%conv_minus(ip_in, pml_dir, field_dir)
514 grad = gradb(pml_dir)%zff_linear(ip, field_dir+3)
515 gradb(pml_dir)%zff_linear(ip, field_dir+3) = pml_c * ((m_one + pml_a)*grad/p_c_ &
516 + rs_sign * pml_b*pml_g)
522 do ip_in = 1, hm%bc%pml%points_number
523 ip = hm%bc%pml%points_map(ip_in)
524 pml_c = hm%bc%pml%c(ip_in, pml_dir)
525 pml_a = hm%bc%pml%a(ip_in, pml_dir)
526 pml_b = hm%bc%pml%b(ip_in, pml_dir)
528 field_dir = field_dirs(pml_dir, ifield)
529 pml_g = hm%bc%pml%conv_plus(ip_in, pml_dir, field_dir)
530 grad = gradb(pml_dir)%zff_pack(field_dir, ip)
531 gradb(pml_dir)%zff_pack(field_dir, ip) = pml_c * ((m_one + pml_a)*grad/p_c_ &
532 + rs_sign * pml_b*pml_g)
533 if (with_medium)
then
534 pml_g = hm%bc%pml%conv_minus(ip_in, pml_dir, field_dir)
535 grad = gradb(pml_dir)%zff_pack(field_dir+3, ip)
536 gradb(pml_dir)%zff_pack(field_dir+3, ip) = pml_c * ((m_one + pml_a)*grad/p_c_ &
537 + rs_sign * pml_b*pml_g)
541 case (batch_device_packed)
542 call accel_kernel_start_call(ker_pml,
'pml.cu',
'pml_apply')
544 if (with_medium)
then
545 call accel_set_kernel_arg(ker_pml, 0, 1_int32)
547 call accel_set_kernel_arg(ker_pml, 0, 0_int32)
549 call accel_set_kernel_arg(ker_pml, 1, hm%bc%pml%points_number)
550 call accel_set_kernel_arg(ker_pml, 2, pml_dir-1)
551 call accel_set_kernel_arg(ker_pml, 3, p_c_)
552 call accel_set_kernel_arg(ker_pml, 4, rs_sign)
553 call accel_set_kernel_arg(ker_pml, 5, hm%bc%pml%buff_map)
554 call accel_set_kernel_arg(ker_pml, 6, gradb(pml_dir)%ff_device)
555 call accel_set_kernel_arg(ker_pml, 7,
log2(int(gradb(pml_dir)%pack_size(1), int32)))
556 call accel_set_kernel_arg(ker_pml, 8, hm%bc%pml%buff_a)
557 call accel_set_kernel_arg(ker_pml, 9, hm%bc%pml%buff_b)
558 call accel_set_kernel_arg(ker_pml, 10, hm%bc%pml%buff_c)
559 call accel_set_kernel_arg(ker_pml, 11, hm%bc%pml%buff_conv_plus)
560 call accel_set_kernel_arg(ker_pml, 12, hm%bc%pml%buff_conv_minus)
563 bsize = accel_max_block_size()
564 call accel_grid_size(hm%bc%pml%points_number, bsize, gsize)
566 call accel_kernel_run(ker_pml, (/gsize/), (/bsize/))
571 if (accel_is_enabled())
then
575 call profiling_out(
"APPLY_PML_BOUNDARY")
581 call profiling_in(
"SCALE_AFTER_CURL")
582 if (.not. hm%cpml_hamiltonian)
then
584 if (with_medium)
then
586 call batch_scal(der%mesh%np, sign_medium * p_c_, hpsib)
588 call batch_scal(der%mesh%np, hm%rs_sign * p_c_, hpsib)
592 if (with_medium)
then
594 call batch_scal(der%mesh%np, sign_medium * m_one, hpsib)
597 call profiling_out(
"SCALE_AFTER_CURL")
603 call profiling_in(
'APPLY_CONSTANT_BC')
604 select case (hpsib%status())
605 case (batch_not_packed)
606 do field_dir = 1, hm%st%dim
607 do ip_in = 1, hm%bc%constant_points_number
608 ip = hm%bc%constant_points_map(ip_in)
609 hpsib%zff_linear(ip, field_dir) = hm%st%rs_state_const(field_dir)
613 do ip_in = 1, hm%bc%constant_points_number
614 ip = hm%bc%constant_points_map(ip_in)
615 do field_dir = 1, hm%st%dim
616 hpsib%zff_pack(field_dir, ip) = hm%st%rs_state_const(field_dir)
619 case (batch_device_packed)
620 call messages_not_implemented(
"Maxwell constant boundary on GPU", namespace=namespace)
622 call profiling_out(
'APPLY_CONSTANT_BC')
627 type(single_medium_box_t),
intent(in) :: medium
632 call profiling_in(
"MEDIUM_BOX")
633 assert(.not. medium%has_mapping)
636 do ip = 1, medium%points_number
637 cc = medium%c(ip)/p_c
638 if (abs(cc) <= m_epsilon) cycle
639 aux_ep(1:3) = medium%aux_ep(ip, 1:3)
640 aux_mu(1:3) = medium%aux_mu(ip, 1:3)
641 sigma_e = medium%sigma_e(ip)
642 sigma_m = medium%sigma_m(ip)
643 select case (hpsib%status())
644 case (batch_not_packed)
645 ff_plus(1:3) = psib%zff_linear(ip, 1:3)
646 ff_minus(1:3) = psib%zff_linear(ip, 4:6)
647 hpsi(1:6) = hpsib%zff_linear(ip, 1:6)
649 ff_plus(1:3) = psib%zff_pack(1:3, ip)
650 ff_minus(1:3) = psib%zff_pack(4:6, ip)
651 hpsi(1:6) = hpsib%zff_pack(1:6, ip)
652 case (batch_device_packed)
653 call messages_not_implemented(
"Maxwell Medium on GPU", namespace=namespace)
655 ff_real = real(ff_plus+ff_minus, real64)
656 ff_imag = aimag(ff_plus-ff_minus)
657 aux_ep = dcross_product(aux_ep, ff_real)
658 aux_mu = dcross_product(aux_mu, ff_imag)
660 coeff_real = - cc * aux_ep(ifield) + sigma_m * ff_imag(ifield)
661 coeff_imag = - cc * aux_mu(ifield) - sigma_e * ff_real(ifield)
662 hpsi(ifield) = cc * hpsi(ifield) + cmplx(coeff_real, coeff_imag, real64)
663 hpsi(ifield+3) = cc * hpsi(ifield+3) + cmplx(-coeff_real, coeff_imag, real64)
665 select case (hpsib%status())
666 case (batch_not_packed)
667 hpsib%zff_linear(ip, 1:6) = hpsi(1:6)
669 hpsib%zff_pack(1:6, ip) = hpsi(1:6)
672 call profiling_out(
"MEDIUM_BOX")
681 type(namespace_t),
intent(in) :: namespace
682 class(mesh_t),
intent(in) :: mesh
683 type(batch_t),
target,
intent(inout) :: psib
684 type(batch_t),
target,
intent(inout) :: hpsib
685 integer,
optional,
intent(in) :: terms
686 logical,
optional,
intent(in) :: set_bc
688 type(batch_t) :: gradb(hm%der%dim)
692 call profiling_in(
"MXLL_HAMILTONIAN_SIMPLE")
694 call zderivatives_batch_grad(hm%der, psib, gradb, set_bc=set_bc)
696 if (any(hm%bc%bc_ab_type == option__maxwellabsorbingboundaries__cpml))
then
702 call zderivatives_batch_curl_from_gradient(hm%der, hpsib, gradb)
704 if (hm%calc_medium_box)
then
709 call batch_scal(mesh%np, p_c*hm%c_factor, hpsib)
712 do idir = 1, hm%der%dim
713 call gradb(idir)%end()
716 call profiling_out(
"MXLL_HAMILTONIAN_SIMPLE")
723 type(batch_t),
intent(inout) :: gradb(1:hm%st%dim)
725 integer :: idir, jdir, ip_in, ip, bsize, gsize
726 type(accel_kernel_t),
save :: ker_pml
730 call profiling_in(
"APPLY_PML_SIMPLE")
733 do jdir = 1, hm%st%dim
734 select case (gradb(1)%status())
735 case (batch_not_packed)
737 do idir = 1, hm%st%dim
738 if (idir == jdir) cycle
739 do ip_in = 1, hm%bc%pml%points_number
740 ip = hm%bc%pml%points_map(ip_in)
742 gradb(jdir)%zff_linear(ip, idir) = gradb(jdir)%zff_linear(ip, idir) * (m_one + hm%bc%pml%c(ip_in, jdir)) + &
743 hm%bc%pml%b(ip_in, jdir) * hm%bc%pml%conv_plus_old(ip_in, jdir, idir)
747 do ip_in = 1, hm%bc%pml%points_number
748 ip = hm%bc%pml%points_map(ip_in)
750 do idir = 1, hm%st%dim
751 if (idir == jdir) cycle
753 gradb(jdir)%zff_pack(idir, ip) = gradb(jdir)%zff_pack(idir, ip) * (m_one + hm%bc%pml%c(ip_in, jdir)) +&
754 hm%bc%pml%b(ip_in, jdir) * hm%bc%pml%conv_plus_old(ip_in, jdir, idir)
757 case (batch_device_packed)
758 call accel_kernel_start_call(ker_pml,
'pml.cu',
'pml_apply_new')
760 call accel_set_kernel_arg(ker_pml, 0, hm%bc%pml%points_number)
761 call accel_set_kernel_arg(ker_pml, 1, jdir-1)
762 call accel_set_kernel_arg(ker_pml, 2, hm%bc%pml%buff_map)
763 call accel_set_kernel_arg(ker_pml, 3, gradb(jdir)%ff_device)
764 call accel_set_kernel_arg(ker_pml, 4,
log2(int(gradb(jdir)%pack_size(1), int32)))
765 call accel_set_kernel_arg(ker_pml, 5, hm%bc%pml%buff_b)
766 call accel_set_kernel_arg(ker_pml, 6, hm%bc%pml%buff_c)
767 call accel_set_kernel_arg(ker_pml, 7, hm%bc%pml%buff_conv_plus_old)
770 bsize = accel_max_block_size()
771 call accel_grid_size(hm%bc%pml%points_number, bsize, gsize)
773 call accel_kernel_run(ker_pml, (/gsize/), (/bsize/))
776 call profiling_out(
"APPLY_PML_SIMPLE")
783 type(batch_t),
intent(inout) :: rs_stateb
785 integer :: idir, jdir, ip, ip_in, bsize, gsize
786 type(accel_kernel_t),
save :: ker_pml_update
787 type(batch_t) :: gradb(1:hm%st%dim)
790 call profiling_in(
"UPDATE_PML_SIMPLE")
792 call zderivatives_batch_grad(hm%der, rs_stateb, gradb)
794 do jdir = 1, hm%st%dim
795 call rs_stateb%check_compatibility_with(gradb(jdir))
796 select case (gradb(jdir)%status())
797 case (batch_not_packed)
799 do idir = 1, hm%st%dim
800 if (idir == jdir) cycle
801 do ip_in = 1, hm%bc%pml%points_number
802 ip = hm%bc%pml%points_map(ip_in)
803 hm%bc%pml%conv_plus(ip_in, jdir, idir) = hm%bc%pml%c(ip_in, jdir) * gradb(jdir)%zff_linear(ip, idir) +&
804 hm%bc%pml%b(ip_in, jdir) * hm%bc%pml%conv_plus_old(ip_in, jdir, idir)
808 do ip_in = 1, hm%bc%pml%points_number
809 ip = hm%bc%pml%points_map(ip_in)
811 do idir = 1, hm%st%dim
812 if (idir == jdir) cycle
813 hm%bc%pml%conv_plus(ip_in, jdir, idir) = hm%bc%pml%c(ip_in, jdir) * gradb(jdir)%zff_pack(idir, ip) +&
814 hm%bc%pml%b(ip_in, jdir) * hm%bc%pml%conv_plus_old(ip_in, jdir, idir)
817 case (batch_device_packed)
818 call accel_kernel_start_call(ker_pml_update,
'pml.cu',
'pml_update_new')
820 call accel_set_kernel_arg(ker_pml_update, 0, hm%bc%pml%points_number)
821 call accel_set_kernel_arg(ker_pml_update, 1, jdir-1)
822 call accel_set_kernel_arg(ker_pml_update, 2, hm%bc%pml%buff_map)
823 call accel_set_kernel_arg(ker_pml_update, 3, gradb(jdir)%ff_device)
824 call accel_set_kernel_arg(ker_pml_update, 4,
log2(int(gradb(jdir)%pack_size(1), int32)))
825 call accel_set_kernel_arg(ker_pml_update, 5, hm%bc%pml%buff_b)
826 call accel_set_kernel_arg(ker_pml_update, 6, hm%bc%pml%buff_c)
827 call accel_set_kernel_arg(ker_pml_update, 7, hm%bc%pml%buff_conv_plus)
828 call accel_set_kernel_arg(ker_pml_update, 8, hm%bc%pml%buff_conv_plus_old)
831 bsize = accel_max_block_size()
832 call accel_grid_size(hm%bc%pml%points_number, bsize, gsize)
834 call accel_kernel_run(ker_pml_update, (/gsize/), (/bsize/))
836 call gradb(jdir)%end()
838 call profiling_out(
"UPDATE_PML_SIMPLE")
845 type(batch_t),
intent(inout) :: rs_stateb
847 integer :: bsize, idir, jdir, ip, ip_in
849 type(accel_kernel_t),
save :: ker_pml_copy
852 call profiling_in(
"COPY_PML_SIMPLE")
854 select case (rs_stateb%status())
855 case (batch_packed, batch_not_packed)
856 do jdir = 1, hm%st%dim
858 do idir = 1, hm%st%dim
859 if (idir == jdir) cycle
860 do ip_in = 1, hm%bc%pml%points_number
861 ip = hm%bc%pml%points_map(ip_in)
862 hm%bc%pml%conv_plus_old(ip_in, jdir, idir) = hm%bc%pml%conv_plus(ip_in, jdir, idir)
866 case (batch_device_packed)
867 call accel_kernel_start_call(ker_pml_copy,
'pml.cu',
'pml_copy')
869 call accel_set_kernel_arg(ker_pml_copy, 0, hm%bc%pml%points_number*9)
870 call accel_set_kernel_arg(ker_pml_copy, 1, hm%bc%pml%buff_conv_plus)
871 call accel_set_kernel_arg(ker_pml_copy, 2, hm%bc%pml%buff_conv_plus_old)
874 bsize = accel_max_block_size()
875 call accel_grid_size(hm%bc%pml%points_number*9, bsize, gsize)
877 call accel_kernel_run(ker_pml_copy, (/gsize/), (/bsize/))
879 call profiling_out(
"COPY_PML_SIMPLE")
886 type(batch_t),
intent(inout) :: rs_stateb
889 logical :: need_to_pack
892 call profiling_in(
"LINEAR_MEDIUM_SIMPLE")
894 do il = 1,
size(hm%medium_boxes)
895 assert(.not. hm%medium_boxes(il)%has_mapping)
896 need_to_pack = .false.
898 if(rs_stateb%status() == batch_device_packed)
then
899 call rs_stateb%do_unpack(force=.
true.)
900 need_to_pack = .
true.
902 select case (rs_stateb%status())
903 case (batch_not_packed)
904 do ip = 1, hm%medium_boxes(il)%points_number
905 if (abs(hm%medium_boxes(il)%c(ip)) <= m_epsilon)
then
907 rs_stateb%zff_linear(ip, 1:3) = p_c*hm%c_factor*rs_stateb%zff_linear(ip, 1:3)
910 rs_stateb%zff_linear(ip, 1:3) = hm%medium_boxes(il)%c(ip)*(rs_stateb%zff_linear(ip, 1:3) + &
912 -dcross_product(hm%medium_boxes(il)%aux_ep(ip, 1:3)*m_two, &
913 real(rs_stateb%zff_linear(ip, 1:3), real64)), &
914 -dcross_product(hm%medium_boxes(il)%aux_mu(ip, 1:3)*m_two, &
915 aimag(rs_stateb%zff_linear(ip, 1:3))), real64)) + &
917 hm%medium_boxes(il)%sigma_m(ip)*aimag(rs_stateb%zff_linear(ip, 1:3)), &
918 -hm%medium_boxes(il)%sigma_e(ip)*
real(rs_stateb%zff_linear(ip, 1:3), real64), real64)
922 do ip = 1, hm%medium_boxes(il)%points_number
923 if (abs(hm%medium_boxes(il)%c(ip)) <= m_epsilon)
then
925 rs_stateb%zff_pack(1:3, ip) = p_c*hm%c_factor*rs_stateb%zff_pack(1:3, ip)
928 rs_stateb%zff_pack(1:3, ip) = hm%medium_boxes(il)%c(ip)*(rs_stateb%zff_pack(1:3, ip) + &
930 -dcross_product(hm%medium_boxes(il)%aux_ep(ip, 1:3)*m_two, &
931 real(rs_stateb%zff_pack(1:3, ip), real64)), &
932 -dcross_product(hm%medium_boxes(il)%aux_mu(ip, 1:3)*m_two, &
933 aimag(rs_stateb%zff_pack(1:3, ip))), real64)) + &
935 hm%medium_boxes(il)%sigma_m(ip)*aimag(rs_stateb%zff_pack(1:3, ip)), &
936 -hm%medium_boxes(il)%sigma_e(ip)*
real(rs_stateb%zff_pack(1:3, ip), real64), real64)
940 if(need_to_pack)
then
941 call rs_stateb%do_pack()
945 call profiling_out(
"LINEAR_MEDIUM_SIMPLE")
953 type(namespace_t),
intent(in) :: namespace
954 class(mesh_t),
intent(in) :: mesh
955 class(batch_t),
target,
intent(inout) :: psib
956 class(batch_t),
target,
intent(inout) :: hpsib
957 integer,
optional,
intent(in) :: terms
958 logical,
optional,
intent(in) :: set_bc
960 write(message(1),
'(a)')
'dhamiltonian_mxll_apply not implemented (states are complex).'
961 call messages_fatal(1, namespace=namespace)
969 type(namespace_t),
intent(in) :: namespace
970 class(mesh_t),
intent(in) :: mesh
971 class(batch_t),
target,
intent(inout) :: psib
972 class(batch_t),
target,
intent(inout) :: hpsib
973 integer,
optional,
intent(in) :: terms
974 logical,
optional,
intent(in) :: set_bc
976 complex(real64),
allocatable :: rs_aux_in(:,:), rs_aux_out(:,:)
982 call profiling_in(
'ZHAMILTONIAN_MXLL_APPLY')
984 on_gpu = psib%status() == batch_device_packed
987 safe_allocate(rs_aux_in(1:mesh%np_part, 1:hm%dim))
988 safe_allocate(rs_aux_out(1:mesh%np, 1:hm%dim))
989 call boundaries_set(hm%der%boundaries, mesh, psib)
991 call batch_get_state(psib, ii,mesh%np_part, rs_aux_in(:, ii))
996 call batch_set_state(hpsib, ii, mesh%np, rs_aux_out(:, ii))
998 safe_deallocate_a(rs_aux_in)
999 safe_deallocate_a(rs_aux_out)
1005 call profiling_out(
'ZHAMILTONIAN_MXLL_APPLY')
1015 type(derivatives_t),
intent(in) :: der
1016 complex(real64),
contiguous,
intent(inout) :: psi(:,:)
1017 complex(real64),
intent(inout) :: oppsi(:,:)
1019 real(real64),
pointer :: mx_rho(:,:)
1020 complex(real64),
allocatable :: tmp(:,:)
1021 complex(real64),
pointer :: kappa_psi(:,:)
1022 integer :: np, np_part, ip, ip_in, rs_sign
1023 real(real64) :: P_c_
1027 call profiling_in(
'MAXWELL_HAMILTONIAN_APPLY_FD')
1030 np_part = der%mesh%np_part
1031 rs_sign = hm%rs_sign
1032 p_c_ = p_c * hm%c_factor
1034 select case (hm%operator)
1040 call profiling_in(
'MXLL_HAM_APPLY_FD_FARADAY_AMP')
1042 safe_allocate(tmp(1:np,1:2))
1045 if (hm%diamag_current)
then
1046 mx_rho => hm%st%grid_rho
1047 kappa_psi => hm%st%kappa_psi
1053 call zderivatives_partial(der, psi(:,3), tmp(:,1), 2, set_bc = .false.)
1054 call zderivatives_partial(der, psi(:,2), tmp(:,2), 3, set_bc = .false.)
1055 tmp = rs_sign * p_c_ * tmp
1058 oppsi(1:np,1) = ( tmp(1:np,1)-tmp(1:np,2))
1063 call zderivatives_partial(der, psi(:,1), tmp(:,1), 3, set_bc = .false.)
1064 call zderivatives_partial(der, psi(:,3), tmp(:,2), 1, set_bc = .false.)
1065 tmp = rs_sign * p_c_ * tmp
1068 oppsi(1:np,2) = ( tmp(1:np,1)-tmp(1:np,2))
1073 call zderivatives_partial(der, psi(:,2), tmp(:,1), 1, set_bc = .false.)
1074 call zderivatives_partial(der, psi(:,1), tmp(:,2), 2, set_bc = .false.)
1075 tmp = rs_sign * p_c_ * tmp
1078 oppsi(1:np,3) = ( tmp(1:np,1)-tmp(1:np,2))
1080 if (hm%bc_constant)
then
1081 do ip_in = 1, hm%bc%constant_points_number
1082 ip = hm%bc%constant_points_map(ip_in)
1083 oppsi(ip,:) = hm%st%rs_state_const(:)
1087 safe_deallocate_a(tmp)
1089 call profiling_out(
'MXLL_HAM_APPLY_FD_FARADAY_AMP')
1094 call profiling_in(
'MXLL_HAM_APP_FAR_AMP_MED')
1096 safe_allocate(tmp(1:np,1:4))
1102 call zderivatives_partial(der, psi(:,3), tmp(:,1), 2, set_bc = .false.)
1103 call zderivatives_partial(der, psi(:,2), tmp(:,3), 3, set_bc = .false.)
1104 call zderivatives_partial(der, psi(:,6), tmp(:,2), 2, set_bc = .false.)
1105 call zderivatives_partial(der, psi(:,5), tmp(:,4), 3, set_bc = .false.)
1109 oppsi(1:np,1) = (tmp(1:np,1)-tmp(1:np,3))
1110 oppsi(1:np,4) = -(tmp(1:np,2)-tmp(1:np,4))
1115 call zderivatives_partial(der, psi(:,1), tmp(:,1), 3, set_bc = .false.)
1116 call zderivatives_partial(der, psi(:,3), tmp(:,3), 1, set_bc = .false.)
1117 call zderivatives_partial(der, psi(:,4), tmp(:,2), 3, set_bc = .false.)
1118 call zderivatives_partial(der, psi(:,6), tmp(:,4), 1, set_bc = .false.)
1122 oppsi(1:np,2) = (tmp(1:np,1)-tmp(1:np,3))
1123 oppsi(1:np,5) = -(tmp(1:np,2)-tmp(1:np,4))
1128 call zderivatives_partial(der, psi(:,2), tmp(:,1), 1, set_bc = .false.)
1129 call zderivatives_partial(der, psi(:,1), tmp(:,3), 2, set_bc = .false.)
1130 call zderivatives_partial(der, psi(:,5), tmp(:,2), 1, set_bc = .false.)
1131 call zderivatives_partial(der, psi(:,4), tmp(:,4), 2, set_bc = .false.)
1135 oppsi(1:np,3) = (tmp(1:np,1)-tmp(1:np,3))
1136 oppsi(1:np,6) = -(tmp(1:np,2)-tmp(1:np,4))
1139 safe_deallocate_a(tmp)
1149 call profiling_out(
'MXLL_HAM_APP_FAR_AMP_MED')
1153 call profiling_out(
'MAXWELL_HAMILTONIAN_APPLY_FD')
1163 type(derivatives_t),
intent(in) :: der
1164 complex(real64),
intent(inout) :: psi(:,:)
1165 integer,
intent(in) :: dir1
1166 integer,
intent(in) :: dir2
1167 complex(real64),
intent(inout) :: tmp(:)
1172 call profiling_in(
'MAXWELL_PML_HAMILTONIAN')
1174 if ((hm%bc%bc_ab_type(dir1) == mxll_ab_cpml) .and. hm%cpml_hamiltonian)
then
1178 call profiling_out(
'MAXWELL_PML_HAMILTONIAN')
1187 type(derivatives_t),
intent(in) :: der
1188 complex(real64),
intent(inout) :: psi(:,:)
1189 integer,
intent(in) :: dir1
1190 integer,
intent(in) :: dir2
1191 complex(real64),
intent(inout) :: tmp(:,:)
1196 call profiling_in(
'MAXWELL_PML_HAMILTONIAN_MEDIUM')
1198 if ((hm%bc%bc_ab_type(dir1) == mxll_ab_cpml) .and. hm%cpml_hamiltonian)
then
1202 call profiling_out(
'MAXWELL_PML_HAMILTONIAN_MEDIUM')
1211 type(derivatives_t),
intent(in) :: der
1212 integer,
intent(in) :: pml_dir
1213 complex(real64),
intent(inout) :: psi(:,:)
1214 integer,
intent(in) :: field_dir
1215 complex(real64),
intent(inout) :: pml(:)
1217 integer :: ip, ip_in, rs_sign
1218 real(real64) :: pml_c
1219 complex(real64),
allocatable :: tmp_partial(:)
1220 complex(real64) :: pml_a, pml_b, pml_g
1224 if (hm%cpml_hamiltonian)
then
1226 rs_sign = hm%rs_sign
1228 safe_allocate(tmp_partial(1:der%mesh%np_part))
1230 call zderivatives_partial(der, psi(:,field_dir), tmp_partial(:), pml_dir, set_bc = .false.)
1231 do ip_in = 1, hm%bc%pml%points_number
1232 ip = hm%bc%pml%points_map(ip_in)
1233 pml_c = hm%bc%pml%c(ip_in, pml_dir)
1234 pml_a = hm%bc%pml%a(ip_in, pml_dir)
1235 pml_b = hm%bc%pml%b(ip_in, pml_dir)
1236 pml_g = hm%bc%pml%conv_plus(ip_in, pml_dir, field_dir)
1237 pml(ip) = rs_sign * pml_c * ( tmp_partial(ip) &
1238 + real(pml_a, real64) * real(tmp_partial(ip), real64) &
1239 + m_zi * aimag(pml_a) * aimag(tmp_partial(ip)) &
1240 + real(pml_b, real64) * real(pml_g, real64) &
1241 + m_zi * aimag(pml_b) * aimag(pml_g))
1244 safe_deallocate_a(tmp_partial)
1256 type(derivatives_t),
intent(in) :: der
1257 integer,
intent(in) :: pml_dir
1258 complex(real64),
intent(inout) :: psi(:,:)
1259 integer,
intent(in) :: field_dir
1260 complex(real64),
intent(inout) :: pml(:,:)
1262 integer :: ip, ip_in, np
1263 real(real64) :: pml_c(3)
1264 complex(real64),
allocatable :: tmp_partial(:,:)
1265 complex(real64) :: pml_a(3), pml_b(3), pml_g_p(3), pml_g_m(3)
1269 if (hm%cpml_hamiltonian)
then
1272 safe_allocate(tmp_partial(1:np,1:2))
1274 call zderivatives_partial(der, psi(:,field_dir ), tmp_partial(:,1), pml_dir, set_bc = .false.)
1275 call zderivatives_partial(der, psi(:,field_dir+3), tmp_partial(:,2), pml_dir, set_bc = .false.)
1276 do ip_in = 1, hm%bc%pml%points_number
1277 ip = hm%bc%pml%points_map(ip_in)
1278 pml_c(1:3) = hm%bc%pml%c(ip_in, 1:3)
1279 pml_a(1:3) = hm%bc%pml%a(ip_in, 1:3)
1280 pml_b(1:3) = hm%bc%pml%b(ip_in, 1:3)
1281 pml_g_p(1:3) = hm%bc%pml%conv_plus(ip_in, pml_dir, 1:3)
1282 pml_g_m(1:3) = hm%bc%pml%conv_minus(ip_in, pml_dir, 1:3)
1283 pml(ip, 1) = pml_c(pml_dir) * tmp_partial(ip, 1) &
1284 + pml_c(pml_dir) * real(pml_a(pml_dir), real64) * real(tmp_partial(ip, 1), real64) &
1285 + m_zi * pml_c(pml_dir) * aimag(pml_a(pml_dir)) * aimag(tmp_partial(ip, 1)) &
1286 + pml_c(pml_dir) * real(pml_b(pml_dir), real64) * real(pml_g_p(field_dir), real64) &
1287 + m_zi * pml_c(pml_dir) * aimag(pml_b(pml_dir)) * aimag(pml_g_p(field_dir))
1288 pml(ip, 2) = pml_c(pml_dir) * tmp_partial(ip, 2) &
1289 + pml_c(pml_dir) * real(pml_a(pml_dir), real64) * real(tmp_partial(ip, 2), real64) &
1290 + m_zi * pml_c(pml_dir) * aimag(pml_a(pml_dir)) * aimag(tmp_partial(ip, 2)) &
1291 + pml_c(pml_dir) * real(pml_b(pml_dir), real64) * real(pml_g_m(field_dir), real64) &
1292 + m_zi * pml_c(pml_dir) * aimag(pml_b(pml_dir)) * aimag(pml_g_m(field_dir))
1297 safe_deallocate_a(tmp_partial)
1307 complex(real64),
intent(in) :: psi(:,:)
1308 complex(real64),
intent(inout) :: oppsi(:,:)
1310 integer :: ip, ip_in, idim
1311 real(real64) :: cc, aux_ep(3), aux_mu(3), sigma_e, sigma_m
1312 complex(real64) :: ff_plus(3), ff_minus(3)
1317 if ((hm%bc%bc_type(idim) == mxll_bc_medium))
then
1318 do ip_in = 1, hm%bc%medium(idim)%points_number
1319 ip = hm%bc%medium(idim)%points_map(ip_in)
1320 cc = hm%bc%medium(idim)%c(ip_in)/p_c
1321 aux_ep(:) = hm%bc%medium(idim)%aux_ep(ip_in, :)
1322 aux_mu(:) = hm%bc%medium(idim)%aux_mu(ip_in, :)
1323 sigma_e = hm%bc%medium(idim)%sigma_e(ip_in)
1324 sigma_m = hm%bc%medium(idim)%sigma_m(ip_in)
1325 ff_plus(1) = psi(ip, 1)
1326 ff_plus(2) = psi(ip, 2)
1327 ff_plus(3) = psi(ip, 3)
1328 ff_minus(1) = psi(ip, 4)
1329 ff_minus(2) = psi(ip, 5)
1330 ff_minus(3) = psi(ip, 6)
1331 aux_ep = dcross_product(aux_ep,real(ff_plus+ff_minus, real64) )
1332 aux_mu = dcross_product(aux_mu,aimag(ff_plus-ff_minus))
1333 oppsi(ip, 1) = oppsi(ip, 1)*cc &
1334 - cc * aux_ep(1) - cc * m_zi * aux_mu(1) &
1335 - m_zi * sigma_e * real(ff_plus(1) + ff_minus(1), real64) &
1336 - m_zi * sigma_m * m_zi * aimag(ff_plus(1) - ff_minus(1))
1337 oppsi(ip, 4) = oppsi(ip, 4)*cc &
1338 + cc * aux_ep(1) - cc * m_zi * aux_mu(1) &
1339 - m_zi * sigma_e * real(ff_plus(1) + ff_minus(1), real64) &
1340 + m_zi * sigma_m * m_zi * aimag(ff_plus(1) - ff_minus(1))
1341 oppsi(ip, 2) = oppsi(ip, 2)*cc &
1342 - cc * aux_ep(2) - cc * m_zi * aux_mu(2) &
1343 - m_zi * sigma_e * real(ff_plus(2) + ff_minus(2), real64) &
1344 - m_zi * sigma_m * m_zi * aimag(ff_plus(2) - ff_minus(2))
1345 oppsi(ip, 5) = oppsi(ip, 5)*cc &
1346 + cc * aux_ep(2) - cc * m_zi * aux_mu(2) &
1347 - m_zi * sigma_e * real(ff_plus(2) + ff_minus(2), real64) &
1348 + m_zi * sigma_m * m_zi * aimag(ff_plus(2) - ff_minus(2))
1349 oppsi(ip, 3) = oppsi(ip, 3)*cc &
1350 - cc * aux_ep(3) - cc * m_zi * aux_mu(3) &
1351 - m_zi * sigma_e * real(ff_plus(3) + ff_minus(3), real64) &
1352 - m_zi * sigma_m * m_zi * aimag(ff_plus(3) - ff_minus(3))
1353 oppsi(ip, 6) = oppsi(ip, 6)*cc &
1354 + cc * aux_ep(3) - cc * m_zi * aux_mu(3) &
1355 - m_zi * sigma_e * real(ff_plus(3) + ff_minus(3), real64) &
1356 + m_zi * sigma_m * m_zi * aimag(ff_plus(3) - ff_minus(3))
1369 type(derivatives_t),
intent(in) :: der
1370 complex(real64),
intent(in) :: psi(:,:)
1371 complex(real64),
intent(inout) :: oppsi(:,:)
1374 real(real64) :: cc, aux_ep(3), aux_mu(3), sigma_e, sigma_m
1375 complex(real64) :: ff_plus(3), ff_minus(3)
1379 if (hm%calc_medium_box)
then
1380 do il = 1,
size(hm%medium_boxes)
1381 assert(.not. hm%medium_boxes(il)%has_mapping)
1382 do ip = 1, hm%medium_boxes(il)%points_number
1383 cc = hm%medium_boxes(il)%c(ip)/p_c
1384 if (abs(cc) <= m_epsilon) cycle
1385 aux_ep(1:3) = hm%medium_boxes(il)%aux_ep(ip, 1:3)
1386 aux_mu(1:3) = hm%medium_boxes(il)%aux_mu(ip, 1:3)
1387 sigma_e = hm%medium_boxes(il)%sigma_e(ip)
1388 sigma_m = hm%medium_boxes(il)%sigma_m(ip)
1389 ff_plus(1) = psi(ip, 1)
1390 ff_plus(2) = psi(ip, 2)
1391 ff_plus(3) = psi(ip, 3)
1392 ff_minus(1) = psi(ip, 4)
1393 ff_minus(2) = psi(ip, 5)
1394 ff_minus(3) = psi(ip, 6)
1395 aux_ep = dcross_product(aux_ep, real(ff_plus+ff_minus, real64) )
1396 aux_mu = dcross_product(aux_mu, aimag(ff_plus-ff_minus))
1397 oppsi(ip, 1) = oppsi(ip,1)*cc &
1398 - cc * aux_ep(1) - cc * m_zi * aux_mu(1) &
1399 - m_zi * sigma_e * real(ff_plus(1) + ff_minus(1), real64) &
1400 - m_zi * sigma_m * m_zi * aimag(ff_plus(1) - ff_minus(1))
1401 oppsi(ip, 4) = oppsi(ip,4)*cc &
1402 + cc * aux_ep(1) - cc * m_zi * aux_mu(1) &
1403 - m_zi * sigma_e * real(ff_plus(1) + ff_minus(1), real64) &
1404 + m_zi * sigma_m * m_zi * aimag(ff_plus(1) - ff_minus(1))
1405 oppsi(ip, 2) = oppsi(ip,2)*cc &
1406 - cc * aux_ep(2) - cc * m_zi * aux_mu(2) &
1407 - m_zi * sigma_e * real(ff_plus(2) + ff_minus(2), real64) &
1408 - m_zi * sigma_m * m_zi * aimag(ff_plus(2) - ff_minus(2))
1409 oppsi(ip, 5) = oppsi(ip,5)*cc &
1410 + cc * aux_ep(2) - cc * m_zi * aux_mu(2) &
1411 - m_zi * sigma_e * real(ff_plus(2) + ff_minus(2), real64) &
1412 + m_zi * sigma_m * m_zi * aimag(ff_plus(2) - ff_minus(2))
1413 oppsi(ip, 3) = oppsi(ip,3)*cc &
1414 - cc * aux_ep(3) - cc * m_zi * aux_mu(3) &
1415 - m_zi * sigma_e * real(ff_plus(3) + ff_minus(3), real64) &
1416 - m_zi * sigma_m * m_zi * aimag(ff_plus(3) - ff_minus(3))
1417 oppsi(ip, 6) = oppsi(ip,6)*cc &
1418 + cc * aux_ep(3) - cc * m_zi * aux_mu(3) &
1419 - m_zi * sigma_e * real(ff_plus(3) + ff_minus(3), real64) &
1420 + m_zi * sigma_m * m_zi * aimag(ff_plus(3) - ff_minus(3))
subroutine apply_medium_box(medium)
subroutine apply_constant_boundary
subroutine apply_pml_boundary(gradb)
subroutine scale_after_curl
double log2(double __x) __attribute__((__nothrow__
double sqrt(double __x) __attribute__((__nothrow__
This module implements batches of mesh functions.
This module implements common operations on batches of mesh functions.
Module implementing boundary conditions in Octopus.
This module calculates the derivatives (gradients, Laplacians, etc.) of a function.
integer, parameter, public der_star
real(real64), parameter, public m_zero
real(real64), parameter, public m_pi
some mathematical constants
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.
This module defines an abstract class for Hamiltonians.
subroutine mxll_apply_pml_simple(hm, gradb)
subroutine, public hamiltonian_mxll_init(hm, namespace, gr, st)
Initializing the Maxwell Hamiltonian.
subroutine, public zhamiltonian_mxll_apply(hm, namespace, mesh, psib, hpsib, terms, set_bc)
Applying the Maxwell Hamiltonian on Maxwell states.
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
real(real64) function, public hamiltonian_mxll_get_time(this)
subroutine, public hamiltonian_mxll_end(hm)
subroutine maxwell_hamiltonian_apply_fd(hm, der, psi, oppsi)
Applying the Maxwell Hamiltonian on Maxwell states with finite difference.
subroutine maxwell_medium_boundaries_calculation(hm, psi, oppsi)
Maxwell Hamiltonian for medium boundaries.
subroutine maxwell_pml_calculation_via_riemann_silberstein_medium(hm, der, psi, pml_dir, field_dir, pml)
Maxwell Hamiltonian is updated for the PML calculation via Riemann-Silberstein vector with medium ins...
subroutine, public hamiltonian_mxll_span(hm, delta, emin, namespace)
integer, parameter, public faraday_ampere_medium
integer, parameter, public mxll_simple
subroutine, public dhamiltonian_mxll_apply(hm, namespace, mesh, psib, hpsib, terms, set_bc)
Apply hamiltonian to real states (not possible)
logical function, public hamiltonian_mxll_hermitian(hm)
subroutine maxwell_pml_calculation_via_riemann_silberstein(hm, der, psi, pml_dir, field_dir, pml)
Maxwell Hamiltonian is updated for the PML calculation via Riemann-Silberstein vector.
subroutine, public hamiltonian_mxll_update(this, time)
Maxwell Hamiltonian update (here only the time is updated, can maybe be added to another routine)
subroutine maxwell_pml_hamiltonian(hm, der, psi, dir1, dir2, tmp)
Maxwell Hamiltonian is updated for the PML calculation.
subroutine, public mxll_copy_pml_simple(hm, rs_stateb)
subroutine, public hamiltonian_mxll_apply_batch(hm, namespace, der, psib, hpsib, time, terms, set_bc)
subroutine, public hamiltonian_mxll_not_adjoint(hm)
subroutine, public hamiltonian_mxll_adjoint(hm)
subroutine maxwell_pml_hamiltonian_medium(hm, der, psi, dir1, dir2, tmp)
Maxwell Hamiltonian is updated for the PML calculation.
logical pure function, public hamiltonian_mxll_apply_packed(this, mesh)
subroutine mxll_linear_medium_terms_simple(hm, rs_stateb)
subroutine maxwell_medium_boxes_calculation(hm, der, psi, oppsi)
subroutine, public single_medium_box_end(medium_box)
Deallocation of medium_box components.
This module is intended to contain "only mathematical" functions and procedures.
subroutine, public bc_mxll_end(bc)
This module defines the meshes, which are used in Octopus.
subroutine, public messages_experimental(name, namespace)
This module defines non-local operators.
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.
The abstract Hamiltonian class defines a skeleton for specific implementations.