Octopus
hamiltonian_mxll.F90
Go to the documentation of this file.
1!! Copyright (C) 2019 R. Jestaedt, H. Appel, F. Bonafe, M. Oliveira, N. Tancogne-Dejean
2!!
3!! This program is free software; you can redistribute it and/or modify
4!! it under the terms of the GNU General Public License as published by
5!! the Free Software Foundation; either version 2, or (at your option)
6!! any later version.
7!!
8!! This program is distributed in the hope that it will be useful,
9!! but WITHOUT ANY WARRANTY; without even the implied warranty of
10!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
11!! GNU General Public License for more details.
12!!
13!! You should have received a copy of the GNU General Public License
14!! along with this program; if not, write to the Free Software
15!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
16!! 02110-1301, USA.
17!!
18
19#include "global.h"
20
22 use accel_oct_m
23 use batch_oct_m
26 use cube_oct_m
27 use debug_oct_m
30 use global_oct_m
31 use grid_oct_m
34 use, intrinsic :: iso_fortran_env
36 use math_oct_m
39 use mesh_oct_m
43 use parser_oct_m
49
50 implicit none
51
52 private
53 public :: &
70
72 integer :: dim
74 logical :: adjoint = .false.
75
76 real(real64) :: current_time
77 logical :: apply_packed
78
79 logical :: time_zero
80
81 type(nl_operator_t), pointer :: operators(:)
82
83 type(bc_mxll_t) :: bc
84 type(derivatives_t), pointer, private :: der
85 type(states_mxll_t), pointer :: st
86
87 integer :: rs_sign
88
89 logical :: propagation_apply = .false.
90
91 integer, pointer :: rs_state_fft_map(:,:,:)
92 integer, pointer :: rs_state_fft_map_inv(:,:)
93
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.
104
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.
112
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.
118
119 logical :: calc_medium_box = .false.
120 type(single_medium_box_t), allocatable :: medium_boxes(:)
121 logical :: medium_boxes_initialized = .false.
122
124 integer :: operator
125 logical :: current_density_ext_flag = .false.
126 logical :: current_density_from_medium = .false.
127
128 type(energy_mxll_t) :: energy
129
130 logical :: cpml_hamiltonian = .false.
131
132 logical :: diamag_current = .false.
133 real(real64) :: c_factor
134 real(real64) :: current_factor
135
136 type(cube_t) :: cube
137 type(mesh_cube_parallel_map_t) :: mesh_cube_map
138
139 contains
140 procedure :: update_span => hamiltonian_mxll_span
141 procedure :: dapply => dhamiltonian_mxll_apply
142 procedure :: zapply => zhamiltonian_mxll_apply
143 procedure :: is_hermitian => hamiltonian_mxll_hermitian
144 end type hamiltonian_mxll_t
145
146 integer, public, parameter :: &
147 FARADAY_AMPERE = 1, &
149 mxll_simple = 3
150
151contains
152
153 ! ---------------------------------------------------------
155 subroutine hamiltonian_mxll_init(hm, namespace, gr, st)
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
160
161
162 push_sub(hamiltonian_mxll_init)
163
164 call profiling_in('HAMILTONIAN_INIT')
165
166 hm%dim = st%dim
167 hm%st => st
168
169 assert(associated(gr%der%lapl))
170
171 hm%operators(1:3) => gr%der%grad(1:3) ! cross product for Maxwell calculation needs dimension >= 2
172 hm%der => gr%der
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.
180 !%Variable MaxwellHamiltonianOperator
181 !%Type integer
182 !%Default faraday_ampere
183 !%Section Maxwell
184 !%Description
185 !% With this variable the the Maxwell Hamiltonian operator can be selected
186 !%Option faraday_ampere 1
187 !% The propagation operation in vacuum with Spin 1 matrices without Gauss law condition.
188 !%Option faraday_ampere_medium 2
189 !% The propagation operation in medium with Spin 1 matrices without Gauss law condition
190 !%Option simple 3
191 !% A simpler implementation of the Hamiltonian, including PML. It does not use an extra
192 !% transformation to (potentially) 6-element vectors, but uses directly the Riemann-Silberstein
193 !% vector as a variable.
194 !%End
195 call parse_variable(namespace, 'MaxwellHamiltonianOperator', faraday_ampere, hm%operator)
197 if (hm%operator == faraday_ampere_medium) then
198 ! TODO: define this operator automatically if there is a linear medium system present
199 hm%dim = 2*hm%dim
200 hm%calc_medium_box = .true.
201 end if
203 !%Variable ExternalCurrent
204 !%Type logical
205 !%Default no
206 !%Section Maxwell
207 !%Description
208 !% If an external current density will be used.
209 !%End
210 call parse_variable(namespace, 'ExternalCurrent', .false., hm%current_density_ext_flag)
212 hm%plane_waves_apply = .false.
213 hm%spatial_constant_apply = .false.
214 hm%spatial_constant_propagate = .true. ! only used if spatially constant field is used
216 hm%propagation_apply = .false.
217
218 !%Variable SpeedOfLightFactor
219 !%Type float
220 !%Default 1.0
221 !%Section Maxwell
222 !%Description
223 !% Fictitous factor to modify the speed of light in vacuum.
224 !% Note: the proper way to handle linear media with a certain refractive index
225 !% is by the user of the linear_medium system.
226 !%End
227 call parse_variable(namespace, 'SpeedOfLightFactor', m_one, hm%c_factor)
229 !%Variable CurrentDensityFactor
230 !%Type float
231 !%Default 1.0
232 !%Section Maxwell
233 !%Description
234 !% Fictitous factor to modify the current density coming from partner systems.
235 !% Note: This factor does not affect the external current density prescribed by the
236 !% <tt>UserDefinedMaxwellExternalCurrent</tt> block.
237 !%End
238 call parse_variable(namespace, 'CurrentDensityFactor', m_one, hm%current_factor)
239
240 hm%rs_state_fft_map => st%rs_state_fft_map
241 hm%rs_state_fft_map_inv => st%rs_state_fft_map_inv
242
243 call parse_variable(namespace, 'StatesPack', .true., hm%apply_packed)
244
245 call parse_variable(namespace, 'TimeZero', .false., hm%time_zero)
246 if (hm%time_zero) call messages_experimental('TimeZero')
247
248 call hm%update_span(gr%spacing(1:gr%der%dim), m_zero, namespace)
249
250 call profiling_out('HAMILTONIAN_INIT')
251 pop_sub(hamiltonian_mxll_init)
252 end subroutine hamiltonian_mxll_init
253
254
255 ! ---------------------------------------------------------
256 subroutine hamiltonian_mxll_end(hm)
257 type(hamiltonian_mxll_t), intent(inout) :: hm
258
259 integer :: il
260
261 push_sub(hamiltonian_mxll_end)
262
263 call profiling_in("HAMILTONIAN_MXLL_END")
264
265 nullify(hm%operators)
266
267 call bc_mxll_end(hm%bc)
268
269 if (allocated(hm%medium_boxes)) then
270 do il = 1, size(hm%medium_boxes)
271 call single_medium_box_end(hm%medium_boxes(il))
272 end do
273 end if
274
275 call profiling_out("HAMILTONIAN_MXLL_END")
276
277 pop_sub(hamiltonian_mxll_end)
278 end subroutine hamiltonian_mxll_end
279
280
281 ! ---------------------------------------------------------
282 logical function hamiltonian_mxll_hermitian(hm)
283 class(hamiltonian_mxll_t), intent(in) :: hm
284
286
287 if (any(hm%bc%bc_ab_type == option__maxwellabsorbingboundaries__cpml)) then
288 ! With PML, the Hamiltonian is not purely Hermitian
290 else
292 end if
293
295 end function hamiltonian_mxll_hermitian
296
297
298 ! ---------------------------------------------------------
299 subroutine hamiltonian_mxll_span(hm, delta, emin, namespace)
300 class(hamiltonian_mxll_t), intent(inout) :: hm
301 real(real64), intent(in) :: delta(:)
302 real(real64), intent(in) :: emin
303 type(namespace_t), intent(in) :: namespace
304
305 integer :: i
306 real(real64) :: emax, fd_factor
307 real(real64), parameter :: fd_factors(4) = [1.0_real64, 1.3723_real64, 1.5861_real64, 1.7307_real64]
308
309 push_sub(hamiltonian_mxll_span)
310
311 ! the discretization of the gradient operator with finite differences of different order
312 ! gives different prefactors, they can be obtained by Fourier transform of the stencil
313 ! and taking the maximum of the resulting expressions
314 if (hm%der%stencil_type == der_star .and. hm%der%order <= 4) then
315 fd_factor = fd_factors(hm%der%order)
316 else
317 ! if we use a different stencil, use pi as an upper bound from the continuous solution
318 fd_factor = m_pi
319 end if
320
321 emax = m_zero
322 do i = 1, size(delta)
323 emax = emax + m_one / delta(i)**2
324 end do
325 emax = p_c * fd_factor * sqrt(emax/m_three)
326
327 hm%spectral_middle_point = m_zero
328 hm%spectral_half_span = emax
329
330 pop_sub(hamiltonian_mxll_span)
331 end subroutine hamiltonian_mxll_span
332
333
334 ! ---------------------------------------------------------
335 subroutine hamiltonian_mxll_adjoint(hm)
336 type(hamiltonian_mxll_t), intent(inout) :: hm
337
339
340 if (.not. hm%adjoint) then
341 hm%adjoint = .true.
342 end if
343
345 end subroutine hamiltonian_mxll_adjoint
346
347
348 ! ---------------------------------------------------------
349 subroutine hamiltonian_mxll_not_adjoint(hm)
350 type(hamiltonian_mxll_t), intent(inout) :: hm
353
354 if (hm%adjoint) then
355 hm%adjoint = .false.
356 end if
357
359 end subroutine hamiltonian_mxll_not_adjoint
360
361
362 ! ---------------------------------------------------------
364 subroutine hamiltonian_mxll_update(this, time)
365 type(hamiltonian_mxll_t), intent(inout) :: this
366 real(real64), optional, intent(in) :: time
367
369
370 this%current_time = m_zero
371 if (present(time)) this%current_time = time
372
374 end subroutine hamiltonian_mxll_update
375
376 ! -----------------------------------------------------------------
378 real(real64) function hamiltonian_mxll_get_time(this) result(time)
379 type(hamiltonian_mxll_t), intent(inout) :: this
380
381 time = this%current_time
382
383 end function hamiltonian_mxll_get_time
384
385 ! -----------------------------------------------------------------
386
387 logical pure function hamiltonian_mxll_apply_packed(this, mesh) result(apply)
388 type(hamiltonian_mxll_t), intent(in) :: this
389 class(mesh_t), intent(in) :: mesh
390
391 apply = this%apply_packed
392 if (mesh%use_curvilinear) apply = .false.
393
395
396 ! ---------------------------------------------------------
397 subroutine hamiltonian_mxll_apply_batch(hm, namespace, der, psib, hpsib, time, terms, set_bc)
398 type(hamiltonian_mxll_t), intent(in) :: hm
399 type(namespace_t), intent(in) :: namespace
400 type(derivatives_t), intent(in) :: der
401 type(batch_t), target, intent(inout) :: psib
402 type(batch_t), target, intent(inout) :: hpsib
403 real(real64), optional, intent(in) :: time
404 integer, optional, intent(in) :: terms
405 logical, optional, intent(in) :: set_bc
406
407 type(batch_t) :: gradb(der%dim)
408 integer :: idir, ifield, field_dir, pml_dir, rs_sign
409 integer :: ip, ip_in, il
410 real(real64) :: pml_a, pml_b, pml_c
411 complex(real64) :: pml_g, grad
412 integer, parameter :: field_dirs(3, 2) = reshape([2, 3, 1, 3, 1, 2], [3, 2])
413 real(real64) :: cc, aux_ep(3), aux_mu(3), sigma_e, sigma_m, ff_real(3), ff_imag(3), coeff_real, coeff_imag
414 complex(real64) :: ff_plus(3), ff_minus(3), hpsi(6)
415 integer, parameter :: sign_medium(6) = [1, 1, 1, -1, -1, -1]
416 logical :: with_medium
417 real(real64) :: p_c_
418
420 call profiling_in("H_MXLL_APPLY_BATCH")
421 assert(psib%status() == hpsib%status())
422
423 assert(psib%nst == hpsib%nst)
424 assert(hm%st%dim == 3)
425 p_c_ = p_c * hm%c_factor
426
427 !Not implemented at the moment
428 assert(.not. present(terms))
429 with_medium = hm%operator == faraday_ampere_medium
431 if (present(time)) then
432 if (abs(time - hm%current_time) > 1e-10_real64) then
433 write(message(1),'(a)') 'hamiltonian_apply_batch time assertion failed.'
434 write(message(2),'(a,f12.6,a,f12.6)') 'time = ', time, '; hm%current_time = ', hm%current_time
435 call messages_fatal(2, namespace=namespace)
436 end if
437 end if
438
439 if (hm%operator == mxll_simple) then
440 call hamiltonian_mxll_apply_simple(hm, namespace, der%mesh, psib, hpsib, terms, set_bc)
441 call profiling_out("H_MXLL_APPLY_BATCH")
443 return
444 end if
445
446 call zderivatives_batch_grad(der, psib, gradb, set_bc=set_bc)
447
448 if (hm%cpml_hamiltonian) then
449 call apply_pml_boundary(gradb)
450 end if
451
452 call zderivatives_batch_curl_from_gradient(der, hpsib, gradb)
453
454 call scale_after_curl()
455
456 if (hm%bc_constant .and. .not. with_medium) then
458 end if
460 if (with_medium) then
461 do idir = 1, 3
462 if ((hm%bc%bc_type(idir) == mxll_bc_medium)) then
463 call apply_medium_box(hm%bc%medium(idir))
464 end if
465 end do
466
467 if (hm%calc_medium_box) then
468 do il = 1, size(hm%medium_boxes)
469 call apply_medium_box(hm%medium_boxes(il))
470 end do
471 end if
472 end if
474 do idir = 1, der%dim
475 call gradb(idir)%end()
476 end do
477
478 call profiling_out("H_MXLL_APPLY_BATCH")
480
481 contains
482 subroutine apply_pml_boundary(gradb)
483 type(batch_t) :: gradb(:)
484 type(accel_kernel_t), save :: ker_pml
485 integer :: bsize, gsize
487 call profiling_in("APPLY_PML_BOUNDARY")
488
489 if (with_medium) then
490 rs_sign = 1
491 else
492 rs_sign = hm%rs_sign
493 end if
494 do idir = 1, der%dim
495 call batch_scal(der%mesh%np, rs_sign * p_c_, gradb(idir))
496 end do
497
498 do pml_dir = 1, hm%st%dim
499 if (hm%bc%bc_ab_type(pml_dir) == mxll_ab_cpml) then
500 select case (gradb(pml_dir)%status())
501 case (batch_not_packed)
502 do ifield = 1, 2
503 field_dir = field_dirs(pml_dir, ifield)
504 !$omp parallel do private(ip, pml_c, pml_a, pml_b, pml_g, grad)
505 do ip_in = 1, hm%bc%pml%points_number
506 ip = hm%bc%pml%points_map(ip_in)
507 pml_c = hm%bc%pml%c(ip_in, pml_dir)
508 pml_a = hm%bc%pml%a(ip_in, pml_dir)
509 pml_b = hm%bc%pml%b(ip_in, pml_dir)
510 pml_g = hm%bc%pml%conv_plus(ip_in, pml_dir, field_dir)
511 grad = gradb(pml_dir)%zff_linear(ip, field_dir)
512 gradb(pml_dir)%zff_linear(ip, field_dir) = pml_c * ((m_one + pml_a)*grad/p_c_ &
513 + rs_sign * pml_b*pml_g)
514 if (with_medium) then
515 pml_g = hm%bc%pml%conv_minus(ip_in, pml_dir, field_dir)
516 grad = gradb(pml_dir)%zff_linear(ip, field_dir+3)
517 gradb(pml_dir)%zff_linear(ip, field_dir+3) = pml_c * ((m_one + pml_a)*grad/p_c_ &
518 + rs_sign * pml_b*pml_g)
519 end if
520 end do
521 end do
522 case (batch_packed)
523 !$omp parallel do private(ip, field_dir, pml_c, pml_a, pml_b, pml_g, grad, ifield)
524 do ip_in = 1, hm%bc%pml%points_number
525 ip = hm%bc%pml%points_map(ip_in)
526 pml_c = hm%bc%pml%c(ip_in, pml_dir)
527 pml_a = hm%bc%pml%a(ip_in, pml_dir)
528 pml_b = hm%bc%pml%b(ip_in, pml_dir)
529 do ifield = 1, 2
530 field_dir = field_dirs(pml_dir, ifield)
531 pml_g = hm%bc%pml%conv_plus(ip_in, pml_dir, field_dir)
532 grad = gradb(pml_dir)%zff_pack(field_dir, ip)
533 gradb(pml_dir)%zff_pack(field_dir, ip) = pml_c * ((m_one + pml_a)*grad/p_c_ &
534 + rs_sign * pml_b*pml_g)
535 if (with_medium) then
536 pml_g = hm%bc%pml%conv_minus(ip_in, pml_dir, field_dir)
537 grad = gradb(pml_dir)%zff_pack(field_dir+3, ip)
538 gradb(pml_dir)%zff_pack(field_dir+3, ip) = pml_c * ((m_one + pml_a)*grad/p_c_ &
539 + rs_sign * pml_b*pml_g)
540 end if
541 end do
542 end do
543 case (batch_device_packed)
544 call accel_kernel_start_call(ker_pml, 'pml.cu', 'pml_apply')
545
546 if (with_medium) then
547 call accel_set_kernel_arg(ker_pml, 0, 1_int32)
548 else
549 call accel_set_kernel_arg(ker_pml, 0, 0_int32)
550 end if
551 call accel_set_kernel_arg(ker_pml, 1, hm%bc%pml%points_number)
552 call accel_set_kernel_arg(ker_pml, 2, pml_dir-1)
553 call accel_set_kernel_arg(ker_pml, 3, p_c_)
554 call accel_set_kernel_arg(ker_pml, 4, rs_sign)
555 call accel_set_kernel_arg(ker_pml, 5, hm%bc%pml%buff_map)
556 call accel_set_kernel_arg(ker_pml, 6, gradb(pml_dir)%ff_device)
557 call accel_set_kernel_arg(ker_pml, 7, log2(int(gradb(pml_dir)%pack_size(1), int32)))
558 call accel_set_kernel_arg(ker_pml, 8, hm%bc%pml%buff_a)
559 call accel_set_kernel_arg(ker_pml, 9, hm%bc%pml%buff_b)
560 call accel_set_kernel_arg(ker_pml, 10, hm%bc%pml%buff_c)
561 call accel_set_kernel_arg(ker_pml, 11, hm%bc%pml%buff_conv_plus)
562 call accel_set_kernel_arg(ker_pml, 12, hm%bc%pml%buff_conv_minus)
563
564 ! Compute the grid size
565 bsize = accel_max_block_size()
566 call accel_grid_size(hm%bc%pml%points_number, bsize, gsize)
567
568 call accel_kernel_run(ker_pml, (/gsize/), (/bsize/))
569 end select
570 end if
571 end do
572
573 if (accel_is_enabled()) then
574 call accel_finish()
575 end if
576
577 call profiling_out("APPLY_PML_BOUNDARY")
579 end subroutine apply_pml_boundary
580
581 subroutine scale_after_curl
583 call profiling_in("SCALE_AFTER_CURL")
584 if (.not. hm%cpml_hamiltonian) then
585 ! if we do not need pml, scale after the curl because it is cheaper
586 if (with_medium) then
587 ! in case of a medium, multiply first 3 components with +, others with -
588 call batch_scal(der%mesh%np, sign_medium * p_c_, hpsib)
589 else
590 call batch_scal(der%mesh%np, hm%rs_sign * p_c_, hpsib)
591 end if
592 else
593 ! this is needed for PML computations with medium
594 if (with_medium) then
595 ! in case of a medium, multiply first 3 components with +, others with -
596 call batch_scal(der%mesh%np, sign_medium * m_one, hpsib)
597 end if
598 end if
599 call profiling_out("SCALE_AFTER_CURL")
601 end subroutine scale_after_curl
602
603 subroutine apply_constant_boundary
605 call profiling_in('APPLY_CONSTANT_BC')
606 select case (hpsib%status())
607 case (batch_not_packed)
608 do field_dir = 1, hm%st%dim
609 do ip_in = 1, hm%bc%constant_points_number
610 ip = hm%bc%constant_points_map(ip_in)
611 hpsib%zff_linear(ip, field_dir) = hm%st%rs_state_const(field_dir)
612 end do
613 end do
614 case (batch_packed)
615 do ip_in = 1, hm%bc%constant_points_number
616 ip = hm%bc%constant_points_map(ip_in)
617 do field_dir = 1, hm%st%dim
618 hpsib%zff_pack(field_dir, ip) = hm%st%rs_state_const(field_dir)
619 end do
620 end do
621 case (batch_device_packed)
622 call messages_not_implemented("Maxwell constant boundary on GPU", namespace=namespace)
623 end select
624 call profiling_out('APPLY_CONSTANT_BC')
626 end subroutine apply_constant_boundary
627
628 subroutine apply_medium_box(medium)
629 type(single_medium_box_t), intent(in) :: medium
630
631 integer :: ifield
632
634 call profiling_in("MEDIUM_BOX")
635 assert(.not. medium%has_mapping)
636 !$omp parallel do private(ip, cc, aux_ep, aux_mu, sigma_e, sigma_m, &
637 !$omp ff_plus, ff_minus, hpsi, ff_real, ff_imag, ifield, coeff_real, coeff_imag)
638 do ip = 1, medium%points_number
639 cc = medium%c(ip)/p_c
640 if (abs(cc) <= m_epsilon) cycle
641 aux_ep(1:3) = medium%aux_ep(ip, 1:3)
642 aux_mu(1:3) = medium%aux_mu(ip, 1:3)
643 sigma_e = medium%sigma_e(ip)
644 sigma_m = medium%sigma_m(ip)
645 select case (hpsib%status())
646 case (batch_not_packed)
647 ff_plus(1:3) = psib%zff_linear(ip, 1:3)
648 ff_minus(1:3) = psib%zff_linear(ip, 4:6)
649 hpsi(1:6) = hpsib%zff_linear(ip, 1:6)
650 case (batch_packed)
651 ff_plus(1:3) = psib%zff_pack(1:3, ip)
652 ff_minus(1:3) = psib%zff_pack(4:6, ip)
653 hpsi(1:6) = hpsib%zff_pack(1:6, ip)
654 case (batch_device_packed)
655 call messages_not_implemented("Maxwell Medium on GPU", namespace=namespace)
656 end select
657 ff_real = real(ff_plus+ff_minus, real64)
658 ff_imag = aimag(ff_plus-ff_minus)
659 aux_ep = dcross_product(aux_ep, ff_real)
660 aux_mu = dcross_product(aux_mu, ff_imag)
661 do ifield = 1, 3
662 coeff_real = - cc * aux_ep(ifield) + sigma_m * ff_imag(ifield)
663 coeff_imag = - cc * aux_mu(ifield) - sigma_e * ff_real(ifield)
664 hpsi(ifield) = cc * hpsi(ifield) + cmplx(coeff_real, coeff_imag, real64)
665 hpsi(ifield+3) = cc * hpsi(ifield+3) + cmplx(-coeff_real, coeff_imag, real64)
666 end do
667 select case (hpsib%status())
668 case (batch_not_packed)
669 hpsib%zff_linear(ip, 1:6) = hpsi(1:6)
670 case (batch_packed)
671 hpsib%zff_pack(1:6, ip) = hpsi(1:6)
672 end select
673 end do
674 call profiling_out("MEDIUM_BOX")
676 end subroutine apply_medium_box
677
678 end subroutine hamiltonian_mxll_apply_batch
679
680 ! ---------------------------------------------------------
681 subroutine hamiltonian_mxll_apply_simple(hm, namespace, mesh, psib, hpsib, terms, set_bc)
682 type(hamiltonian_mxll_t), intent(in) :: hm
683 type(namespace_t), intent(in) :: namespace
684 class(mesh_t), intent(in) :: mesh
685 type(batch_t), target, intent(inout) :: psib
686 type(batch_t), target, intent(inout) :: hpsib
687 integer, optional, intent(in) :: terms
688 logical, optional, intent(in) :: set_bc
689
690 type(batch_t) :: gradb(hm%der%dim)
691 integer :: idir
692
694 call profiling_in("MXLL_HAMILTONIAN_SIMPLE")
695
696 call zderivatives_batch_grad(hm%der, psib, gradb, set_bc=set_bc)
697
698 if (any(hm%bc%bc_ab_type == option__maxwellabsorbingboundaries__cpml)) then
699 ! in the pml layer, the gradient is modified
700 call mxll_apply_pml_simple(hm, gradb)
701 end if
702
703 ! compute the curl from the gradient
704 call zderivatives_batch_curl_from_gradient(hm%der, hpsib, gradb)
705
706 if (hm%calc_medium_box) then
707 ! as the speed of light depends on space, hpsib is already scaled correctly
708 call mxll_linear_medium_terms_simple(hm, hpsib)
709 else
710 ! scale rs_state_tmpb (prefactor of curl)
711 call batch_scal(mesh%np, p_c*hm%c_factor, hpsib)
712 end if
713
714 do idir = 1, hm%der%dim
715 call gradb(idir)%end()
716 end do
717
718 call profiling_out("MXLL_HAMILTONIAN_SIMPLE")
720 end subroutine hamiltonian_mxll_apply_simple
721
722 ! ---------------------------------------------------------
723 subroutine mxll_apply_pml_simple(hm, gradb)
724 type(hamiltonian_mxll_t), target, intent(in) :: hm
725 type(batch_t), intent(inout) :: gradb(1:hm%st%dim)
726
727 integer :: idir, jdir, ip_in, ip, bsize, gsize
728 type(accel_kernel_t), save :: ker_pml
729
730
731 push_sub(mxll_apply_pml_simple)
732 call profiling_in("APPLY_PML_SIMPLE")
733 ! update pml values
734 ! loop over gradient direction
735 do jdir = 1, hm%st%dim
736 select case (gradb(1)%status())
737 case (batch_not_packed)
738 ! loop over space direction
739 do idir = 1, hm%st%dim
740 if (idir == jdir) cycle
741 do ip_in = 1, hm%bc%pml%points_number
742 ip = hm%bc%pml%points_map(ip_in)
743 ! update gradient to take pml into account when computing the curl
744 gradb(jdir)%zff_linear(ip, idir) = gradb(jdir)%zff_linear(ip, idir) * (m_one + hm%bc%pml%c(ip_in, jdir)) + &
745 hm%bc%pml%b(ip_in, jdir) * hm%bc%pml%conv_plus_old(ip_in, jdir, idir)
746 end do
747 end do
748 case (batch_packed)
749 do ip_in = 1, hm%bc%pml%points_number
750 ip = hm%bc%pml%points_map(ip_in)
751 ! loop over space direction
752 do idir = 1, hm%st%dim
753 if (idir == jdir) cycle
754 ! update gradient to take pml into account when computing the curl
755 gradb(jdir)%zff_pack(idir, ip) = gradb(jdir)%zff_pack(idir, ip) * (m_one + hm%bc%pml%c(ip_in, jdir)) +&
756 hm%bc%pml%b(ip_in, jdir) * hm%bc%pml%conv_plus_old(ip_in, jdir, idir)
757 end do
758 end do
759 case (batch_device_packed)
760 call accel_kernel_start_call(ker_pml, 'pml.cu', 'pml_apply_new')
761
762 call accel_set_kernel_arg(ker_pml, 0, hm%bc%pml%points_number)
763 call accel_set_kernel_arg(ker_pml, 1, jdir-1)
764 call accel_set_kernel_arg(ker_pml, 2, hm%bc%pml%buff_map)
765 call accel_set_kernel_arg(ker_pml, 3, gradb(jdir)%ff_device)
766 call accel_set_kernel_arg(ker_pml, 4, log2(int(gradb(jdir)%pack_size(1), int32)))
767 call accel_set_kernel_arg(ker_pml, 5, hm%bc%pml%buff_b)
768 call accel_set_kernel_arg(ker_pml, 6, hm%bc%pml%buff_c)
769 call accel_set_kernel_arg(ker_pml, 7, hm%bc%pml%buff_conv_plus_old)
770
771 ! Compute the grid size
772 bsize = accel_max_block_size()
773 call accel_grid_size(hm%bc%pml%points_number, bsize, gsize)
774
775 call accel_kernel_run(ker_pml, (/gsize/), (/bsize/))
776 end select
777 end do
778 call profiling_out("APPLY_PML_SIMPLE")
779 pop_sub(mxll_apply_pml_simple)
780 end subroutine mxll_apply_pml_simple
781
782 ! ---------------------------------------------------------
783 subroutine mxll_update_pml_simple(hm, rs_stateb)
784 type(hamiltonian_mxll_t),intent(inout) :: hm
785 type(batch_t), intent(inout) :: rs_stateb
786
787 integer :: idir, jdir, ip, ip_in, bsize, gsize
788 type(accel_kernel_t), save :: ker_pml_update
789 type(batch_t) :: gradb(1:hm%st%dim)
790
791 push_sub(mxll_update_pml_simple)
792 call profiling_in("UPDATE_PML_SIMPLE")
793
794 call zderivatives_batch_grad(hm%der, rs_stateb, gradb)
795
796 do jdir = 1, hm%st%dim
797 call rs_stateb%check_compatibility_with(gradb(jdir))
798 select case (gradb(jdir)%status())
799 case (batch_not_packed)
800 ! loop over space direction
801 do idir = 1, hm%st%dim
802 if (idir == jdir) cycle
803 do ip_in = 1, hm%bc%pml%points_number
804 ip = hm%bc%pml%points_map(ip_in)
805 hm%bc%pml%conv_plus(ip_in, jdir, idir) = hm%bc%pml%c(ip_in, jdir) * gradb(jdir)%zff_linear(ip, idir) +&
806 hm%bc%pml%b(ip_in, jdir) * hm%bc%pml%conv_plus_old(ip_in, jdir, idir)
807 end do
808 end do
809 case (batch_packed)
810 do ip_in = 1, hm%bc%pml%points_number
811 ip = hm%bc%pml%points_map(ip_in)
812 ! loop over space direction
813 do idir = 1, hm%st%dim
814 if (idir == jdir) cycle
815 hm%bc%pml%conv_plus(ip_in, jdir, idir) = hm%bc%pml%c(ip_in, jdir) * gradb(jdir)%zff_pack(idir, ip) +&
816 hm%bc%pml%b(ip_in, jdir) * hm%bc%pml%conv_plus_old(ip_in, jdir, idir)
817 end do
818 end do
819 case (batch_device_packed)
820 call accel_kernel_start_call(ker_pml_update, 'pml.cu', 'pml_update_new')
821
822 call accel_set_kernel_arg(ker_pml_update, 0, hm%bc%pml%points_number)
823 call accel_set_kernel_arg(ker_pml_update, 1, jdir-1)
824 call accel_set_kernel_arg(ker_pml_update, 2, hm%bc%pml%buff_map)
825 call accel_set_kernel_arg(ker_pml_update, 3, gradb(jdir)%ff_device)
826 call accel_set_kernel_arg(ker_pml_update, 4, log2(int(gradb(jdir)%pack_size(1), int32)))
827 call accel_set_kernel_arg(ker_pml_update, 5, hm%bc%pml%buff_b)
828 call accel_set_kernel_arg(ker_pml_update, 6, hm%bc%pml%buff_c)
829 call accel_set_kernel_arg(ker_pml_update, 7, hm%bc%pml%buff_conv_plus)
830 call accel_set_kernel_arg(ker_pml_update, 8, hm%bc%pml%buff_conv_plus_old)
831
832 ! Compute the grid size
833 bsize = accel_max_block_size()
834 call accel_grid_size(hm%bc%pml%points_number, bsize, gsize)
835
836 call accel_kernel_run(ker_pml_update, (/gsize/), (/bsize/))
837 end select
838 call gradb(jdir)%end()
839 end do
840 call profiling_out("UPDATE_PML_SIMPLE")
842 end subroutine mxll_update_pml_simple
843
844 ! ---------------------------------------------------------
845 subroutine mxll_copy_pml_simple(hm, rs_stateb)
846 type(hamiltonian_mxll_t),intent(inout) :: hm
847 type(batch_t), intent(inout) :: rs_stateb
848
849 integer :: bsize, idir, jdir, ip, ip_in
850 integer :: gsize
851 type(accel_kernel_t), save :: ker_pml_copy
852
853 push_sub(mxll_copy_pml_simple)
854 call profiling_in("COPY_PML_SIMPLE")
855
856 select case (rs_stateb%status())
857 case (batch_packed, batch_not_packed)
858 do jdir = 1, hm%st%dim
859 ! loop over space direction
860 do idir = 1, hm%st%dim
861 if (idir == jdir) cycle
862 do ip_in = 1, hm%bc%pml%points_number
863 ip = hm%bc%pml%points_map(ip_in)
864 hm%bc%pml%conv_plus_old(ip_in, jdir, idir) = hm%bc%pml%conv_plus(ip_in, jdir, idir)
865 end do
866 end do
867 end do
868 case (batch_device_packed)
869 call accel_kernel_start_call(ker_pml_copy, 'pml.cu', 'pml_copy')
870
871 call accel_set_kernel_arg(ker_pml_copy, 0, hm%bc%pml%points_number*9)
872 call accel_set_kernel_arg(ker_pml_copy, 1, hm%bc%pml%buff_conv_plus)
873 call accel_set_kernel_arg(ker_pml_copy, 2, hm%bc%pml%buff_conv_plus_old)
874
875 ! Compute the grid size
876 bsize = accel_max_block_size()
877 call accel_grid_size(hm%bc%pml%points_number*9, bsize, gsize)
879 call accel_kernel_run(ker_pml_copy, (/gsize/), (/bsize/))
880 end select
881 call profiling_out("COPY_PML_SIMPLE")
882 pop_sub(mxll_copy_pml_simple)
883 end subroutine mxll_copy_pml_simple
884
885 ! ---------------------------------------------------------
886 subroutine mxll_linear_medium_terms_simple(hm, rs_stateb)
887 type(hamiltonian_mxll_t),intent(in) :: hm
888 type(batch_t), intent(inout) :: rs_stateb
889
890 integer :: il, ip
891 logical :: need_to_pack
892
894 call profiling_in("LINEAR_MEDIUM_SIMPLE")
895
896 do il = 1, size(hm%medium_boxes)
897 assert(.not. hm%medium_boxes(il)%has_mapping)
898 need_to_pack = .false.
899 ! copy to the CPU for GPU runs
900 if(rs_stateb%status() == batch_device_packed) then
901 call rs_stateb%do_unpack(force=.true.)
902 need_to_pack = .true.
903 end if
904 select case (rs_stateb%status())
905 case (batch_not_packed)
906 do ip = 1, hm%medium_boxes(il)%points_number
907 if (abs(hm%medium_boxes(il)%c(ip)) <= m_epsilon) then
908 ! Hamiltonian without medium
909 rs_stateb%zff_linear(ip, 1:3) = p_c*hm%c_factor*rs_stateb%zff_linear(ip, 1:3)
910 else
911 ! Hamiltonian with medium terms
912 rs_stateb%zff_linear(ip, 1:3) = hm%medium_boxes(il)%c(ip)*(rs_stateb%zff_linear(ip, 1:3) + &
913 cmplx( &
914 -dcross_product(hm%medium_boxes(il)%aux_ep(ip, 1:3)*m_two, &
915 real(rs_stateb%zff_linear(ip, 1:3), real64)), &
916 -dcross_product(hm%medium_boxes(il)%aux_mu(ip, 1:3)*m_two, &
917 aimag(rs_stateb%zff_linear(ip, 1:3))), real64)) + &
918 cmplx(&
919 hm%medium_boxes(il)%sigma_m(ip)*aimag(rs_stateb%zff_linear(ip, 1:3)), &
920 -hm%medium_boxes(il)%sigma_e(ip)*real(rs_stateb%zff_linear(ip, 1:3), real64), real64)
921 end if
922 end do
923 case (batch_packed)
924 do ip = 1, hm%medium_boxes(il)%points_number
925 if (abs(hm%medium_boxes(il)%c(ip)) <= m_epsilon) then
926 ! Hamiltonian without medium
927 rs_stateb%zff_pack(1:3, ip) = p_c*hm%c_factor*rs_stateb%zff_pack(1:3, ip)
928 else
929 ! Hamiltonian with medium terms
930 rs_stateb%zff_pack(1:3, ip) = hm%medium_boxes(il)%c(ip)*(rs_stateb%zff_pack(1:3, ip) + &
931 cmplx( &
932 -dcross_product(hm%medium_boxes(il)%aux_ep(ip, 1:3)*m_two, &
933 real(rs_stateb%zff_pack(1:3, ip), real64)), &
934 -dcross_product(hm%medium_boxes(il)%aux_mu(ip, 1:3)*m_two, &
935 aimag(rs_stateb%zff_pack(1:3, ip))), real64)) + &
936 cmplx(&
937 hm%medium_boxes(il)%sigma_m(ip)*aimag(rs_stateb%zff_pack(1:3, ip)), &
938 -hm%medium_boxes(il)%sigma_e(ip)*real(rs_stateb%zff_pack(1:3, ip), real64), real64)
939 end if
940 end do
941 end select
942 if(need_to_pack) then
943 call rs_stateb%do_pack()
944 end if
945 end do
946
947 call profiling_out("LINEAR_MEDIUM_SIMPLE")
950
951 ! --------------------------------------------------------
953 subroutine dhamiltonian_mxll_apply(hm, namespace, mesh, psib, hpsib, terms, set_bc)
954 class(hamiltonian_mxll_t), intent(in) :: hm
955 type(namespace_t), intent(in) :: namespace
956 class(mesh_t), intent(in) :: mesh
957 class(batch_t), target, intent(inout) :: psib
958 class(batch_t), target, intent(inout) :: hpsib
959 integer, optional, intent(in) :: terms
960 logical, optional, intent(in) :: set_bc
961
962 write(message(1),'(a)') 'dhamiltonian_mxll_apply not implemented (states are complex).'
963 call messages_fatal(1, namespace=namespace)
964
965 end subroutine dhamiltonian_mxll_apply
966
967 ! ---------------------------------------------------------
969 subroutine zhamiltonian_mxll_apply(hm, namespace, mesh, psib, hpsib, terms, set_bc)
970 class(hamiltonian_mxll_t), intent(in) :: hm
971 type(namespace_t), intent(in) :: namespace
972 class(mesh_t), intent(in) :: mesh
973 class(batch_t), target, intent(inout) :: psib
974 class(batch_t), target, intent(inout) :: hpsib
975 integer, optional, intent(in) :: terms
976 logical, optional, intent(in) :: set_bc
977
978 complex(real64), allocatable :: rs_aux_in(:,:), rs_aux_out(:,:)
979 integer :: ii
980 logical :: on_gpu
983
984 call profiling_in('ZHAMILTONIAN_MXLL_APPLY')
985
986 on_gpu = psib%status() == batch_device_packed
987 if (hm%operator == faraday_ampere_medium .and. on_gpu) then
988 ! legacy code, keep for the moment
989 safe_allocate(rs_aux_in(1:mesh%np_part, 1:hm%dim))
990 safe_allocate(rs_aux_out(1:mesh%np, 1:hm%dim))
991 call boundaries_set(hm%der%boundaries, mesh, psib)
992 do ii = 1, hm%dim
993 call batch_get_state(psib, ii,mesh%np_part, rs_aux_in(:, ii))
994 end do
995 ! This uses the old non-batch implementation
996 call maxwell_hamiltonian_apply_fd(hm, hm%der, rs_aux_in, rs_aux_out)
997 do ii = 1, hm%dim
998 call batch_set_state(hpsib, ii, mesh%np, rs_aux_out(:, ii))
999 end do
1000 safe_deallocate_a(rs_aux_in)
1001 safe_deallocate_a(rs_aux_out)
1002 else
1003 ! default branch, should be the only one at some point
1004 call hamiltonian_mxll_apply_batch(hm, namespace, hm%der, psib, hpsib, set_bc=set_bc)
1005 end if
1006
1007 call profiling_out('ZHAMILTONIAN_MXLL_APPLY')
1008
1010 end subroutine zhamiltonian_mxll_apply
1011
1012
1013 ! ---------------------------------------------------------
1015 subroutine maxwell_hamiltonian_apply_fd(hm, der, psi, oppsi)
1016 type(hamiltonian_mxll_t), intent(in) :: hm
1017 type(derivatives_t), intent(in) :: der
1018 complex(real64), contiguous, intent(inout) :: psi(:,:)
1019 complex(real64), intent(inout) :: oppsi(:,:)
1020
1021 real(real64), pointer :: mx_rho(:,:)
1022 complex(real64), allocatable :: tmp(:,:)
1023 complex(real64), pointer :: kappa_psi(:,:)
1024 integer :: np, np_part, ip, ip_in, rs_sign
1025 real(real64) :: P_c_
1026
1028
1029 call profiling_in('MAXWELL_HAMILTONIAN_APPLY_FD')
1030
1031 np = der%mesh%np
1032 np_part = der%mesh%np_part
1033 rs_sign = hm%rs_sign
1034 p_c_ = p_c * hm%c_factor
1035
1036 select case (hm%operator)
1037
1038 !=================================================================================================
1039 ! Maxwell Hamiltonian - Hamiltonian operation in vacuum via partial derivatives:
1040
1041 case (faraday_ampere)
1042 call profiling_in('MXLL_HAM_APPLY_FD_FARADAY_AMP')
1043
1044 safe_allocate(tmp(1:np,1:2))
1045 oppsi = m_z0
1046
1047 if (hm%diamag_current) then
1048 mx_rho => hm%st%grid_rho
1049 kappa_psi => hm%st%kappa_psi
1050 end if
1051
1052 !-----------------------------------------------------------------------------------------------
1053 ! Riemann-Silberstein vector component 1 calculation:
1054 tmp = m_z0
1055 call zderivatives_partial(der, psi(:,3), tmp(:,1), 2, set_bc = .false.)
1056 call zderivatives_partial(der, psi(:,2), tmp(:,2), 3, set_bc = .false.)
1057 tmp = rs_sign * p_c_ * tmp
1058 call maxwell_pml_hamiltonian(hm, der, psi, 2, 3, tmp(:,1))
1059 call maxwell_pml_hamiltonian(hm, der, psi, 3, 2, tmp(:,2))
1060 oppsi(1:np,1) = ( tmp(1:np,1)-tmp(1:np,2))
1061
1062 !-----------------------------------------------------------------------------------------------
1063 ! Riemann-Silberstein vector component 2 calculation:
1064 tmp = m_z0
1065 call zderivatives_partial(der, psi(:,1), tmp(:,1), 3, set_bc = .false.)
1066 call zderivatives_partial(der, psi(:,3), tmp(:,2), 1, set_bc = .false.)
1067 tmp = rs_sign * p_c_ * tmp
1068 call maxwell_pml_hamiltonian(hm, der, psi, 3, 1, tmp(:,1))
1069 call maxwell_pml_hamiltonian(hm, der, psi, 1, 3, tmp(:,2))
1070 oppsi(1:np,2) = ( tmp(1:np,1)-tmp(1:np,2))
1071
1072 !-----------------------------------------------------------------------------------------------
1073 ! Riemann-Silberstein vector component 3 calculation:
1074 tmp = m_z0
1075 call zderivatives_partial(der, psi(:,2), tmp(:,1), 1, set_bc = .false.)
1076 call zderivatives_partial(der, psi(:,1), tmp(:,2), 2, set_bc = .false.)
1077 tmp = rs_sign * p_c_ * tmp
1078 call maxwell_pml_hamiltonian(hm, der, psi, 1, 2, tmp(:,1))
1079 call maxwell_pml_hamiltonian(hm, der, psi, 2, 1, tmp(:,2))
1080 oppsi(1:np,3) = ( tmp(1:np,1)-tmp(1:np,2))
1081
1082 if (hm%bc_constant) then
1083 do ip_in = 1, hm%bc%constant_points_number
1084 ip = hm%bc%constant_points_map(ip_in)
1085 oppsi(ip,:) = hm%st%rs_state_const(:)
1086 end do
1087 end if
1088
1089 safe_deallocate_a(tmp)
1090
1091 call profiling_out('MXLL_HAM_APPLY_FD_FARADAY_AMP')
1092 !=================================================================================================
1093 ! Maxwell Hamiltonian - Hamiltonian operation in medium via partial derivatives:
1094
1096 call profiling_in('MXLL_HAM_APP_FAR_AMP_MED')
1097
1098 safe_allocate(tmp(1:np,1:4))
1099 oppsi = m_z0
1100
1101 !-----------------------------------------------------------------------------------------------
1102 ! Riemann-Silberstein vector component 1 calculation:
1103 tmp = m_z0
1104 call zderivatives_partial(der, psi(:,3), tmp(:,1), 2, set_bc = .false.)
1105 call zderivatives_partial(der, psi(:,2), tmp(:,3), 3, set_bc = .false.)
1106 call zderivatives_partial(der, psi(:,6), tmp(:,2), 2, set_bc = .false.)
1107 call zderivatives_partial(der, psi(:,5), tmp(:,4), 3, set_bc = .false.)
1108 tmp = p_c_ * tmp
1109 call maxwell_pml_hamiltonian_medium(hm, der, psi, 2, 3, tmp(:,1:2))
1110 call maxwell_pml_hamiltonian_medium(hm, der, psi, 3, 2, tmp(:,3:4))
1111 oppsi(1:np,1) = (tmp(1:np,1)-tmp(1:np,3))
1112 oppsi(1:np,4) = -(tmp(1:np,2)-tmp(1:np,4))
1113
1114 !-----------------------------------------------------------------------------------------------
1115 ! Riemann-Silberstein vector component 2 calculation:
1116 tmp = m_z0
1117 call zderivatives_partial(der, psi(:,1), tmp(:,1), 3, set_bc = .false.)
1118 call zderivatives_partial(der, psi(:,3), tmp(:,3), 1, set_bc = .false.)
1119 call zderivatives_partial(der, psi(:,4), tmp(:,2), 3, set_bc = .false.)
1120 call zderivatives_partial(der, psi(:,6), tmp(:,4), 1, set_bc = .false.)
1121 tmp = p_c_ * tmp
1122 call maxwell_pml_hamiltonian_medium(hm, der, psi, 3, 1, tmp(:,1:2))
1123 call maxwell_pml_hamiltonian_medium(hm, der, psi, 1, 3, tmp(:,3:4))
1124 oppsi(1:np,2) = (tmp(1:np,1)-tmp(1:np,3))
1125 oppsi(1:np,5) = -(tmp(1:np,2)-tmp(1:np,4))
1126
1127 !-----------------------------------------------------------------------------------------------
1128 ! Riemann-Silberstein vector component 3 calculation:
1129 tmp = m_z0
1130 call zderivatives_partial(der, psi(:,2), tmp(:,1), 1, set_bc = .false.)
1131 call zderivatives_partial(der, psi(:,1), tmp(:,3), 2, set_bc = .false.)
1132 call zderivatives_partial(der, psi(:,5), tmp(:,2), 1, set_bc = .false.)
1133 call zderivatives_partial(der, psi(:,4), tmp(:,4), 2, set_bc = .false.)
1134 tmp = p_c_ * tmp
1135 call maxwell_pml_hamiltonian_medium(hm, der, psi, 1, 2, tmp(:,1:2))
1136 call maxwell_pml_hamiltonian_medium(hm, der, psi, 2, 1, tmp(:,3:4))
1137 oppsi(1:np,3) = (tmp(1:np,1)-tmp(1:np,3))
1138 oppsi(1:np,6) = -(tmp(1:np,2)-tmp(1:np,4))
1139
1140
1141 safe_deallocate_a(tmp)
1142
1143 !-----------------------------------------------------------------------------------------------
1144 ! Riemann-Silberstein vector calculation if medium boundaries is set:
1145 call maxwell_medium_boundaries_calculation(hm, psi, oppsi)
1146
1147 !-----------------------------------------------------------------------------------------------
1148 ! Riemann-Silberstein vector calculation for medium boxes:
1149 call maxwell_medium_boxes_calculation(hm, der, psi, oppsi)
1150
1151 call profiling_out('MXLL_HAM_APP_FAR_AMP_MED')
1152
1153 end select
1154
1155 call profiling_out('MAXWELL_HAMILTONIAN_APPLY_FD')
1156
1158 end subroutine maxwell_hamiltonian_apply_fd
1159
1160
1161 ! ---------------------------------------------------------
1163 subroutine maxwell_pml_hamiltonian(hm, der, psi, dir1, dir2, tmp)
1164 type(hamiltonian_mxll_t), intent(in) :: hm
1165 type(derivatives_t), intent(in) :: der
1166 complex(real64), intent(inout) :: psi(:,:)
1167 integer, intent(in) :: dir1
1168 integer, intent(in) :: dir2
1169 complex(real64), intent(inout) :: tmp(:)
1170
1171
1172 push_sub(maxwell_pml_hamiltonian)
1173
1174 call profiling_in('MAXWELL_PML_HAMILTONIAN')
1175
1176 if ((hm%bc%bc_ab_type(dir1) == mxll_ab_cpml) .and. hm%cpml_hamiltonian) then
1177 call maxwell_pml_calculation_via_riemann_silberstein(hm, der, psi, dir1, dir2, tmp(:))
1178 end if
1179
1180 call profiling_out('MAXWELL_PML_HAMILTONIAN')
1181
1183 end subroutine maxwell_pml_hamiltonian
1184
1185 ! ---------------------------------------------------------
1187 subroutine maxwell_pml_hamiltonian_medium(hm, der, psi, dir1, dir2, tmp)
1188 type(hamiltonian_mxll_t), intent(in) :: hm
1189 type(derivatives_t), intent(in) :: der
1190 complex(real64), intent(inout) :: psi(:,:)
1191 integer, intent(in) :: dir1
1192 integer, intent(in) :: dir2
1193 complex(real64), intent(inout) :: tmp(:,:)
1194
1195
1197
1198 call profiling_in('MAXWELL_PML_HAMILTONIAN_MEDIUM')
1199
1200 if ((hm%bc%bc_ab_type(dir1) == mxll_ab_cpml) .and. hm%cpml_hamiltonian) then
1201 call maxwell_pml_calculation_via_riemann_silberstein_medium(hm, der, psi, dir1, dir2, tmp(:,:))
1202 end if
1203
1204 call profiling_out('MAXWELL_PML_HAMILTONIAN_MEDIUM')
1205
1207 end subroutine maxwell_pml_hamiltonian_medium
1208
1209 ! ---------------------------------------------------------
1211 subroutine maxwell_pml_calculation_via_riemann_silberstein(hm, der, psi, pml_dir, field_dir, pml)
1212 type(hamiltonian_mxll_t), intent(in) :: hm
1213 type(derivatives_t), intent(in) :: der
1214 integer, intent(in) :: pml_dir
1215 complex(real64), intent(inout) :: psi(:,:)
1216 integer, intent(in) :: field_dir
1217 complex(real64), intent(inout) :: pml(:)
1218
1219 integer :: ip, ip_in, rs_sign
1220 real(real64) :: pml_c
1221 complex(real64), allocatable :: tmp_partial(:)
1222 complex(real64) :: pml_a, pml_b, pml_g
1223
1225
1226 if (hm%cpml_hamiltonian) then
1227
1228 rs_sign = hm%rs_sign
1229
1230 safe_allocate(tmp_partial(1:der%mesh%np_part))
1231
1232 call zderivatives_partial(der, psi(:,field_dir), tmp_partial(:), pml_dir, set_bc = .false.)
1233 do ip_in = 1, hm%bc%pml%points_number
1234 ip = hm%bc%pml%points_map(ip_in)
1235 pml_c = hm%bc%pml%c(ip_in, pml_dir)
1236 pml_a = hm%bc%pml%a(ip_in, pml_dir)
1237 pml_b = hm%bc%pml%b(ip_in, pml_dir)
1238 pml_g = hm%bc%pml%conv_plus(ip_in, pml_dir, field_dir)
1239 pml(ip) = rs_sign * pml_c * ( tmp_partial(ip) &
1240 + real(pml_a, real64) * real(tmp_partial(ip), real64) &
1241 + m_zi * aimag(pml_a) * aimag(tmp_partial(ip)) &
1242 + real(pml_b, real64) * real(pml_g, real64) &
1243 + m_zi * aimag(pml_b) * aimag(pml_g))
1244 end do
1245
1246 safe_deallocate_a(tmp_partial)
1247 end if
1248
1251
1252
1253 ! ---------------------------------------------------------
1256 subroutine maxwell_pml_calculation_via_riemann_silberstein_medium(hm, der, psi, pml_dir, field_dir, pml)
1257 type(hamiltonian_mxll_t), intent(in) :: hm
1258 type(derivatives_t), intent(in) :: der
1259 integer, intent(in) :: pml_dir
1260 complex(real64), intent(inout) :: psi(:,:)
1261 integer, intent(in) :: field_dir
1262 complex(real64), intent(inout) :: pml(:,:)
1263
1264 integer :: ip, ip_in, np
1265 real(real64) :: pml_c(3)
1266 complex(real64), allocatable :: tmp_partial(:,:)
1267 complex(real64) :: pml_a(3), pml_b(3), pml_g_p(3), pml_g_m(3)
1268
1270
1271 if (hm%cpml_hamiltonian) then
1272
1273 np = der%mesh%np
1274 safe_allocate(tmp_partial(1:np,1:2))
1275
1276 call zderivatives_partial(der, psi(:,field_dir ), tmp_partial(:,1), pml_dir, set_bc = .false.)
1277 call zderivatives_partial(der, psi(:,field_dir+3), tmp_partial(:,2), pml_dir, set_bc = .false.)
1278 do ip_in = 1, hm%bc%pml%points_number
1279 ip = hm%bc%pml%points_map(ip_in)
1280 pml_c(1:3) = hm%bc%pml%c(ip_in, 1:3)
1281 pml_a(1:3) = hm%bc%pml%a(ip_in, 1:3)
1282 pml_b(1:3) = hm%bc%pml%b(ip_in, 1:3)
1283 pml_g_p(1:3) = hm%bc%pml%conv_plus(ip_in, pml_dir, 1:3)
1284 pml_g_m(1:3) = hm%bc%pml%conv_minus(ip_in, pml_dir, 1:3)
1285 pml(ip, 1) = pml_c(pml_dir) * tmp_partial(ip, 1) &
1286 + pml_c(pml_dir) * real(pml_a(pml_dir), real64) * real(tmp_partial(ip, 1), real64) &
1287 + m_zi * pml_c(pml_dir) * aimag(pml_a(pml_dir)) * aimag(tmp_partial(ip, 1)) &
1288 + pml_c(pml_dir) * real(pml_b(pml_dir), real64) * real(pml_g_p(field_dir), real64) &
1289 + m_zi * pml_c(pml_dir) * aimag(pml_b(pml_dir)) * aimag(pml_g_p(field_dir))
1290 pml(ip, 2) = pml_c(pml_dir) * tmp_partial(ip, 2) &
1291 + pml_c(pml_dir) * real(pml_a(pml_dir), real64) * real(tmp_partial(ip, 2), real64) &
1292 + m_zi * pml_c(pml_dir) * aimag(pml_a(pml_dir)) * aimag(tmp_partial(ip, 2)) &
1293 + pml_c(pml_dir) * real(pml_b(pml_dir), real64) * real(pml_g_m(field_dir), real64) &
1294 + m_zi * pml_c(pml_dir) * aimag(pml_b(pml_dir)) * aimag(pml_g_m(field_dir))
1295 end do
1296
1297 end if
1298
1299 safe_deallocate_a(tmp_partial)
1300
1303
1304
1305 ! ---------------------------------------------------------
1307 subroutine maxwell_medium_boundaries_calculation(hm, psi, oppsi)
1308 type(hamiltonian_mxll_t), intent(in) :: hm
1309 complex(real64), intent(in) :: psi(:,:)
1310 complex(real64), intent(inout) :: oppsi(:,:)
1311
1312 integer :: ip, ip_in, idim
1313 real(real64) :: cc, aux_ep(3), aux_mu(3), sigma_e, sigma_m
1314 complex(real64) :: ff_plus(3), ff_minus(3)
1315
1317
1318 do idim = 1, 3
1319 if ((hm%bc%bc_type(idim) == mxll_bc_medium)) then
1320 do ip_in = 1, hm%bc%medium(idim)%points_number
1321 ip = hm%bc%medium(idim)%points_map(ip_in)
1322 cc = hm%bc%medium(idim)%c(ip_in)/p_c
1323 aux_ep(:) = hm%bc%medium(idim)%aux_ep(ip_in, :)
1324 aux_mu(:) = hm%bc%medium(idim)%aux_mu(ip_in, :)
1325 sigma_e = hm%bc%medium(idim)%sigma_e(ip_in)
1326 sigma_m = hm%bc%medium(idim)%sigma_m(ip_in)
1327 ff_plus(1) = psi(ip, 1)
1328 ff_plus(2) = psi(ip, 2)
1329 ff_plus(3) = psi(ip, 3)
1330 ff_minus(1) = psi(ip, 4)
1331 ff_minus(2) = psi(ip, 5)
1332 ff_minus(3) = psi(ip, 6)
1333 aux_ep = dcross_product(aux_ep,real(ff_plus+ff_minus, real64) )
1334 aux_mu = dcross_product(aux_mu,aimag(ff_plus-ff_minus))
1335 oppsi(ip, 1) = oppsi(ip, 1)*cc &
1336 - cc * aux_ep(1) - cc * m_zi * aux_mu(1) &
1337 - m_zi * sigma_e * real(ff_plus(1) + ff_minus(1), real64) &
1338 - m_zi * sigma_m * m_zi * aimag(ff_plus(1) - ff_minus(1))
1339 oppsi(ip, 4) = oppsi(ip, 4)*cc &
1340 + cc * aux_ep(1) - cc * m_zi * aux_mu(1) &
1341 - m_zi * sigma_e * real(ff_plus(1) + ff_minus(1), real64) &
1342 + m_zi * sigma_m * m_zi * aimag(ff_plus(1) - ff_minus(1))
1343 oppsi(ip, 2) = oppsi(ip, 2)*cc &
1344 - cc * aux_ep(2) - cc * m_zi * aux_mu(2) &
1345 - m_zi * sigma_e * real(ff_plus(2) + ff_minus(2), real64) &
1346 - m_zi * sigma_m * m_zi * aimag(ff_plus(2) - ff_minus(2))
1347 oppsi(ip, 5) = oppsi(ip, 5)*cc &
1348 + cc * aux_ep(2) - cc * m_zi * aux_mu(2) &
1349 - m_zi * sigma_e * real(ff_plus(2) + ff_minus(2), real64) &
1350 + m_zi * sigma_m * m_zi * aimag(ff_plus(2) - ff_minus(2))
1351 oppsi(ip, 3) = oppsi(ip, 3)*cc &
1352 - cc * aux_ep(3) - cc * m_zi * aux_mu(3) &
1353 - m_zi * sigma_e * real(ff_plus(3) + ff_minus(3), real64) &
1354 - m_zi * sigma_m * m_zi * aimag(ff_plus(3) - ff_minus(3))
1355 oppsi(ip, 6) = oppsi(ip, 6)*cc &
1356 + cc * aux_ep(3) - cc * m_zi * aux_mu(3) &
1357 - m_zi * sigma_e * real(ff_plus(3) + ff_minus(3), real64) &
1358 + m_zi * sigma_m * m_zi * aimag(ff_plus(3) - ff_minus(3))
1359 end do
1360 end if
1361 end do
1362
1365
1366
1367 ! ---------------------------------------------------------
1368 ! > Maxwell Hamiltonian including medium boxes
1369 subroutine maxwell_medium_boxes_calculation(hm, der, psi, oppsi)
1370 type(hamiltonian_mxll_t), intent(in) :: hm
1371 type(derivatives_t), intent(in) :: der
1372 complex(real64), intent(in) :: psi(:,:)
1373 complex(real64), intent(inout) :: oppsi(:,:)
1374
1375 integer :: ip, il
1376 real(real64) :: cc, aux_ep(3), aux_mu(3), sigma_e, sigma_m
1377 complex(real64) :: ff_plus(3), ff_minus(3)
1378
1380
1381 if (hm%calc_medium_box) then
1382 do il = 1, size(hm%medium_boxes)
1383 assert(.not. hm%medium_boxes(il)%has_mapping)
1384 do ip = 1, hm%medium_boxes(il)%points_number
1385 cc = hm%medium_boxes(il)%c(ip)/p_c
1386 if (abs(cc) <= m_epsilon) cycle
1387 aux_ep(1:3) = hm%medium_boxes(il)%aux_ep(ip, 1:3)
1388 aux_mu(1:3) = hm%medium_boxes(il)%aux_mu(ip, 1:3)
1389 sigma_e = hm%medium_boxes(il)%sigma_e(ip)
1390 sigma_m = hm%medium_boxes(il)%sigma_m(ip)
1391 ff_plus(1) = psi(ip, 1)
1392 ff_plus(2) = psi(ip, 2)
1393 ff_plus(3) = psi(ip, 3)
1394 ff_minus(1) = psi(ip, 4)
1395 ff_minus(2) = psi(ip, 5)
1396 ff_minus(3) = psi(ip, 6)
1397 aux_ep = dcross_product(aux_ep, real(ff_plus+ff_minus, real64) )
1398 aux_mu = dcross_product(aux_mu, aimag(ff_plus-ff_minus))
1399 oppsi(ip, 1) = oppsi(ip,1)*cc &
1400 - cc * aux_ep(1) - cc * m_zi * aux_mu(1) &
1401 - m_zi * sigma_e * real(ff_plus(1) + ff_minus(1), real64) &
1402 - m_zi * sigma_m * m_zi * aimag(ff_plus(1) - ff_minus(1))
1403 oppsi(ip, 4) = oppsi(ip,4)*cc &
1404 + cc * aux_ep(1) - cc * m_zi * aux_mu(1) &
1405 - m_zi * sigma_e * real(ff_plus(1) + ff_minus(1), real64) &
1406 + m_zi * sigma_m * m_zi * aimag(ff_plus(1) - ff_minus(1))
1407 oppsi(ip, 2) = oppsi(ip,2)*cc &
1408 - cc * aux_ep(2) - cc * m_zi * aux_mu(2) &
1409 - m_zi * sigma_e * real(ff_plus(2) + ff_minus(2), real64) &
1410 - m_zi * sigma_m * m_zi * aimag(ff_plus(2) - ff_minus(2))
1411 oppsi(ip, 5) = oppsi(ip,5)*cc &
1412 + cc * aux_ep(2) - cc * m_zi * aux_mu(2) &
1413 - m_zi * sigma_e * real(ff_plus(2) + ff_minus(2), real64) &
1414 + m_zi * sigma_m * m_zi * aimag(ff_plus(2) - ff_minus(2))
1415 oppsi(ip, 3) = oppsi(ip,3)*cc &
1416 - cc * aux_ep(3) - cc * m_zi * aux_mu(3) &
1417 - m_zi * sigma_e * real(ff_plus(3) + ff_minus(3), real64) &
1418 - m_zi * sigma_m * m_zi * aimag(ff_plus(3) - ff_minus(3))
1419 oppsi(ip, 6) = oppsi(ip,6)*cc &
1420 + cc * aux_ep(3) - cc * m_zi * aux_mu(3) &
1421 - m_zi * sigma_e * real(ff_plus(3) + ff_minus(3), real64) &
1422 + m_zi * sigma_m * m_zi * aimag(ff_plus(3) - ff_minus(3))
1423 end do
1424 end do
1425 end if
1426
1429
1430end module hamiltonian_mxll_oct_m
1431
1432!! Local Variables:
1433!! mode: f90
1434!! coding: utf-8
1435!! End:
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.
Definition: batch.F90:135
This module implements common operations on batches of mesh functions.
Definition: batch_ops.F90:118
Module implementing boundary conditions in Octopus.
Definition: boundaries.F90:124
This module calculates the derivatives (gradients, Laplacians, etc.) of a function.
integer, parameter, public der_star
real(real64), parameter, public m_zero
Definition: global.F90:191
real(real64), parameter, public m_pi
some mathematical constants
Definition: global.F90:189
real(real64), parameter, public p_c
Electron gyromagnetic ratio, see Phys. Rev. Lett. 130, 071801 (2023)
Definition: global.F90:233
real(real64), parameter, public m_one
Definition: global.F90:192
real(real64), parameter, public m_three
Definition: global.F90:194
This module implements the underlying real-space grid.
Definition: grid.F90:119
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.
Definition: math.F90:117
subroutine, public bc_mxll_end(bc)
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:120
subroutine, public messages_experimental(name, namespace)
Definition: messages.F90:1063
This module defines non-local operators.
subroutine, public profiling_out(label)
Increment out counter and sum up difference between entry and exit time.
Definition: profiling.F90:631
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
Definition: profiling.F90:554
This module handles spin dimensions of the states and the k-point distribution.
The abstract Hamiltonian class defines a skeleton for specific implementations.
int true(void)