Octopus
maxwell.F90
Go to the documentation of this file.
1!! Copyright (C) 2019-2020 Franco Bonafe, Heiko Appel, Rene Jestaedt
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
21module maxwell_oct_m
22 use accel_oct_m
23 use batch_oct_m
28 use debug_oct_m
35 use global_oct_m
36 use grid_oct_m
39 use index_oct_m
43 use io_oct_m
48 use loct_oct_m
50 use math_oct_m
52 use mesh_oct_m
55 use mpi_oct_m
61 use output_oct_m
64 use parser_oct_m
76 use sort_oct_m
77 use space_oct_m
78 use system_oct_m
83 use types_oct_m
84 use unit_oct_m
86
87
88 implicit none
89
90 private
91 public :: &
92 maxwell_t, &
94
95 integer, parameter, public :: &
96 MULTIGRID_MX_TO_MA_EQUAL = 1, &
98
101 type, extends(system_t) :: maxwell_t
102 type(space_t) :: space
103 type(states_mxll_t) :: st
104 type(hamiltonian_mxll_t) :: hm
105 type(grid_t) :: gr
106 type(output_t) :: outp
107 type(multicomm_t) :: mc
108
109 type(mesh_interpolation_t) :: mesh_interpolate
110
111 type(propagator_mxll_t) :: tr_mxll
112 type(td_write_t) :: write_handler
113 type(c_ptr) :: output_handle
115 complex(real64), allocatable :: ff_rs_inhom_t1(:,:), ff_rs_inhom_t2(:,:)
116 complex(real64), allocatable :: rs_state_init(:,:)
117 type(time_interpolation_t), pointer :: current_interpolation
118 real(real64) :: bc_bounds(2, 3), dt_bounds(2, 3)
119 integer :: energy_update_iter
120 type(restart_t) :: restart_dump
121 type(restart_t) :: restart
122
123 type(lattice_vectors_t) :: latt
124
125 type(helmholtz_decomposition_t) :: helmholtz
126
127 logical :: write_previous_state = .false.
128
129 contains
130 procedure :: init_interaction => maxwell_init_interaction
131 procedure :: init_parallelization => maxwell_init_parallelization
132 procedure :: initialize => maxwell_initialize
133 procedure :: do_algorithmic_operation => maxwell_do_algorithmic_operation
134 procedure :: is_tolerance_reached => maxwell_is_tolerance_reached
135 procedure :: init_interaction_as_partner => maxwell_init_interaction_as_partner
136 procedure :: copy_quantities_to_interaction => maxwell_copy_quantities_to_interaction
137 procedure :: output_start => maxwell_output_start
138 procedure :: output_write => maxwell_output_write
139 procedure :: output_finish => maxwell_output_finish
140 procedure :: update_interactions_start => maxwell_update_interactions_start
141 procedure :: update_interactions_finish => maxwell_update_interactions_finish
142 procedure :: restart_write_data => maxwell_restart_write_data
143 procedure :: restart_read_data => maxwell_restart_read_data
144 procedure :: update_kinetic_energy => maxwell_update_kinetic_energy
145 procedure :: get_current => maxwell_get_current
146 final :: maxwell_finalize
147 end type maxwell_t
148
149 interface maxwell_t
150 procedure maxwell_constructor
151 end interface maxwell_t
152
153contains
154
155 ! ---------------------------------------------------------
156 function maxwell_constructor(namespace) result(sys)
157 class(maxwell_t), pointer :: sys
158 type(namespace_t), intent(in) :: namespace
159
160 push_sub(maxwell_constructor)
161
162 allocate(sys)
163
164 call maxwell_init(sys, namespace)
165
166 pop_sub(maxwell_constructor)
167 end function maxwell_constructor
168
169
170 ! ---------------------------------------------------------
171 subroutine maxwell_init(this, namespace)
172 class(maxwell_t), intent(inout) :: this
173 type(namespace_t), intent(in) :: namespace
174
175
176 push_sub(maxwell_init)
177
178 call profiling_in("MAXWELL_INIT")
179
180 this%namespace = namespace
181
182 call messages_obsolete_variable(this%namespace, 'SystemName')
183 call messages_print_with_emphasis(msg='Maxwell Simulation Start', namespace=this%namespace)
184 this%space = space_t(this%namespace)
185 if (this%space%is_periodic()) then
186 call messages_not_implemented('Maxwell for periodic systems', namespace=namespace)
187 end if
189 call grid_init_stage_1(this%gr, this%namespace, this%space)
190 call states_mxll_init(this%st, this%namespace, this%space)
191 this%latt = lattice_vectors_t(this%namespace, this%space)
192
193 call this%quantities%add(quantity_t("E field", updated_on_demand = .false.))
194 call this%quantities%add(quantity_t("B field", updated_on_demand = .false.))
195 call this%quantities%add(quantity_t("vector potential", updated_on_demand = .false.))
197 call mesh_interpolation_init(this%mesh_interpolate, this%gr)
199 this%supported_interactions = [linear_medium_to_em_field, current_to_mxll_field]
201
202 call profiling_out("MAXWELL_INIT")
203
205 end subroutine maxwell_init
207 ! ---------------------------------------------------------
208 subroutine maxwell_init_interaction(this, interaction)
209 class(maxwell_t), target, intent(inout) :: this
210 class(interaction_t), intent(inout) :: interaction
212 integer :: depth
215
216 select type (interaction)
218 call interaction%init(this%gr)
220 call interaction%init(this%gr, this%st%dim)
221 call interaction%init_space_latt(this%space, this%latt)
222
223 ! set interpolation depth for interaction
224 ! interpolation depth depends on the propagator
225 select type (prop => this%algo)
227 depth = 2
229 depth = 2
231 depth = 2
233 depth = 4
234 class default
235 message(1) = "The chosen propagator does not yet support interaction interpolation"
236 call messages_fatal(1, namespace=this%namespace)
237 end select
238 call interaction%init_interpolation(depth, interaction%label, cmplx=.true.)
240 this%hm%current_density_from_medium = .true.
241 class default
242 message(1) = "Trying to initialize an unsupported interaction by Maxwell."
243 call messages_fatal(1, namespace=this%namespace)
244 end select
245
247 end subroutine maxwell_init_interaction
248
249 ! ---------------------------------------------------------
250 subroutine maxwell_init_parallelization(this, grp)
251 class(maxwell_t), intent(inout) :: this
252 type(mpi_grp_t), intent(in) :: grp
253
254 integer(int64) :: index_range(4)
255 integer :: ierr, ip, pos_index, rankmin
256 real(real64) :: dmin
257
259
260 call system_init_parallelization(this, grp)
261
262 ! store the ranges for these two indices (serves as initial guess
263 ! for parallelization strategy)
264 index_range(1) = this%gr%np_global ! Number of points in mesh
265 index_range(2) = this%st%nst ! Number of states
266 index_range(3) = 1 ! Number of k-points
267 index_range(4) = 100000 ! Some large number
268
269 ! create index and domain communicators
270 call multicomm_init(this%mc, this%namespace, mpi_world, calc_mode_par, mpi_world%size, &
271 index_range, (/ 5000, 1, 1, 1 /))
272
273 call grid_init_stage_2(this%gr, this%namespace, this%space, this%mc)
274 call output_mxll_init(this%outp, this%namespace, this%space)
275 call hamiltonian_mxll_init(this%hm, this%namespace, this%gr, this%st)
276
277 this%st%energy_rate = m_zero
278 this%st%delta_energy = m_zero
279 this%st%energy_via_flux_calc = m_zero
280 this%st%trans_energy_rate = m_zero
281 this%st%trans_delta_energy = m_zero
282 this%st%trans_energy_via_flux_calc = m_zero
283 this%st%plane_waves_energy_rate = m_zero
284 this%st%plane_waves_delta_energy = m_zero
285 this%st%plane_waves_energy_via_flux_calc = m_zero
286
287 safe_allocate(this%rs_state_init(1:this%gr%np_part, 1:this%st%dim))
288 this%rs_state_init(:,:) = m_z0
289
290 this%energy_update_iter = 1
291
292 call poisson_init(this%st%poisson, this%namespace, this%space, this%gr%der, this%mc, this%gr%stencil)
293
294 call propagator_mxll_init(this%gr, this%namespace, this%st, this%hm, this%tr_mxll)
295 call states_mxll_allocate(this%st, this%gr)
296 call external_current_init(this%st, this%space, this%namespace, this%gr)
297 this%hm%propagation_apply = .true.
298
299 if (parse_is_defined(this%namespace, 'MaxwellIncidentWaves') .and. (this%tr_mxll%bc_plane_waves)) then
300 this%st%rs_state_plane_waves(:,:) = m_z0
301 end if
302
303 ! set map for selected points
304 do ip = 1, this%st%selected_points_number
305 pos_index = mesh_nearest_point(this%gr, this%st%selected_points_coordinate(:,ip), dmin, rankmin)
306 if (this%gr%mpi_grp%rank == rankmin) then
307 this%st%selected_points_map(ip) = pos_index
308 else
309 this%st%selected_points_map(ip) = -1
310 end if
311 end do
312 if (accel_is_enabled()) then
313 call accel_create_buffer(this%st%buff_selected_points_map, accel_mem_read_only, type_integer, &
314 this%st%selected_points_number)
315 call accel_write_buffer(this%st%buff_selected_points_map, this%st%selected_points_number, &
316 this%st%selected_points_map)
317 end if
318
319 this%hm%plane_waves_apply = .true.
320 this%hm%spatial_constant_apply = .true.
321
322 call bc_mxll_init(this%hm%bc, this%namespace, this%space, this%gr, this%st)
323 this%bc_bounds(:,1:3) = this%hm%bc%bc_bounds(:,1:3)
324 call inner_and_outer_points_mapping(this%gr, this%st, this%bc_bounds)
325 this%dt_bounds(2, 1:3) = this%bc_bounds(1, 1:3)
326 this%dt_bounds(1, 1:3) = this%bc_bounds(1, 1:3) - this%gr%der%order * this%gr%spacing(1:3)
327 call surface_grid_points_mapping(this%gr, this%st, this%dt_bounds)
328
329 call restart_init(this%restart, this%namespace, restart_td, restart_type_load, this%mc, ierr, mesh=this%gr)
330 call restart_init(this%restart_dump, this%namespace, restart_td, restart_type_dump, this%mc, ierr, mesh=this%gr)
331
332 ! initialize batches
333 call zbatch_init(this%st%rs_stateb, 1, 1, this%st%dim, this%gr%np_part)
334 if (this%st%pack_states) call this%st%rs_stateb%do_pack()
335 call this%st%rs_stateb%copy_to(this%st%rs_state_prevb)
336 call this%st%rs_stateb%copy_to(this%st%inhomogeneousb)
337 if (this%tr_mxll%bc_plane_waves) then
338 call this%st%rs_stateb%copy_to(this%st%rs_state_plane_wavesb)
339 end if
340
341 ! the order should depend on the propagator
342 !this%current_interpolation => time_interpolation_t(this%gr%np, this%st%dim, 2, .true., "current")
343 ! Initialize Helmholtz decomposition
344 call this%helmholtz%init(this%namespace, this%gr, this%mc, this%space)
345
347 end subroutine maxwell_init_parallelization
348
349 ! ---------------------------------------------------------
350 subroutine maxwell_initialize(this)
351 class(maxwell_t), intent(inout) :: this
352
353 real(real64) :: courant
354
355
356 push_sub(maxwell_initialize)
357
358 call profiling_in("MAXWELL_INIT_CONDITIONS")
359
360 select type (algo => this%algo)
361 class is (propagator_t)
362
363 ! The reference Courant–Friedrichs–Lewy condition for stability here
364 ! is calculated with respect to R. Courant, K. Friedrichs, and H. Lewy
365 ! On the Partial Difference Equations of
366 ! Mathematical Physics, IBM Journal of Research and Development, vol. 11
367 ! pp. 215-234, Mar. 1967.
368 courant = m_one/(p_c * sqrt(m_one/this%gr%spacing(1)**2 + m_one/this%gr%spacing(2)**2 + &
369 m_one/this%gr%spacing(3)**2))
370
371 if (algo%dt > (m_one + m_fourth) * courant) then
372 write(message(1),'(a)') 'The timestep is too large compared to the Courant-Friedrichs-Lewy'
373 write(message(2),'(a)') 'stability condition. Time propagation might not be stable.'
374 call messages_warning(2, namespace=this%namespace)
375 end if
376
377 if (parse_is_defined(this%namespace, 'UserDefinedInitialMaxwellStates')) then
378 call states_mxll_read_user_def(this%namespace, this%space, this%gr, this%gr%der, this%st, this%hm%bc,&
379 this%rs_state_init)
380 call messages_print_with_emphasis(msg="Setting initial EM field inside box", namespace=this%namespace)
381 ! TODO: add consistency check that initial state fulfills Gauss laws
382 this%st%rs_state(:,:) = this%st%rs_state + this%rs_state_init
383 if (this%tr_mxll%bc_plane_waves) then
384 this%st%rs_state_plane_waves(:,:) = this%rs_state_init
385 end if
386 end if
387
388 ! initialize the spatial constant field according to the conditions set in the
389 ! UserDefinedConstantSpatialMaxwellField block
390 if (this%tr_mxll%bc_constant) then
391 call spatial_constant_calculation(this%tr_mxll%bc_constant, this%st, this%gr, this%hm, m_zero, &
392 algo%dt, m_zero, this%st%rs_state, set_initial_state = .true.)
393 ! for mesh parallelization, this needs communication!
394 this%st%rs_state_const(:) = this%st%rs_state(mesh_global_index_from_coords(this%gr, [0,0,0]),:)
395 end if
396
397 if (parse_is_defined(this%namespace, 'UserDefinedInitialMaxwellStates')) then
398 safe_deallocate_a(this%rs_state_init)
399 end if
400
401 call hamiltonian_mxll_update(this%hm, time = m_zero)
402
403 ! calculate Maxwell energy
404 call energy_mxll_calc(this%gr, this%st, this%hm, this%hm%energy, this%st%rs_state, &
405 this%st%rs_state_plane_waves)
406
407 ! Get RS states values for selected points
408 call get_rs_state_at_point(this%st%selected_points_rs_state(:,:), this%st%rs_state, &
409 this%st%selected_points_coordinate(:,:), this%st, this%gr)
410
411 call mxll_set_batch(this%st%rs_stateb, this%st%rs_state, this%gr%np, this%st%dim)
412 call batch_set_zero(this%st%inhomogeneousb)
413 if (this%tr_mxll%bc_plane_waves) then
414 call mxll_set_batch(this%st%rs_state_plane_wavesb, this%st%rs_state_plane_waves, this%gr%np, this%st%dim)
415 end if
416 end select
417
418 call profiling_out("MAXWELL_INIT_CONDITIONS")
419
420 pop_sub(maxwell_initialize)
421 end subroutine maxwell_initialize
422
423 ! ---------------------------------------------------------
424 logical function maxwell_do_algorithmic_operation(this, operation, updated_quantities) result(done)
425 class(maxwell_t), intent(inout) :: this
426 class(algorithmic_operation_t), intent(in) :: operation
427 character(len=:), allocatable, intent(out) :: updated_quantities(:)
428
429 complex(real64), allocatable :: charge_density_ext(:)
430 real(real64) :: current_time
431
433 call profiling_in(trim(this%namespace%get())//":"//trim(operation%id))
434
435 current_time = this%iteration%value()
436
437 select type (algo => this%algo)
438 class is (propagator_t)
439
440 done = .true.
441 select case (operation%id)
443 ! For the moment we do nothing
444
445 case (expmid_start)
446 safe_allocate(this%ff_rs_inhom_t1(1:this%gr%np_part, 1:this%hm%dim))
447 safe_allocate(this%ff_rs_inhom_t2(1:this%gr%np_part, 1:this%hm%dim))
448
449 case (expmid_finish)
450 safe_deallocate_a(this%ff_rs_inhom_t1)
451 safe_deallocate_a(this%ff_rs_inhom_t2)
452
453 case (expmid_extrapolate)
454 if (this%hm%current_density_ext_flag .or. this%hm%current_density_from_medium) then
455 call this%get_current(current_time, this%st%rs_current_density_t1)
456 call this%get_current(current_time+algo%dt, this%st%rs_current_density_t2)
457 end if
458
459 case (expmid_propagate)
460 ! Propagation
461
462 !We first compute three external charge and current densities and we convert them as RS vectors
463 safe_allocate(charge_density_ext(1:this%gr%np))
464
465 !No charge density at the moment
466 charge_density_ext = m_z0
467
468 call transform_rs_densities(this%hm, this%gr, charge_density_ext, &
469 this%st%rs_current_density_t1, this%ff_rs_inhom_t1, rs_trans_forward)
470 call transform_rs_densities(this%hm, this%gr, charge_density_ext, &
471 this%st%rs_current_density_t2, this%ff_rs_inhom_t2, rs_trans_forward)
472
473 safe_deallocate_a(charge_density_ext)
474
475 ! Propagation dt with H_maxwell
476 call mxll_propagation_step(this%hm, this%namespace, this%gr, this%space, this%st, this%tr_mxll,&
477 this%st%rs_stateb, this%ff_rs_inhom_t1, this%ff_rs_inhom_t2, current_time, algo%dt)
478
479 updated_quantities = [character(16) :: "E field", "B field", "vector potential"]
480
481 case (leapfrog_start)
482 if (any(this%hm%bc%bc_ab_type == option__maxwellabsorbingboundaries__cpml)) then
483 call bc_mxll_initialize_pml_simple(this%hm%bc, this%space, this%gr, this%hm%c_factor, algo%dt)
484 end if
485
486 case (leapfrog_finish)
487
488 case (leapfrog_propagate)
489 if (this%hm%current_density_ext_flag .or. this%hm%current_density_from_medium) then
490 call this%get_current(current_time, this%st%rs_current_density_t1)
491 call mxll_set_batch(this%st%inhomogeneousb, this%st%rs_current_density_t1, this%gr%np, this%st%dim)
492 call batch_scal(this%gr%np, -m_one, this%st%inhomogeneousb)
493 end if
494
495 call mxll_propagate_leapfrog(this%hm, this%namespace, this%gr, this%space, this%st, this%tr_mxll, &
496 current_time, algo%dt, this%iteration%counter())
497
498 updated_quantities = [character(16) :: "E field", "B field", "vector potential"]
499
500 case (exp_gauss1_start)
501 safe_allocate(this%ff_rs_inhom_t1(1:this%gr%np_part, 1:this%hm%dim))
502 if (any(this%hm%bc%bc_ab_type == option__maxwellabsorbingboundaries__cpml)) then
503 call bc_mxll_initialize_pml_simple(this%hm%bc, this%space, this%gr, this%hm%c_factor, algo%dt)
504 end if
505
506 case (exp_gauss1_finish)
507 safe_deallocate_a(this%ff_rs_inhom_t1)
508
510 if (this%hm%current_density_ext_flag .or. this%hm%current_density_from_medium) then
511 call this%get_current(current_time+algo%dt*m_half, this%st%rs_current_density_t1)
512 end if
513
515 ! the propagation step is
516 ! F_{n+1} = F_n + dt * phi_1(-i*dt*H) [-i*H F_n - J_1]
517 ! where J_1 = J(t_n + dt/2)
518 call mxll_propagate_expgauss1(this%hm, this%namespace, this%gr, this%space, this%st, this%tr_mxll, &
519 current_time, algo%dt)
520
521 updated_quantities = [character(16) :: "E field", "B field", "vector potential"]
522
523 case (exp_gauss2_start)
524 safe_allocate(this%ff_rs_inhom_t1(1:this%gr%np_part, 1:this%hm%dim))
525 safe_allocate(this%ff_rs_inhom_t2(1:this%gr%np_part, 1:this%hm%dim))
526 if (any(this%hm%bc%bc_ab_type == option__maxwellabsorbingboundaries__cpml)) then
527 call bc_mxll_initialize_pml_simple(this%hm%bc, this%space, this%gr, this%hm%c_factor, algo%dt)
528 end if
529
530 case (exp_gauss2_finish)
531 safe_deallocate_a(this%ff_rs_inhom_t1)
532 safe_deallocate_a(this%ff_rs_inhom_t2)
533
535 if (this%hm%current_density_ext_flag .or. this%hm%current_density_from_medium) then
536 call this%get_current(current_time+algo%dt*(m_half-sqrt(m_three)/6.0_real64), this%st%rs_current_density_t1)
537 call this%get_current(current_time+algo%dt*(m_half+sqrt(m_three)/6.0_real64), this%st%rs_current_density_t2)
538 end if
539
541 ! the propagation step is
542 ! F_{n+1} = F_n + dt * phi_1(-i*dt*H) [-i*H F_n - a_1 J_1 - a_2 J_2]
543 ! + dt * phi_2(-i*dt*H) [-b_1 J_1 - b_2 J_2]
544 ! where J_1 = J(t_n + (1/2 - sqrt(3)/6)*dt),
545 ! J_2 = J(t_n + (1/2 + sqrt(3)/6)*dt),
546 ! a_1 = 1/2*(1+sqrt(3)),
547 ! a_2 = 1/2*(1-sqrt(3)),
548 ! b_1 = -sqrt(3),
549 ! b_2 = sqrt(3)
550 call mxll_propagate_expgauss2(this%hm, this%namespace, this%gr, this%space, this%st, this%tr_mxll, &
551 current_time, algo%dt)
552
553 updated_quantities = [character(16) :: "E field", "B field", "vector potential"]
554
555 case (iteration_done)
557 done = .false.
558
559 case default
560 done = .false.
561 end select
562
563 end select
564
565 call profiling_out(trim(this%namespace%get())//":"//trim(operation%id))
568
569 ! ---------------------------------------------------------
570 logical function maxwell_is_tolerance_reached(this, tol) result(converged)
571 class(maxwell_t), intent(in) :: this
572 real(real64), intent(in) :: tol
573
575
576 converged = .false.
577
580
581 ! ---------------------------------------------------------
582 subroutine maxwell_init_interaction_as_partner(partner, interaction)
583 class(maxwell_t), intent(in) :: partner
584 class(interaction_surrogate_t), intent(inout) :: interaction
585
587
588 select type (interaction)
589 type is (lorentz_force_t)
590 ! Nothing to be initialized for the Lorentz force.
592 call interaction%init_from_partner(partner%gr, partner%space, partner%namespace)
594 call interaction%init_from_partner(partner%gr, partner%space, partner%namespace)
596 call interaction%init_from_partner(partner%gr, partner%space, partner%namespace)
597 class default
598 message(1) = "Unsupported interaction."
599 call messages_fatal(1, namespace=partner%namespace)
600 end select
601
604
605 ! ---------------------------------------------------------
606 subroutine maxwell_copy_quantities_to_interaction(partner, interaction)
607 class(maxwell_t), intent(inout) :: partner
608 class(interaction_surrogate_t), intent(inout) :: interaction
609
610 integer :: ip
611 complex(real64) :: interpolated_value(3)
612 real(real64), allocatable :: b_field(:,:), vec_pot(:,:)
613
615 call profiling_in(trim(partner%namespace%get())//":"//"COPY_QUANTITY_INTER")
616
617 select type (interaction)
618 type is (lorentz_force_t)
619 call mxll_get_batch(partner%st%rs_stateb, partner%st%rs_state, partner%gr%np, partner%st%dim)
620
621 do ip = 1, interaction%system_np
622 call mesh_interpolation_evaluate(partner%mesh_interpolate, partner%st%rs_state(:,1), &
623 interaction%system_pos(:, ip), interpolated_value(1))
624 call mesh_interpolation_evaluate(partner%mesh_interpolate, partner%st%rs_state(:,2), &
625 interaction%system_pos(:, ip), interpolated_value(2))
626 call mesh_interpolation_evaluate(partner%mesh_interpolate, partner%st%rs_state(:,3), &
627 interaction%system_pos(:, ip), interpolated_value(3))
628 call get_electric_field_vector(interpolated_value, interaction%partner_e_field(:, ip))
629 call get_magnetic_field_vector(interpolated_value, 1, interaction%partner_b_field(:, ip))
630 end do
631
633 call mxll_get_batch(partner%st%rs_stateb, partner%st%rs_state, partner%gr%np, partner%st%dim)
634 select case (interaction%type)
635
636 case (mxll_field_total)
637 call get_electric_field_state(partner%st%rs_state, partner%gr, interaction%partner_field)
638
639 case (mxll_field_trans)
640 call get_transverse_rs_state(partner%helmholtz, partner%st, partner%namespace)
641 call get_electric_field_state(partner%st%rs_state_trans, partner%gr, interaction%partner_field)
642
643 case (mxll_field_long)
644 call partner%helmholtz%get_long_field(partner%namespace, partner%st%rs_state_long, total_field=partner%st%rs_state)
645 call get_electric_field_state(partner%st%rs_state_long, partner%gr, interaction%partner_field)
646
647 case default
648 message(1) = "Unknown type of field requested by interaction."
649 call messages_fatal(1, namespace=partner%namespace)
650 end select
651 call interaction%do_mapping()
652
654 call mxll_get_batch(partner%st%rs_stateb, partner%st%rs_state, partner%gr%np, partner%st%dim)
655 safe_allocate(b_field(1:partner%gr%np_part, 1:partner%gr%box%dim))
656 safe_allocate(vec_pot(1:partner%gr%np_part, 1:partner%gr%box%dim))
657 ! magnetic field is always transverse
658 call get_magnetic_field_state(partner%st%rs_state, partner%gr, partner%st%rs_sign, b_field, &
659 partner%st%mu(1:partner%gr%np), partner%gr%np)
660 ! vector potential stored in partner_field
661 call partner%helmholtz%get_vector_potential(partner%namespace, vec_pot, trans_field=b_field)
662 ! in the convention used by the electronic system, the -1/c factor is included in the vector potential
663 ! but we do not divide by it, because the Maxwell code is in SI units, so the vector potential units
664 ! are suitable for a (p + q.A) minimal coupling interaction
665 interaction%partner_field(1:partner%gr%np,1:partner%gr%box%dim) = &
666 vec_pot(1:partner%gr%np,1:partner%gr%box%dim)
667 safe_deallocate_a(b_field)
668 safe_deallocate_a(vec_pot)
669 call interaction%do_mapping()
670
672 call mxll_get_batch(partner%st%rs_stateb, partner%st%rs_state, partner%gr%np, partner%st%dim)
673 call get_magnetic_field_state(partner%st%rs_state, partner%gr, partner%st%rs_sign, &
674 interaction%partner_field, partner%st%mu(1:partner%gr%np), partner%gr%np)
675 call interaction%do_mapping()
676
677 class default
678 message(1) = "Unsupported interaction."
679 call messages_fatal(1, namespace=partner%namespace)
680 end select
681
682 call profiling_out(trim(partner%namespace%get())//":"//"COPY_QUANTITY_INTER")
685
686 ! ---------------------------------------------------------
687 subroutine maxwell_output_start(this)
688 class(maxwell_t), intent(inout) :: this
689
690 real(real64) :: itr_value
691
692 push_sub(maxwell_output_start)
693 call profiling_in(trim(this%namespace%get())//":"//"OUTPUT_START")
694
695 select type (algo => this%algo)
696 class is (propagator_t)
697 call mxll_get_batch(this%st%rs_stateb, this%st%rs_state, this%gr%np, this%st%dim)
698
699 call td_write_mxll_init(this%write_handler, this%namespace, this%iteration%counter(), algo%dt)
700 if (this%st%fromScratch) then
701 call td_write_mxll_iter(this%write_handler, this%space, this%gr, this%st, this%hm, this%helmholtz, algo%dt, &
702 this%iteration%counter(), this%namespace)
703 itr_value = this%iteration%value()
704 call td_write_mxll_free_data(this%write_handler, this%namespace, this%space, this%gr, this%st, this%hm, this%helmholtz, &
705 this%outp, this%iteration%counter(), itr_value)
706 end if
707
708 end select
709
710 call profiling_out(trim(this%namespace%get())//":"//"OUTPUT_START")
711 pop_sub(maxwell_output_start)
712 end subroutine maxwell_output_start
713
714 ! ---------------------------------------------------------
715 subroutine maxwell_output_write(this)
716 class(maxwell_t), intent(inout) :: this
717
718 logical :: stopping, reached_output_interval
719 real(real64) :: itr_value
720 integer :: iout
721
722 push_sub(maxwell_output_write)
723 call profiling_in(trim(this%namespace%get())//":"//"OUTPUT_WRITE")
724
725 select type (algo => this%algo)
726 class is (propagator_t)
727 stopping = clean_stop(this%mc%master_comm)
728
729 call td_write_mxll_iter(this%write_handler, this%space, this%gr, this%st, this%hm, this%helmholtz, algo%dt, &
730 this%iteration%counter(), this%namespace)
731
732 reached_output_interval = .false.
733 do iout = 1, out_maxwell_max
734 if (this%outp%output_interval(iout) > 0) then
735 if (mod(this%iteration%counter(), this%outp%output_interval(iout)) == 0) then
736 reached_output_interval = .true.
737 exit
738 end if
739 end if
740 end do
741
742 if (reached_output_interval .or. stopping) then
743 call mxll_get_batch(this%st%rs_stateb, this%st%rs_state, this%gr%np, this%st%dim)
744
745 itr_value = this%iteration%value()
746 call td_write_mxll_free_data(this%write_handler, this%namespace, this%space, this%gr, this%st, this%hm, this%helmholtz, &
747 this%outp, this%iteration%counter(), itr_value)
748 end if
749
750 end select
751
752 call profiling_out(trim(this%namespace%get())//":"//"OUTPUT_WRITE")
753 pop_sub(maxwell_output_write)
754 end subroutine maxwell_output_write
755
756 ! ---------------------------------------------------------
757 subroutine maxwell_output_finish(this)
758 class(maxwell_t), intent(inout) :: this
759
760
761 push_sub(maxwell_output_finish)
762
763 call profiling_in("MAXWELL_OUTPUT_FINISH")
764
765 select type (algo => this%algo)
766 class is (propagator_t)
767 call td_write_mxll_end(this%write_handler)
768 end select
769
770 call profiling_out("MAXWELL_OUTPUT_FINISH")
771
772 pop_sub(maxwell_output_finish)
773 end subroutine maxwell_output_finish
774
775 ! ---------------------------------------------------------
776 subroutine maxwell_update_interactions_start(this)
777 class(maxwell_t), intent(inout) :: this
778
779 type(interaction_iterator_t) :: iter
780 integer :: int_counter
781
783
784 int_counter = 0
785 call iter%start(this%interactions)
786 do while (iter%has_next())
787 select type (interaction => iter%get_next())
789 int_counter = int_counter + 1
790 end select
791 end do
792
793 if (int_counter /= 0 .and. .not. allocated(this%hm%medium_boxes)) then
794 safe_allocate(this%hm%medium_boxes(1:int_counter))
795 this%hm%calc_medium_box = .true.
796 end if
797
800
801 ! ---------------------------------------------------------
803 class(maxwell_t), intent(inout) :: this
804
805 type(interaction_iterator_t) :: iter
806 integer :: iint
807
809
810 iint = 0
811
812 call iter%start(this%interactions)
813 do while (iter%has_next())
814 select type (interaction => iter%get_next())
816 if (allocated(this%hm%medium_boxes) .and. .not. this%hm%medium_boxes_initialized) then
817 iint = iint + 1
818 this%hm%medium_boxes(iint) = interaction%medium_box
819 end if
820 end select
821 end do
822
823 if (allocated(this%hm%medium_boxes) .and. .not. this%hm%medium_boxes_initialized) then
824 call set_medium_rs_state(this%st, this%gr, this%hm)
825 this%hm%medium_boxes_initialized = .true.
826 end if
827
828 if (this%hm%medium_boxes_initialized .and. this%hm%operator == faraday_ampere) then
829 message(1) = "A linear medium has been defined in the input file but the Hamiltonian"
830 message(2) = "type you specified is not capable of dealing with the medium."
831 message(3) = "Please use MaxwellHamiltonianOperator = faraday_ampere_medium or simple to enable"
832 message(4) = "the medium propagation."
833 call messages_fatal(4, namespace=this%namespace)
834 end if
835
836 if (.not. this%hm%medium_boxes_initialized .and. this%hm%operator == faraday_ampere_medium) then
837 message(1) = "The variable MaxwellHamiltonianOperator has been defined as faraday_ampere_medium"
838 message(2) = "in the input file but no linear medium has been defined in the system block."
839 message(3) = "Please either use a different option for MaxwellHamiltonianOperator or add"
840 message(4) = "a linear medium to the system block."
841 call messages_fatal(4, namespace=this%namespace)
842 end if
843
846
847 ! ---------------------------------------------------------
848 subroutine maxwell_restart_write_data(this)
849 class(maxwell_t), intent(inout) :: this
851 integer :: ierr, err, zff_dim, id, id1, id2, ip_in, offset, iout
852 logical :: pml_check, write_previous_state
853 complex(real64), allocatable :: zff(:,:)
854 type(interaction_iterator_t) :: iter
855 class(interaction_t), pointer :: interaction
856
858 call profiling_in(trim(this%namespace%get())//":"//"RESTART_WRITE")
859
860 ierr = 0
861
862 if (mpi_grp_is_root(mpi_world)) then
863 do iout = 1, out_maxwell_max
864 if (this%write_handler%out(iout)%write) then
865 call write_iter_flush(this%write_handler%out(iout)%handle)
866 end if
867 end do
868 end if
870 if (.not. restart_skip(this%restart_dump)) then
871 pml_check = any(this%hm%bc%bc_ab_type(1:3) == option__maxwellabsorbingboundaries__cpml)
872
873 message(1) = "Debug: Writing td_maxwell restart."
874 call messages_info(1, namespace=this%namespace, debug_only=.true.)
875
876 if (this%tr_mxll%bc_plane_waves) then
877 zff_dim = 2 * this%st%dim
878 else
879 zff_dim = 1 * this%st%dim
880 end if
881 if (pml_check) then
882 zff_dim = zff_dim + 18
883 end if
884 select type (prop => this%algo)
885 type is (propagator_leapfrog_t)
886 write_previous_state = .true.
887 zff_dim = zff_dim + this%st%dim
888 class default
889 write_previous_state = .false.
890 end select
891 if (pml_check .and. accel_is_enabled()) then
892 call accel_read_buffer(this%hm%bc%pml%buff_conv_plus, &
893 int(this%hm%bc%pml%points_number, int64)*3*3, this%hm%bc%pml%conv_plus)
894 call accel_read_buffer(this%hm%bc%pml%buff_conv_minus, &
895 int(this%hm%bc%pml%points_number, int64)*3*3, this%hm%bc%pml%conv_minus)
896 end if
897
898 call mxll_get_batch(this%st%rs_stateb, this%st%rs_state, this%gr%np, this%st%dim)
899 if (write_previous_state) then
900 call mxll_get_batch(this%st%rs_state_prevb, this%st%rs_state_prev, this%gr%np, this%st%dim)
901 end if
902
903 safe_allocate(zff(1:this%gr%np,1:zff_dim))
904 zff = m_z0
905
906 zff(1:this%gr%np, 1:this%st%dim) = this%st%rs_state(1:this%gr%np, 1:this%st%dim)
907 if (this%tr_mxll%bc_plane_waves) then
908 call mxll_get_batch(this%st%rs_state_plane_wavesb, this%st%rs_state_plane_waves, &
909 this%gr%np, this%st%dim)
910 zff(1:this%gr%np, this%st%dim+1:this%st%dim+this%st%dim) = &
911 this%st%rs_state_plane_waves(1:this%gr%np, 1:this%st%dim)
912 offset = 2*this%st%dim
913 else
914 offset = this%st%dim
915 end if
916 if (pml_check) then
917 id = 0
918 do id1 = 1, 3
919 do id2 = 1, 3
920 id = id + 1
921 do ip_in = 1, this%hm%bc%pml%points_number
922 zff(ip_in, offset+id) = this%hm%bc%pml%conv_plus(ip_in, id1, id2)
923 zff(ip_in, offset+9+id) = this%hm%bc%pml%conv_minus(ip_in, id1, id2)
924 end do
925 end do
926 end do
927 offset = offset + 18
928 end if
929 if (write_previous_state) then
930 zff(1:this%gr%np, offset+1:offset+this%st%dim) = &
931 this%st%rs_state_prev(1:this%gr%np, 1:this%st%dim)
932 end if
933
934 call states_mxll_dump(this%restart_dump, this%st, this%space, this%gr, zff, zff_dim, err, this%iteration%counter())
935 if (err /= 0) ierr = ierr + 1
936
937 if (this%hm%current_density_from_medium) then
938 !call this%current_interpolation%write_restart(this%gr, this%space, this%restart_dump, err)
939 call iter%start(this%interactions)
940 do while (iter%has_next())
941 interaction => iter%get_next()
942 select type (interaction)
943 class is (current_to_mxll_field_t)
944 call interaction%write_restart(this%gr, this%space, this%restart_dump, err)
945 end select
946 end do
947 if (err /= 0) ierr = ierr + 1
948 end if
949
950 message(1) = "Debug: Writing td_maxwell restart done."
951 call messages_info(1, namespace=this%namespace, debug_only=.true.)
952
953 safe_deallocate_a(zff)
954
955 if (ierr /=0) then
956 message(1) = "Unable to write time-dependent Maxwell restart information."
957 call messages_warning(1, namespace=this%namespace)
958 end if
959 end if
960
961 call profiling_out(trim(this%namespace%get())//":"//"RESTART_WRITE")
963 end subroutine maxwell_restart_write_data
964
965 ! ---------------------------------------------------------
966 ! this function returns true if restart data could be read
967 logical function maxwell_restart_read_data(this)
968 class(maxwell_t), intent(inout) :: this
969
970 integer :: ierr, err, zff_dim, id, id1, id2, ip_in, offset
971 logical :: pml_check, read_previous_state
972 complex(real64), allocatable :: zff(:,:)
973 type(interaction_iterator_t) :: iter
974 class(interaction_t), pointer :: interaction
975
977 call profiling_in(trim(this%namespace%get())//":"//"RESTART_READ")
978
979 if (.not. restart_skip(this%restart)) then
980 ierr = 0
981 pml_check = any(this%hm%bc%bc_ab_type(1:3) == option__maxwellabsorbingboundaries__cpml)
982
983 if (restart_skip(this%restart)) then
984 ierr = -1
985 call profiling_out(trim(this%namespace%get())//":"//"RESTART_READ")
987 return
988 end if
989
990 message(1) = "Debug: Reading td_maxwell restart."
991 call messages_info(1, namespace=this%namespace, debug_only=.true.)
992
993 if (this%tr_mxll%bc_plane_waves) then
994 zff_dim = 2 * this%st%dim
995 else
996 zff_dim = 1 * this%st%dim
997 end if
998 if (pml_check) then
999 zff_dim = zff_dim + 18
1000 end if
1001 select type (prop => this%algo)
1002 type is (propagator_leapfrog_t)
1003 read_previous_state = .true.
1004 zff_dim = zff_dim + this%st%dim
1005 class default
1006 read_previous_state = .false.
1007 end select
1008
1009 safe_allocate(zff(1:this%gr%np,1:zff_dim))
1010
1011 call states_mxll_load(this%restart, this%st, this%gr, this%namespace, this%space, zff, &
1012 zff_dim, err, label = ": td_maxwell")
1013 this%st%rs_current_density_restart = .true.
1014
1015 this%st%rs_state(1:this%gr%np,1:this%st%dim) = zff(1:this%gr%np, 1:this%st%dim)
1016 if (this%tr_mxll%bc_plane_waves) then
1017 this%st%rs_state_plane_waves(1:this%gr%np,1:this%st%dim) = &
1018 zff(1:this%gr%np,this%st%dim+1:this%st%dim+3)
1019 offset = 2*this%st%dim
1020 else
1021 offset = this%st%dim
1022 end if
1023 if (pml_check) then
1024 id = 0
1025 do id1 = 1, 3
1026 do id2 = 1, 3
1027 id = id+1
1028 do ip_in = 1, this%hm%bc%pml%points_number
1029 this%hm%bc%pml%conv_plus(ip_in,id1,id2) = zff(ip_in, offset+ id)
1030 this%hm%bc%pml%conv_minus(ip_in,id1,id2) = zff(ip_in, offset+9+id)
1031 end do
1032 end do
1033 end do
1034 this%hm%bc%pml%conv_plus_old = this%hm%bc%pml%conv_plus
1035 this%hm%bc%pml%conv_minus_old = this%hm%bc%pml%conv_minus
1036 offset = offset + 18
1037 end if
1038 if (read_previous_state) then
1039 this%st%rs_state_prev(1:this%gr%np, 1:this%st%dim) = &
1040 zff(1:this%gr%np, offset+1:offset+this%st%dim)
1041 end if
1042
1043 if (err /= 0) then
1044 ierr = ierr + 1
1045 end if
1046
1047 message(1) = "Debug: Reading td restart done."
1048 call messages_info(1, namespace=this%namespace, debug_only=.true.)
1049
1050 safe_deallocate_a(zff)
1051
1052 if (pml_check .and. accel_is_enabled()) then
1053 call accel_write_buffer(this%hm%bc%pml%buff_conv_plus, &
1054 int(this%hm%bc%pml%points_number, int64)*3*3, this%hm%bc%pml%conv_plus)
1055 call accel_write_buffer(this%hm%bc%pml%buff_conv_minus, &
1056 int(this%hm%bc%pml%points_number, int64)*3*3, this%hm%bc%pml%conv_minus)
1057 call accel_write_buffer(this%hm%bc%pml%buff_conv_plus_old, &
1058 int(this%hm%bc%pml%points_number, int64)*3*3, this%hm%bc%pml%conv_plus_old)
1059 end if
1061 if (this%hm%current_density_from_medium) then
1062 !call this%current_interpolation%read_restart(this%gr, this%space, this%restart, err)
1063 call iter%start(this%interactions)
1064 do while (iter%has_next())
1065 interaction => iter%get_next()
1066 select type (interaction)
1067 class is (current_to_mxll_field_t)
1068 call interaction%read_restart(this%gr, this%space, this%restart, err)
1069 end select
1070 end do
1071 if (err /= 0) then
1072 ierr = ierr + 1
1073 end if
1074 end if
1075
1076 ! set batches
1077 call mxll_set_batch(this%st%rs_stateb, this%st%rs_state, this%gr%np, this%st%dim)
1078 if (read_previous_state) then
1079 call mxll_set_batch(this%st%rs_state_prevb, this%st%rs_state_prev, this%gr%np, this%st%dim)
1080 end if
1081 call batch_set_zero(this%st%inhomogeneousb)
1082 if (this%tr_mxll%bc_plane_waves) then
1083 call mxll_set_batch(this%st%rs_state_plane_wavesb, this%st%rs_state_plane_waves, this%gr%np, this%st%dim)
1084 end if
1085
1086 this%st%fromScratch = .false.
1088 else
1089 message(1) = "Unable to read time-dependent Maxwell restart information: Starting from scratch"
1090 call messages_warning(1, namespace=this%namespace)
1092 end if
1093
1094 call profiling_out(trim(this%namespace%get())//":"//"RESTART_READ")
1096 end function maxwell_restart_read_data
1097
1098 subroutine maxwell_update_kinetic_energy(this)
1099 class(maxwell_t), intent(inout) :: this
1100
1102
1103 ! the energy has already been computed at the end of the timestep
1104
1105 ! the energy of the EM wave is computed and stored in energy_mxll%energy;
1106 ! energy_mxll%energy = energy_mxll%e_energy + energy_mxll%b_energy
1107 ! here this%hm%energy is 'energy_mxll'
1108 this%kinetic_energy = this%hm%energy%energy
1109
1111 end subroutine maxwell_update_kinetic_energy
1112
1113 ! ---------------------------------------------------------
1114 subroutine maxwell_exec_end_of_timestep_tasks(this)
1115 class(maxwell_t), intent(inout) :: this
1116
1118 call profiling_in(trim(this%namespace%get())//":"//"END_TIMESTEP")
1119
1120 ! calculate Maxwell energy
1121 call energy_mxll_calc_batch(this%gr, this%st, this%hm, this%hm%energy, this%st%rs_stateb, this%st%rs_state_plane_wavesb)
1122
1123 ! get RS state values for selected points
1124 call get_rs_state_batch_selected_points(this%st%selected_points_rs_state(:,:), this%st%rs_stateb, &
1125 this%st, this%gr)
1126
1127 call profiling_out(trim(this%namespace%get())//":"//"END_TIMESTEP")
1130
1131 ! ---------------------------------------------------------
1132 subroutine maxwell_get_current(this, time, current)
1133 class(maxwell_t), intent(inout) :: this
1134 real(real64), intent(in) :: time
1135 complex(real64), contiguous, intent(inout) :: current(:, :)
1136
1137 type(interaction_iterator_t) :: iter
1138 complex(real64), allocatable :: current_density_ext(:, :)
1139
1140 push_sub(maxwell_get_current)
1141
1142 safe_allocate(current_density_ext(1:this%gr%np, 1:this%st%dim))
1143 current = m_z0
1144 if (this%hm%current_density_from_medium) then
1145 ! interpolate current from interaction
1146 call iter%start(this%interactions)
1147 do while (iter%has_next())
1148 select type (interaction => iter%get_next())
1149 class is (current_to_mxll_field_t)
1150 call interaction%interpolate(time, current_density_ext)
1151 call lalg_axpy(this%gr%np, 3, m_one, current_density_ext * this%hm%current_factor, current)
1152 end select
1153 end do
1154 end if
1155 ! calculation of external RS density
1156 if (this%hm%current_density_ext_flag) then
1157 call get_rs_density_ext(this%st, this%space, this%gr, time, current_density_ext)
1158 call lalg_axpy(this%gr%np, 3, m_one, current_density_ext, current)
1159 end if
1160 safe_deallocate_a(current_density_ext)
1161
1162 pop_sub(maxwell_get_current)
1163 end subroutine maxwell_get_current
1164
1165 ! ---------------------------------------------------------
1166 subroutine maxwell_finalize(this)
1167 type(maxwell_t), intent(inout) :: this
1168
1169
1170 push_sub(maxwell_finalize)
1171
1172 call profiling_in("MAXWELL_FINALIZE")
1173
1174 call system_end(this)
1175
1176 ! free memory
1177 safe_deallocate_a(this%rs_state_init)
1178 call hamiltonian_mxll_end(this%hm)
1179
1180 call multicomm_end(this%mc)
1181
1182 call this%st%rs_stateb%end()
1183 call this%st%rs_state_prevb%end()
1184 call this%st%inhomogeneousb%end()
1185 if (this%tr_mxll%bc_plane_waves) then
1186 call this%st%rs_state_plane_wavesb%end()
1187 end if
1188
1189 call states_mxll_end(this%st)
1190
1191 call grid_end(this%gr)
1192
1193 call restart_end(this%restart)
1194 call restart_end(this%restart_dump)
1195
1196 call poisson_end(this%st%poisson)
1197
1198 call profiling_out("MAXWELL_FINALIZE")
1199
1200 pop_sub(maxwell_finalize)
1201 end subroutine maxwell_finalize
1202
1203end module maxwell_oct_m
1204
1205!! Local Variables:
1206!! mode: f90
1207!! coding: utf-8
1208!! End:
scale a batch by a constant or vector
Definition: batch_ops.F90:162
constant times a vector plus a vector
Definition: lalg_basic.F90:171
pure logical function, public accel_is_enabled()
Definition: accel.F90:401
integer, parameter, public accel_mem_read_only
Definition: accel.F90:183
This module defines the abstract interfact for algorithm factories.
This module implements the basic elements defining algorithms.
Definition: algorithm.F90:141
character(len=algo_label_len), parameter, public iteration_done
Definition: algorithm.F90:172
This module implements batches of mesh functions.
Definition: batch.F90:133
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
Definition: batch.F90:1780
This module implements common operations on batches of mesh functions.
Definition: batch_ops.F90:116
subroutine, public batch_set_zero(this, np, async)
fill all mesh functions of the batch with zero
Definition: batch_ops.F90:242
This module handles the calculation mode.
type(calc_mode_par_t), public calc_mode_par
Singleton instance of parallel calculation mode.
This module calculates the derivatives (gradients, Laplacians, etc.) of a function.
subroutine, public get_rs_density_ext(st, space, mesh, time, rs_current_density_ext)
subroutine, public external_current_init(st, space, namespace, mesh)
This module implements the field transfer.
real(real64), parameter, public m_zero
Definition: global.F90:188
real(real64), parameter, public m_fourth
Definition: global.F90:197
complex(real64), parameter, public m_z0
Definition: global.F90:198
real(real64), parameter, public m_half
Definition: global.F90:194
real(real64), parameter, public p_c
Electron gyromagnetic ratio, see Phys. Rev. Lett. 130, 071801 (2023)
Definition: global.F90:224
real(real64), parameter, public m_one
Definition: global.F90:189
real(real64), parameter, public m_three
Definition: global.F90:191
This module implements the underlying real-space grid.
Definition: grid.F90:117
subroutine, public grid_init_stage_1(gr, namespace, space, symm, latt, n_sites, site_position)
First stage of the grid initialization.
Definition: grid.F90:194
subroutine, public grid_init_stage_2(gr, namespace, space, mc, qvector)
Second stage of the grid initialization.
Definition: grid.F90:465
subroutine, public grid_end(gr)
finalize a grid object
Definition: grid.F90:493
subroutine, public hamiltonian_mxll_init(hm, namespace, gr, st)
Initializing the Maxwell Hamiltonian.
integer, parameter, public faraday_ampere
subroutine, public hamiltonian_mxll_end(hm)
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)
The Helmholtz decomposition is intended to contain "only mathematical" functions and procedures to co...
This module implements the index, used for the mesh points.
Definition: index.F90:122
integer, parameter, public mxll_vec_pot_to_matter
integer, parameter, public linear_medium_to_em_field
integer, parameter, public lorentz_force
integer, parameter, public mxll_b_field_to_matter
integer, parameter, public mxll_e_field_to_matter
integer, parameter, public current_to_mxll_field
This module defines the abstract interaction_t class, and some auxiliary classes for interactions.
Definition: io.F90:114
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:115
subroutine, public bc_mxll_init(bc, namespace, space, gr, st)
subroutine, public surface_grid_points_mapping(mesh, st, bounds)
subroutine, public bc_mxll_initialize_pml_simple(bc, space, gr, c_factor, dt)
subroutine, public inner_and_outer_points_mapping(mesh, st, bounds)
subroutine maxwell_update_interactions_start(this)
Definition: maxwell.F90:870
subroutine maxwell_exec_end_of_timestep_tasks(this)
Definition: maxwell.F90:1208
subroutine maxwell_init_interaction(this, interaction)
Definition: maxwell.F90:302
subroutine maxwell_output_start(this)
Definition: maxwell.F90:781
subroutine maxwell_update_interactions_finish(this)
Definition: maxwell.F90:896
subroutine maxwell_init_interaction_as_partner(partner, interaction)
Definition: maxwell.F90:676
subroutine maxwell_output_finish(this)
Definition: maxwell.F90:851
subroutine maxwell_restart_write_data(this)
Definition: maxwell.F90:942
logical function maxwell_do_algorithmic_operation(this, operation, updated_quantities)
Definition: maxwell.F90:518
integer, parameter, public multigrid_mx_to_ma_large
Definition: maxwell.F90:188
subroutine maxwell_init_parallelization(this, grp)
Definition: maxwell.F90:344
subroutine maxwell_get_current(this, time, current)
Definition: maxwell.F90:1226
subroutine, public maxwell_init(this, namespace)
Definition: maxwell.F90:265
subroutine maxwell_initialize(this)
Definition: maxwell.F90:444
subroutine maxwell_output_write(this)
Definition: maxwell.F90:809
logical function maxwell_is_tolerance_reached(this, tol)
Definition: maxwell.F90:664
class(maxwell_t) function, pointer maxwell_constructor(namespace)
Definition: maxwell.F90:250
subroutine maxwell_update_kinetic_energy(this)
Definition: maxwell.F90:1192
logical function maxwell_restart_read_data(this)
Definition: maxwell.F90:1061
subroutine maxwell_finalize(this)
Definition: maxwell.F90:1260
subroutine maxwell_copy_quantities_to_interaction(partner, interaction)
Definition: maxwell.F90:700
subroutine, public mesh_interpolation_init(this, mesh)
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
integer(int64) function, public mesh_global_index_from_coords(mesh, ix)
This function returns the true global index of the point for a given vector of integer coordinates.
Definition: mesh.F90:920
integer function, public mesh_nearest_point(mesh, pos, dmin, rankmin)
Returns the index of the point which is nearest to a given vector position pos.
Definition: mesh.F90:380
subroutine, public messages_print_with_emphasis(msg, iunit, namespace)
Definition: messages.F90:921
subroutine, public messages_not_implemented(feature, namespace)
Definition: messages.F90:1114
character(len=512), private msg
Definition: messages.F90:166
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:538
subroutine, public messages_obsolete_variable(namespace, name, rep)
Definition: messages.F90:1046
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:161
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:415
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:617
logical function mpi_grp_is_root(grp)
Is the current MPI process of grpcomm, root.
Definition: mpi.F90:430
type(mpi_grp_t), public mpi_world
Definition: mpi.F90:266
This module handles the communicators for the various parallelization strategies.
Definition: multicomm.F90:145
subroutine, public multicomm_end(mc)
Definition: multicomm.F90:897
subroutine, public multicomm_init(mc, namespace, base_grp, mode_para, n_node, index_range, min_range)
create index and domain communicators
Definition: multicomm.F90:268
integer, parameter, public mxll_field_total
integer, parameter, public mxll_field_trans
integer, parameter, public mxll_field_long
this module contains the low-level part of the output system
Definition: output_low.F90:115
subroutine, public output_mxll_init(outp, namespace, space)
this module contains the output system
Definition: output.F90:115
logical function, public parse_is_defined(namespace, name)
Definition: parser.F90:502
subroutine, public poisson_init(this, namespace, space, der, mc, stencil, qtot, label, solver, verbose, force_serial, force_cmplx)
Definition: poisson.F90:233
subroutine, public poisson_end(this)
Definition: poisson.F90:709
subroutine, public profiling_out(label)
Increment out counter and sum up difference between entry and exit time.
Definition: profiling.F90:623
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
Definition: profiling.F90:552
character(len=algo_label_len), parameter, public exp_gauss1_start
character(len=algo_label_len), parameter, public exp_gauss1_finish
character(len=algo_label_len), parameter, public exp_gauss1_extrapolate
character(len=algo_label_len), parameter, public exp_gauss1_propagate
character(len=algo_label_len), parameter, public exp_gauss2_finish
character(len=algo_label_len), parameter, public exp_gauss2_propagate
character(len=algo_label_len), parameter, public exp_gauss2_extrapolate
character(len=algo_label_len), parameter, public exp_gauss2_start
character(len=algo_label_len), parameter, public expmid_extrapolate
character(len=algo_label_len), parameter, public expmid_finish
character(len=algo_label_len), parameter, public expmid_start
character(len=algo_label_len), parameter, public expmid_propagate
character(len=algo_label_len), parameter, public leapfrog_propagate
character(len=algo_label_len), parameter, public leapfrog_finish
character(len=algo_label_len), parameter, public leapfrog_start
subroutine, public mxll_propagate_expgauss2(hm, namespace, gr, space, st, tr, time, dt)
Exponential propagation scheme with Gauss collocation points, s=2.
subroutine, public energy_mxll_calc_batch(gr, st, hm, energy_mxll, rs_fieldb, rs_field_plane_wavesb)
integer, parameter, public rs_trans_forward
subroutine, public mxll_propagate_leapfrog(hm, namespace, gr, space, st, tr, time, dt, counter)
subroutine, public mxll_propagate_expgauss1(hm, namespace, gr, space, st, tr, time, dt)
Exponential propagation scheme with Gauss collocation points, s=1.
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 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, public transform_rs_densities(hm, mesh, rs_charge_density, rs_current_density, ff_density, sign)
subroutine, public energy_mxll_calc(gr, st, hm, energy_mxll, rs_field, rs_field_plane_waves)
subroutine, public propagator_mxll_init(gr, namespace, st, hm, tr)
This module implements the basic propagator framework.
Definition: propagator.F90:117
character(len=algo_label_len), parameter, public store_current_status
Definition: propagator.F90:169
This module defines the quantity_t class and the IDs for quantities, which can be exposed by a system...
Definition: quantity.F90:138
Implementation details for regridding.
Definition: regridding.F90:170
logical function, public clean_stop(comm)
returns true if a file named stop exists
Definition: restart.F90:276
subroutine, public restart_init(restart, namespace, data_type, type, mc, ierr, mesh, dir, exact)
Initializes a restart object.
Definition: restart.F90:516
integer, parameter, public restart_type_dump
Definition: restart.F90:245
logical pure function, public restart_skip(restart)
Returns true if the restart information should neither be read nor written. This might happen because...
Definition: restart.F90:969
integer, parameter, public restart_td
Definition: restart.F90:200
integer, parameter, public restart_type_load
Definition: restart.F90:245
subroutine, public restart_end(restart)
Definition: restart.F90:722
This module is intended to contain "only mathematical" functions and procedures.
Definition: sort.F90:117
subroutine, public mxll_set_batch(rs_stateb, rs_state, np, dim, offset)
subroutine, public get_rs_state_at_point(rs_state_point, rs_state, pos, st, mesh)
subroutine, public get_electric_field_vector(rs_state_vector, electric_field_vector, ep_element)
subroutine, public get_transverse_rs_state(helmholtz, st, namespace)
subroutine, public mxll_get_batch(rs_stateb, rs_state, np, dim, offset)
subroutine, public get_rs_state_batch_selected_points(rs_state_point, rs_stateb, st, mesh)
subroutine, public states_mxll_end(st)
subroutine, public states_mxll_allocate(st, mesh)
Allocates the Maxwell states defined within a states_mxll_t structure.
subroutine, public states_mxll_init(st, namespace, space)
subroutine, public get_electric_field_state(rs_state, mesh, electric_field, ep_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 states_mxll_load(restart, st, mesh, namespace, space, zff, zff_dim, ierr, iter, lowest_missing, label, verbose)
subroutine, public states_mxll_dump(restart, st, space, mesh, zff, zff_dim, ierr, iter, st_start_writing, verbose)
subroutine, public states_mxll_read_user_def(namespace, space, mesh, der, st, bc, user_def_rs_state)
This module implements the abstract system type.
Definition: system.F90:118
subroutine, public system_init_parallelization(this, grp)
Basic functionality: copy the MPI group. This function needs to be implemented by extended types that...
Definition: system.F90:1236
subroutine, public system_end(this)
Definition: system.F90:1145
integer, parameter, public out_maxwell_max
Definition: td_write.F90:291
subroutine, public td_write_mxll_init(writ, namespace, iter, dt)
Definition: td_write.F90:3628
subroutine, public td_write_mxll_end(writ)
Definition: td_write.F90:3816
subroutine, public td_write_mxll_free_data(writ, namespace, space, gr, st, hm, helmholtz, outp, iter, time)
Definition: td_write.F90:4274
subroutine, public td_write_mxll_iter(writ, space, gr, st, hm, helmholtz, dt, iter, namespace)
Definition: td_write.F90:3832
type(type_t), parameter, public type_integer
Definition: types.F90:135
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:132
This module defines the unit system, used for input and output.
Descriptor of one algorithmic operation.
Definition: algorithm.F90:163
Class to transfer a current to a Maxwell field.
These class extend the list and list iterator to make an interaction list.
abstract interaction class
surrogate interaction class to avoid circular dependencies between modules.
Lorenz force between a systems of particles and an electromagnetic field.
Class describing Maxwell systems.
Definition: maxwell.F90:194
This is defined even when running serial.
Definition: mpi.F90:142
class to transfer a Maxwell B field to a matter system
class to transfer a Maxwell field to a medium
class to transfer a Maxwell vector potential to a medium
Implements the an exponential RK scheme with Gauss collocation points, s=1 see also Hochbruck,...
Implements the an exponential RK scheme with Gauss collocation points, s=2 see also Hochbruck,...
Implements the explicit exponential midpoint propagator (without predictor-corrector)
Implements a propagator for the leap frog algorithm.
Abstract class implementing propagators.
Definition: propagator.F90:138
Systems (system_t) can expose quantities that can be used to calculate interactions with other system...
Definition: quantity.F90:171
Abstract class for systems.
Definition: system.F90:172
int true(void)