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
114
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 :: initialize => maxwell_initialize
132 procedure :: do_algorithmic_operation => maxwell_do_algorithmic_operation
133 procedure :: is_tolerance_reached => maxwell_is_tolerance_reached
134 procedure :: init_interaction_as_partner => maxwell_init_interaction_as_partner
135 procedure :: copy_quantities_to_interaction => maxwell_copy_quantities_to_interaction
136 procedure :: output_start => maxwell_output_start
137 procedure :: output_write => maxwell_output_write
138 procedure :: output_finish => maxwell_output_finish
139 procedure :: update_interactions_start => maxwell_update_interactions_start
140 procedure :: update_interactions_finish => maxwell_update_interactions_finish
141 procedure :: restart_write_data => maxwell_restart_write_data
142 procedure :: restart_read_data => maxwell_restart_read_data
143 procedure :: update_kinetic_energy => maxwell_update_kinetic_energy
144 procedure :: get_current => maxwell_get_current
145 final :: maxwell_finalize
146 end type maxwell_t
147
148 interface maxwell_t
149 procedure maxwell_constructor
150 end interface maxwell_t
151
152contains
153
154 ! ---------------------------------------------------------
155 function maxwell_constructor(namespace, grp) result(sys)
156 class(maxwell_t), pointer :: sys
157 type(namespace_t), intent(in) :: namespace
158 type(mpi_grp_t), intent(in) :: grp
159
160 push_sub(maxwell_constructor)
161
162 allocate(sys)
163
164 call maxwell_init(sys, namespace, grp)
165
166 pop_sub(maxwell_constructor)
167 end function maxwell_constructor
168
169
170 ! ---------------------------------------------------------
171 subroutine maxwell_init(this, namespace, grp)
172 class(maxwell_t), intent(inout) :: this
173 type(namespace_t), intent(in) :: namespace
174 type(mpi_grp_t), intent(in) :: grp
175
176
177 push_sub(maxwell_init)
178
179 call profiling_in("MAXWELL_INIT")
180
181 this%namespace = namespace
182 call this%init_parallelization(grp)
183
184 call messages_obsolete_variable(this%namespace, 'SystemName')
185 call messages_print_with_emphasis(msg='Maxwell Simulation Start', namespace=this%namespace)
186 this%space = space_t(this%namespace)
187
188 call grid_init_stage_1(this%gr, this%namespace, this%space, this%grp, latt=lattice_vectors_t(namespace, this%space))
189 call states_mxll_init(this%st, this%namespace, this%space)
190 this%latt = lattice_vectors_t(this%namespace, this%space)
191
192 call this%quantities%add(quantity_t("E field", updated_on_demand = .false.))
193 call this%quantities%add(quantity_t("B field", updated_on_demand = .false.))
194 call this%quantities%add(quantity_t("vector potential", updated_on_demand = .false.))
195
196 call mesh_interpolation_init(this%mesh_interpolate, this%gr)
198 this%supported_interactions = [linear_medium_to_em_field, current_to_mxll_field]
201 call profiling_out("MAXWELL_INIT")
205 pop_sub(maxwell_init)
206 end subroutine maxwell_init
208 ! ---------------------------------------------------------
209 subroutine maxwell_init_interaction(this, interaction)
210 class(maxwell_t), target, intent(inout) :: this
211 class(interaction_t), intent(inout) :: interaction
213 integer :: depth
217 select type (interaction)
219 call interaction%init(this%gr)
221 call interaction%init(this%gr, this%st%dim)
222 call interaction%init_space_latt(this%space, this%latt)
223
224 ! set interpolation depth for interaction
225 ! interpolation depth depends on the propagator
226 select type (prop => this%algo)
228 depth = 2
230 depth = 2
232 depth = 2
234 depth = 4
235 class default
236 message(1) = "The chosen propagator does not yet support interaction interpolation"
237 call messages_fatal(1, namespace=this%namespace)
238 end select
239 call interaction%init_interpolation(depth, interaction%label, cmplx=.true.)
241 this%hm%current_density_from_medium = .true.
242 class default
243 message(1) = "Trying to initialize an unsupported interaction by Maxwell."
244 call messages_fatal(1, namespace=this%namespace)
245 end select
246
248 end subroutine maxwell_init_interaction
249
250 ! ---------------------------------------------------------
251 subroutine maxwell_init_parallelization(this)
252 class(maxwell_t), intent(inout) :: this
253
254 integer(int64) :: index_range(4)
255 integer :: ierr, ip, pos_index, rankmin
256 real(real64) :: dmin
257
259
260
261 ! store the ranges for these two indices (serves as initial guess
262 ! for parallelization strategy)
263 index_range(1) = this%gr%np_global ! Number of points in mesh
264 index_range(2) = this%st%nst ! Number of states
265 index_range(3) = 1 ! Number of k-points
266 index_range(4) = 100000 ! Some large number
267
268 ! create index and domain communicators
269 call multicomm_init(this%mc, this%namespace, mpi_world, calc_mode_par, mpi_world%size, &
270 index_range, (/ 5000, 1, 1, 1 /))
271
272 call mpi_grp_init(this%st%system_grp, this%mc%master_comm)
273
274 call grid_init_stage_2(this%gr, this%namespace, this%space, this%mc)
275 call output_mxll_init(this%outp, this%namespace, this%space)
276 call hamiltonian_mxll_init(this%hm, this%namespace, this%gr, this%st)
277
278 this%st%energy_rate = m_zero
279 this%st%delta_energy = m_zero
280 this%st%energy_via_flux_calc = m_zero
281 this%st%trans_energy_rate = m_zero
282 this%st%trans_delta_energy = m_zero
283 this%st%trans_energy_via_flux_calc = m_zero
284 this%st%plane_waves_energy_rate = m_zero
285 this%st%plane_waves_delta_energy = m_zero
286 this%st%plane_waves_energy_via_flux_calc = m_zero
287
288 safe_allocate(this%rs_state_init(1:this%gr%np_part, 1:this%st%dim))
289 this%rs_state_init(:,:) = m_z0
290
291 this%energy_update_iter = 1
292
293 call poisson_init(this%st%poisson, this%namespace, this%space, this%gr%der, this%mc, this%gr%stencil)
294
295 call bc_mxll_init(this%hm%bc, this%namespace, this%space, this%gr, this%st)
296 this%bc_bounds(:,1:3) = this%hm%bc%bc_bounds(:,1:3)
297 call inner_and_outer_points_mapping(this%gr, this%st, this%bc_bounds)
298 this%dt_bounds(2, 1:3) = this%bc_bounds(1, 1:3)
299 this%dt_bounds(1, 1:3) = this%bc_bounds(1, 1:3) - this%gr%der%order * this%gr%spacing(1:3)
300 call surface_grid_points_mapping(this%gr, this%st, this%dt_bounds)
301
302 call propagator_mxll_init(this%gr, this%namespace, this%st, this%hm, this%tr_mxll)
303 call states_mxll_allocate(this%st, this%gr)
304 call external_current_init(this%st, this%space, this%namespace, this%gr)
305 this%hm%propagation_apply = .true.
306
307 ! set map for selected points
308 do ip = 1, this%st%selected_points_number
309 pos_index = mesh_nearest_point(this%gr, this%st%selected_points_coordinate(:,ip), dmin, rankmin)
310 if (this%gr%mpi_grp%rank == rankmin) then
311 this%st%selected_points_map(ip) = pos_index
312 else
313 this%st%selected_points_map(ip) = -1
314 end if
315 end do
316 if (accel_is_enabled()) then
317 call accel_create_buffer(this%st%buff_selected_points_map, accel_mem_read_only, type_integer, &
318 this%st%selected_points_number)
319 call accel_write_buffer(this%st%buff_selected_points_map, this%st%selected_points_number, &
320 this%st%selected_points_map)
321 end if
322
323 this%hm%plane_waves_apply = .true.
324 this%hm%spatial_constant_apply = .true.
325
326 call this%restart%init(this%namespace, restart_td, restart_type_load, this%mc, ierr, mesh=this%gr)
327 call this%restart_dump%init(this%namespace, restart_td, restart_type_dump, this%mc, ierr, mesh=this%gr)
328
329 ! initialize batches
330 call zbatch_init(this%st%rs_stateb, 1, 1, this%st%dim, this%gr%np_part)
331 if (this%st%pack_states) call this%st%rs_stateb%do_pack()
332 call this%st%rs_stateb%copy_to(this%st%rs_state_prevb)
333 call this%st%rs_stateb%copy_to(this%st%inhomogeneousb)
334 if (this%tr_mxll%bc_plane_waves) then
335 call this%st%rs_stateb%copy_to(this%st%rs_state_plane_wavesb)
336 end if
337
338 ! the order should depend on the propagator
339 !this%current_interpolation => time_interpolation_t(this%gr%np, this%st%dim, 2, .true., "current")
340 ! Initialize Helmholtz decomposition
341 call this%helmholtz%init(this%namespace, this%gr, this%mc, this%space)
342
344 end subroutine maxwell_init_parallelization
345
346 ! ---------------------------------------------------------
347 subroutine maxwell_initialize(this)
348 class(maxwell_t), intent(inout) :: this
349
350 real(real64) :: courant
351
352
353 push_sub(maxwell_initialize)
354
355 call profiling_in("MAXWELL_INIT_CONDITIONS")
356
357 select type (algo => this%algo)
358 class is (propagator_t)
359
360 ! The reference Courant–Friedrichs–Lewy condition for stability here
361 ! is calculated with respect to R. Courant, K. Friedrichs, and H. Lewy
362 ! On the Partial Difference Equations of
363 ! Mathematical Physics, IBM Journal of Research and Development, vol. 11
364 ! pp. 215-234, Mar. 1967.
365 courant = m_one/(p_c * sqrt(m_one/this%gr%spacing(1)**2 + m_one/this%gr%spacing(2)**2 + &
366 m_one/this%gr%spacing(3)**2))
367
368 if (algo%dt > (m_one + m_fourth) * courant) then
369 write(message(1),'(a)') 'The timestep is too large compared to the Courant-Friedrichs-Lewy'
370 write(message(2),'(a)') 'stability condition. Time propagation might not be stable.'
371 call messages_warning(2, namespace=this%namespace)
372 end if
373
374 if (parse_is_defined(this%namespace, 'UserDefinedInitialMaxwellStates')) then
375 call states_mxll_read_user_def(this%namespace, this%space, this%gr, this%gr%der, this%st, this%hm%bc,&
376 this%rs_state_init)
377 call messages_print_with_emphasis(msg="Setting initial EM field inside box", namespace=this%namespace)
378 ! TODO: add consistency check that initial state fulfills Gauss laws
379 this%st%rs_state(:,:) = this%st%rs_state + this%rs_state_init
380 if (this%tr_mxll%bc_plane_waves) then
381 this%st%rs_state_plane_waves(:,:) = this%rs_state_init
382 end if
383 end if
384
385 ! initialize the spatial constant field according to the conditions set in the
386 ! UserDefinedConstantSpatialMaxwellField block
387 if (this%tr_mxll%bc_constant) then
388 call spatial_constant_calculation(this%tr_mxll%bc_constant, this%st, this%gr, this%hm, m_zero, &
389 algo%dt, m_zero, this%st%rs_state, set_initial_state = .true.)
390 ! for mesh parallelization, this needs communication!
391 this%st%rs_state_const(:) = this%st%rs_state(mesh_global_index_from_coords(this%gr, [0,0,0]),:)
392 end if
393
394 if (parse_is_defined(this%namespace, 'UserDefinedInitialMaxwellStates')) then
395 safe_deallocate_a(this%rs_state_init)
396 end if
397
398 call hamiltonian_mxll_update(this%hm, time = m_zero)
399
400 ! calculate Maxwell energy
401 call energy_mxll_calc(this%gr, this%st, this%hm, this%hm%energy, this%st%rs_state, &
402 this%st%rs_state_plane_waves)
403
404 ! Get RS states values for selected points
405 call get_rs_state_at_point(this%st%selected_points_rs_state(:,:), this%st%rs_state, &
406 this%st%selected_points_coordinate(:,:), this%st, this%gr)
407
408 call mxll_set_batch(this%st%rs_stateb, this%st%rs_state, this%gr%np, this%st%dim)
409 call batch_set_zero(this%st%inhomogeneousb)
410 if (this%tr_mxll%bc_plane_waves) then
411 call mxll_set_batch(this%st%rs_state_plane_wavesb, this%st%rs_state_plane_waves, this%gr%np, this%st%dim)
412 end if
413 end select
414
415 call profiling_out("MAXWELL_INIT_CONDITIONS")
416
417 pop_sub(maxwell_initialize)
418 end subroutine maxwell_initialize
419
420 ! ---------------------------------------------------------
421 logical function maxwell_do_algorithmic_operation(this, operation, updated_quantities) result(done)
422 class(maxwell_t), intent(inout) :: this
423 class(algorithmic_operation_t), intent(in) :: operation
424 character(len=:), allocatable, intent(out) :: updated_quantities(:)
425
426 complex(real64), allocatable :: charge_density_ext(:)
427 real(real64) :: current_time
428
430 call profiling_in(trim(this%namespace%get())//":"//trim(operation%id))
431
432 current_time = this%iteration%value()
433
434 select type (algo => this%algo)
435 class is (propagator_t)
436
437 done = .true.
438 select case (operation%id)
440 ! For the moment we do nothing
441
443 safe_allocate(this%ff_rs_inhom_t1(1:this%gr%np_part, 1:this%hm%dim))
444 safe_allocate(this%ff_rs_inhom_t2(1:this%gr%np_part, 1:this%hm%dim))
445
446 case (expmid_finish)
447 safe_deallocate_a(this%ff_rs_inhom_t1)
448 safe_deallocate_a(this%ff_rs_inhom_t2)
449
450 case (expmid_extrapolate)
451 if (this%hm%current_density_ext_flag .or. this%hm%current_density_from_medium) then
452 call this%get_current(current_time, this%st%rs_current_density_t1)
453 call this%get_current(current_time+algo%dt, this%st%rs_current_density_t2)
454 end if
455
456 case (expmid_propagate)
457 ! Propagation
458
459 !We first compute three external charge and current densities and we convert them as RS vectors
460 safe_allocate(charge_density_ext(1:this%gr%np))
461
462 !No charge density at the moment
463 charge_density_ext = m_z0
464
465 call transform_rs_densities(this%hm, this%gr, &
466 this%st%rs_current_density_t1, this%ff_rs_inhom_t1, rs_trans_forward)
467 call transform_rs_densities(this%hm, this%gr, &
468 this%st%rs_current_density_t2, this%ff_rs_inhom_t2, rs_trans_forward)
469
470 safe_deallocate_a(charge_density_ext)
471
472 ! Propagation dt with H_maxwell
473 call mxll_propagation_step(this%hm, this%namespace, this%gr, this%space, this%st, this%tr_mxll,&
474 this%st%rs_stateb, this%ff_rs_inhom_t1, this%ff_rs_inhom_t2, current_time, algo%dt)
475
476 updated_quantities = [character(16) :: "E field", "B field", "vector potential"]
477
478 case (leapfrog_start)
479 if (any(this%hm%bc%bc_ab_type == option__maxwellabsorbingboundaries__cpml)) then
480 call bc_mxll_initialize_pml_simple(this%hm%bc, this%space, this%gr, this%hm%c_factor, algo%dt)
481 end if
482
483 case (leapfrog_finish)
484
485 case (leapfrog_propagate)
486 if (this%hm%current_density_ext_flag .or. this%hm%current_density_from_medium) then
487 call this%get_current(current_time, this%st%rs_current_density_t1)
488 call mxll_set_batch(this%st%inhomogeneousb, this%st%rs_current_density_t1, this%gr%np, this%st%dim)
489 call batch_scal(this%gr%np, -m_one, this%st%inhomogeneousb)
490 end if
491
492 call mxll_propagate_leapfrog(this%hm, this%namespace, this%gr, this%st, this%tr_mxll, &
493 current_time, algo%dt, this%iteration%counter())
494
495 updated_quantities = [character(16) :: "E field", "B field", "vector potential"]
496
497 case (exp_gauss1_start)
498 safe_allocate(this%ff_rs_inhom_t1(1:this%gr%np_part, 1:this%hm%dim))
499 if (any(this%hm%bc%bc_ab_type == option__maxwellabsorbingboundaries__cpml)) then
500 call bc_mxll_initialize_pml_simple(this%hm%bc, this%space, this%gr, this%hm%c_factor, algo%dt)
501 end if
502
503 case (exp_gauss1_finish)
504 safe_deallocate_a(this%ff_rs_inhom_t1)
505
507 if (this%hm%current_density_ext_flag .or. this%hm%current_density_from_medium) then
508 call this%get_current(current_time+algo%dt*m_half, this%st%rs_current_density_t1)
509 end if
510
512 ! the propagation step is
513 ! F_{n+1} = F_n + dt * phi_1(-i*dt*H) [-i*H F_n - J_1]
514 ! where J_1 = J(t_n + dt/2)
515 call mxll_propagate_expgauss1(this%hm, this%namespace, this%gr, this%st, this%tr_mxll, &
516 current_time, algo%dt)
517
518 updated_quantities = [character(16) :: "E field", "B field", "vector potential"]
519
520 case (exp_gauss2_start)
521 safe_allocate(this%ff_rs_inhom_t1(1:this%gr%np_part, 1:this%hm%dim))
522 safe_allocate(this%ff_rs_inhom_t2(1:this%gr%np_part, 1:this%hm%dim))
523 if (any(this%hm%bc%bc_ab_type == option__maxwellabsorbingboundaries__cpml)) then
524 call bc_mxll_initialize_pml_simple(this%hm%bc, this%space, this%gr, this%hm%c_factor, algo%dt)
525 end if
526
527 case (exp_gauss2_finish)
528 safe_deallocate_a(this%ff_rs_inhom_t1)
529 safe_deallocate_a(this%ff_rs_inhom_t2)
530
532 if (this%hm%current_density_ext_flag .or. this%hm%current_density_from_medium) then
533 call this%get_current(current_time+algo%dt*(m_half-sqrt(m_three)/6.0_real64), this%st%rs_current_density_t1)
534 call this%get_current(current_time+algo%dt*(m_half+sqrt(m_three)/6.0_real64), this%st%rs_current_density_t2)
535 end if
536
538 ! the propagation step is
539 ! F_{n+1} = F_n + dt * phi_1(-i*dt*H) [-i*H F_n - a_1 J_1 - a_2 J_2]
540 ! + dt * phi_2(-i*dt*H) [-b_1 J_1 - b_2 J_2]
541 ! where J_1 = J(t_n + (1/2 - sqrt(3)/6)*dt),
542 ! J_2 = J(t_n + (1/2 + sqrt(3)/6)*dt),
543 ! a_1 = 1/2*(1+sqrt(3)),
544 ! a_2 = 1/2*(1-sqrt(3)),
545 ! b_1 = -sqrt(3),
546 ! b_2 = sqrt(3)
547 call mxll_propagate_expgauss2(this%hm, this%namespace, this%gr, this%st, this%tr_mxll, &
548 current_time, algo%dt)
549
550 updated_quantities = [character(16) :: "E field", "B field", "vector potential"]
551
552 case (iteration_done)
554 done = .false.
555
556 case default
557 done = .false.
558 end select
559
560 end select
561
562 call profiling_out(trim(this%namespace%get())//":"//trim(operation%id))
565
566 ! ---------------------------------------------------------
567 logical function maxwell_is_tolerance_reached(this, tol) result(converged)
568 class(maxwell_t), intent(in) :: this
569 real(real64), intent(in) :: tol
570
572
573 converged = .false.
574
577
578 ! ---------------------------------------------------------
579 subroutine maxwell_init_interaction_as_partner(partner, interaction)
580 class(maxwell_t), intent(in) :: partner
581 class(interaction_surrogate_t), intent(inout) :: interaction
582
584
585 select type (interaction)
586 type is (lorentz_force_t)
587 ! Nothing to be initialized for the Lorentz force.
589 call interaction%init_from_partner(partner%gr, partner%space, partner%namespace)
591 call interaction%init_from_partner(partner%gr, partner%space, partner%namespace)
593 call interaction%init_from_partner(partner%gr, partner%space, partner%namespace)
594 class default
595 message(1) = "Unsupported interaction."
596 call messages_fatal(1, namespace=partner%namespace)
597 end select
598
601
602 ! ---------------------------------------------------------
603 subroutine maxwell_copy_quantities_to_interaction(partner, interaction)
604 class(maxwell_t), intent(inout) :: partner
605 class(interaction_surrogate_t), intent(inout) :: interaction
606
607 integer :: ip
608 complex(real64) :: interpolated_value(3)
609 real(real64), allocatable :: b_field(:,:), vec_pot(:,:)
610
612 call profiling_in(trim(partner%namespace%get())//":"//"COPY_QUANTITY_INTER")
613
614 select type (interaction)
615 type is (lorentz_force_t)
616 call mxll_get_batch(partner%st%rs_stateb, partner%st%rs_state, partner%gr%np, partner%st%dim)
617
618 do ip = 1, interaction%system_np
619 call mesh_interpolation_evaluate(partner%mesh_interpolate, partner%st%rs_state(:,1), &
620 interaction%system_pos(:, ip), interpolated_value(1))
621 call mesh_interpolation_evaluate(partner%mesh_interpolate, partner%st%rs_state(:,2), &
622 interaction%system_pos(:, ip), interpolated_value(2))
623 call mesh_interpolation_evaluate(partner%mesh_interpolate, partner%st%rs_state(:,3), &
624 interaction%system_pos(:, ip), interpolated_value(3))
625 call get_electric_field_vector(interpolated_value, interaction%partner_e_field(:, ip))
626 call get_magnetic_field_vector(interpolated_value, 1, interaction%partner_b_field(:, ip))
627 end do
628
630 call mxll_get_batch(partner%st%rs_stateb, partner%st%rs_state, partner%gr%np, partner%st%dim)
631 select case (interaction%type)
632
633 case (mxll_field_total)
634 call get_electric_field_state(partner%st%rs_state, partner%gr, interaction%partner_field)
635
636 case (mxll_field_trans)
637 call get_transverse_rs_state(partner%helmholtz, partner%st, partner%namespace)
638 call get_electric_field_state(partner%st%rs_state_trans, partner%gr, interaction%partner_field)
639
640 case (mxll_field_long)
641 call partner%helmholtz%get_long_field(partner%namespace, partner%st%rs_state_long, total_field=partner%st%rs_state)
642 call get_electric_field_state(partner%st%rs_state_long, partner%gr, interaction%partner_field)
643
644 case default
645 message(1) = "Unknown type of field requested by interaction."
646 call messages_fatal(1, namespace=partner%namespace)
647 end select
648 call interaction%do_mapping()
649
651 call mxll_get_batch(partner%st%rs_stateb, partner%st%rs_state, partner%gr%np, partner%st%dim)
652 safe_allocate(b_field(1:partner%gr%np_part, 1:partner%gr%box%dim))
653 safe_allocate(vec_pot(1:partner%gr%np_part, 1:partner%gr%box%dim))
654 ! magnetic field is always transverse
655 call get_magnetic_field_state(partner%st%rs_state, partner%gr, partner%st%rs_sign, b_field, &
656 partner%st%mu(1:partner%gr%np), partner%gr%np)
657 ! vector potential stored in partner_field
658 call partner%helmholtz%get_vector_potential(partner%namespace, vec_pot, trans_field=b_field)
659 ! in the convention used by the electronic system, the -1/c factor is included in the vector potential
660 ! but we do not divide by it, because the Maxwell code is in SI units, so the vector potential units
661 ! are suitable for a (p + q.A) minimal coupling interaction
662 interaction%partner_field(1:partner%gr%np,1:partner%gr%box%dim) = &
663 vec_pot(1:partner%gr%np,1:partner%gr%box%dim)
664 safe_deallocate_a(b_field)
665 safe_deallocate_a(vec_pot)
666 call interaction%do_mapping()
667
669 call mxll_get_batch(partner%st%rs_stateb, partner%st%rs_state, partner%gr%np, partner%st%dim)
670 call get_magnetic_field_state(partner%st%rs_state, partner%gr, partner%st%rs_sign, &
671 interaction%partner_field, partner%st%mu(1:partner%gr%np), partner%gr%np)
672 call interaction%do_mapping()
673
674 class default
675 message(1) = "Unsupported interaction."
676 call messages_fatal(1, namespace=partner%namespace)
677 end select
678
679 call profiling_out(trim(partner%namespace%get())//":"//"COPY_QUANTITY_INTER")
682
683 ! ---------------------------------------------------------
684 subroutine maxwell_output_start(this)
685 class(maxwell_t), intent(inout) :: this
686
687 real(real64) :: itr_value
688
689 push_sub(maxwell_output_start)
690 call profiling_in(trim(this%namespace%get())//":"//"OUTPUT_START")
691
692 select type (algo => this%algo)
693 class is (propagator_t)
694 call mxll_get_batch(this%st%rs_stateb, this%st%rs_state, this%gr%np, this%st%dim)
695
696 call td_write_mxll_init(this%write_handler, this%namespace, this%iteration%counter(), algo%dt, &
697 this%st%system_grp)
698 if (this%st%fromScratch) then
699 call td_write_mxll_iter(this%write_handler, this%space, this%gr, this%st, this%hm, this%helmholtz, algo%dt, &
700 this%iteration%counter(), this%namespace)
701 itr_value = this%iteration%value()
702 call td_write_mxll_free_data(this%namespace, this%space, this%gr, this%st, this%hm, this%helmholtz, &
703 this%outp, this%iteration%counter(), itr_value)
704 end if
705
706 end select
707
708 call profiling_out(trim(this%namespace%get())//":"//"OUTPUT_START")
709 pop_sub(maxwell_output_start)
710 end subroutine maxwell_output_start
711
712 ! ---------------------------------------------------------
713 subroutine maxwell_output_write(this)
714 class(maxwell_t), intent(inout) :: this
715
716 logical :: stopping, reached_output_interval
717 real(real64) :: itr_value
718 integer :: iout
719
720 push_sub(maxwell_output_write)
721 call profiling_in(trim(this%namespace%get())//":"//"OUTPUT_WRITE")
722
723 select type (algo => this%algo)
724 class is (propagator_t)
725 stopping = clean_stop(this%mc%master_comm)
726
727 call td_write_mxll_iter(this%write_handler, this%space, this%gr, this%st, this%hm, this%helmholtz, algo%dt, &
728 this%iteration%counter(), this%namespace)
729
730 reached_output_interval = .false.
731 do iout = 1, out_maxwell_max
732 if (this%outp%output_interval(iout) > 0) then
733 if (mod(this%iteration%counter(), this%outp%output_interval(iout)) == 0) then
734 reached_output_interval = .true.
735 exit
736 end if
737 end if
738 end do
739
740 if (reached_output_interval .or. stopping) then
741 call mxll_get_batch(this%st%rs_stateb, this%st%rs_state, this%gr%np, this%st%dim)
742
743 itr_value = this%iteration%value()
744 call td_write_mxll_free_data(this%namespace, this%space, this%gr, this%st, this%hm, this%helmholtz, &
745 this%outp, this%iteration%counter(), itr_value)
746 end if
747
748 end select
749
750 call profiling_out(trim(this%namespace%get())//":"//"OUTPUT_WRITE")
751 pop_sub(maxwell_output_write)
752 end subroutine maxwell_output_write
753
754 ! ---------------------------------------------------------
755 subroutine maxwell_output_finish(this)
756 class(maxwell_t), intent(inout) :: this
757
758
759 push_sub(maxwell_output_finish)
760
761 call profiling_in("MAXWELL_OUTPUT_FINISH")
762
763 select type (algo => this%algo)
764 class is (propagator_t)
765 call td_write_mxll_end(this%write_handler)
766 end select
767
768 call profiling_out("MAXWELL_OUTPUT_FINISH")
769
770 pop_sub(maxwell_output_finish)
771 end subroutine maxwell_output_finish
772
773 ! ---------------------------------------------------------
774 subroutine maxwell_update_interactions_start(this)
775 class(maxwell_t), intent(inout) :: this
776
777 type(interaction_iterator_t) :: iter
778 integer :: int_counter
781
782 int_counter = 0
783 call iter%start(this%interactions)
784 do while (iter%has_next())
785 select type (interaction => iter%get_next())
787 int_counter = int_counter + 1
788 end select
789 end do
790
791 if (int_counter /= 0 .and. .not. allocated(this%hm%medium_boxes)) then
792 safe_allocate(this%hm%medium_boxes(1:int_counter))
793 this%hm%calc_medium_box = .true.
794 end if
795
798
799 ! ---------------------------------------------------------
801 class(maxwell_t), intent(inout) :: this
802
803 type(interaction_iterator_t) :: iter
804 integer :: iint
805
807
808 iint = 0
809
810 call iter%start(this%interactions)
811 do while (iter%has_next())
812 select type (interaction => iter%get_next())
814 if (allocated(this%hm%medium_boxes) .and. .not. this%hm%medium_boxes_initialized) then
815 iint = iint + 1
816 this%hm%medium_boxes(iint) = interaction%medium_box
817 end if
818 end select
819 end do
820
821 if (allocated(this%hm%medium_boxes) .and. .not. this%hm%medium_boxes_initialized) then
822 call set_medium_rs_state(this%st, this%gr, this%hm)
823 this%hm%medium_boxes_initialized = .true.
824 end if
825
826 if (this%hm%medium_boxes_initialized .and. this%hm%operator == faraday_ampere) then
827 message(1) = "A linear medium has been defined in the input file but the Hamiltonian"
828 message(2) = "type you specified is not capable of dealing with the medium."
829 message(3) = "Please use MaxwellHamiltonianOperator = faraday_ampere_medium or simple to enable"
830 message(4) = "the medium propagation."
831 call messages_fatal(4, namespace=this%namespace)
832 end if
833
834 if (.not. this%hm%medium_boxes_initialized .and. &
835 .not. any(this%hm%bc%bc_type == mxll_bc_medium) .and. &
836 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 neither a medium boundary nor a linear medium system"
839 message(3) = "has been defined. Please either use a different option for"
840 message(4) = "MaxwellHamiltonianOperator or define a Maxwell medium."
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_world%is_root()) 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. this%restart_dump%skip()) 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 this%hm%bc%pml%points_number, 3, 3, this%hm%bc%pml%conv_plus)
894 call accel_read_buffer(this%hm%bc%pml%buff_conv_minus, &
895 this%hm%bc%pml%points_number, 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)
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. this%restart%skip()) then
980 ierr = 0
981 pml_check = any(this%hm%bc%bc_ab_type(1:3) == option__maxwellabsorbingboundaries__cpml)
982
983 if (this%restart%skip()) 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, this%hm%bc%pml%points_number, 3, 3, &
1054 this%hm%bc%pml%conv_plus)
1055 call accel_write_buffer(this%hm%bc%pml%buff_conv_minus, this%hm%bc%pml%points_number, 3, 3, &
1056 this%hm%bc%pml%conv_minus)
1057 call accel_write_buffer(this%hm%bc%pml%buff_conv_plus_old, this%hm%bc%pml%points_number, 3, 3, &
1058 this%hm%bc%pml%conv_plus_old)
1059 end if
1060
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 this%restart%end()
1194 call this%restart_dump%end()
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:167
constant times a vector plus a vector
Definition: lalg_basic.F90:173
pure logical function, public accel_is_enabled()
Definition: accel.F90:403
integer, parameter, public accel_mem_read_only
Definition: accel.F90:186
This module defines the abstract interfact for algorithm factories.
This module implements the basic elements defining algorithms.
Definition: algorithm.F90:143
character(len=algo_label_len), parameter, public iteration_done
Definition: algorithm.F90:174
This module implements batches of mesh functions.
Definition: batch.F90:135
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:1924
This module implements common operations on batches of mesh functions.
Definition: batch_ops.F90:118
subroutine, public batch_set_zero(this, np, async)
fill all mesh functions of the batch with zero
Definition: batch_ops.F90:265
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:200
real(real64), parameter, public m_fourth
Definition: global.F90:209
complex(real64), parameter, public m_z0
Definition: global.F90:210
real(real64), parameter, public m_half
Definition: global.F90:206
real(real64), parameter, public p_c
Electron gyromagnetic ratio, see Phys. Rev. Lett. 130, 071801 (2023)
Definition: global.F90:242
real(real64), parameter, public m_one
Definition: global.F90:201
real(real64), parameter, public m_three
Definition: global.F90:203
This module implements the underlying real-space grid.
Definition: grid.F90:119
subroutine, public grid_init_stage_1(gr, namespace, space, grp, symm, latt, n_sites, site_position)
First stage of the grid initialization.
Definition: grid.F90:197
subroutine, public grid_init_stage_2(gr, namespace, space, mc, qvector)
Second stage of the grid initialization.
Definition: grid.F90:476
subroutine, public grid_end(gr)
finalize a grid object
Definition: grid.F90:503
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:124
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:116
System information (time, memory, sysname)
Definition: loct.F90:117
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:117
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)
integer, parameter, public mxll_bc_medium
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:1210
subroutine maxwell_init_interaction(this, interaction)
Definition: maxwell.F90:305
subroutine maxwell_output_start(this)
Definition: maxwell.F90:780
subroutine maxwell_update_interactions_finish(this)
Definition: maxwell.F90:896
class(maxwell_t) function, pointer maxwell_constructor(namespace, grp)
Definition: maxwell.F90:251
subroutine maxwell_init_interaction_as_partner(partner, interaction)
Definition: maxwell.F90:675
subroutine maxwell_output_finish(this)
Definition: maxwell.F90:851
subroutine maxwell_restart_write_data(this)
Definition: maxwell.F90:944
logical function maxwell_do_algorithmic_operation(this, operation, updated_quantities)
Definition: maxwell.F90:517
integer, parameter, public multigrid_mx_to_ma_large
Definition: maxwell.F90:190
subroutine maxwell_get_current(this, time, current)
Definition: maxwell.F90:1228
subroutine, public maxwell_init(this, namespace, grp)
Definition: maxwell.F90:267
subroutine maxwell_initialize(this)
Definition: maxwell.F90:443
subroutine maxwell_output_write(this)
Definition: maxwell.F90:809
logical function maxwell_is_tolerance_reached(this, tol)
Definition: maxwell.F90:663
subroutine maxwell_init_parallelization(this)
Definition: maxwell.F90:347
subroutine maxwell_update_kinetic_energy(this)
Definition: maxwell.F90:1194
logical function maxwell_restart_read_data(this)
Definition: maxwell.F90:1063
subroutine maxwell_finalize(this)
Definition: maxwell.F90:1262
subroutine maxwell_copy_quantities_to_interaction(partner, interaction)
Definition: maxwell.F90:699
subroutine, public mesh_interpolation_init(this, mesh)
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:120
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:919
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:386
subroutine, public messages_print_with_emphasis(msg, iunit, namespace)
Definition: messages.F90:898
character(len=512), private msg
Definition: messages.F90:167
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:525
subroutine, public messages_obsolete_variable(namespace, name, rep)
Definition: messages.F90:1000
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:162
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:410
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:594
type(mpi_grp_t), public mpi_world
Definition: mpi.F90:272
subroutine mpi_grp_init(grp, comm)
Initialize MPI group instance.
Definition: mpi.F90:341
This module handles the communicators for the various parallelization strategies.
Definition: multicomm.F90:147
subroutine, public multicomm_end(mc)
Definition: multicomm.F90:706
subroutine, public multicomm_init(mc, namespace, base_grp, mode_para, n_node, index_range, min_range)
create index and domain communicators
Definition: multicomm.F90:273
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:117
subroutine, public output_mxll_init(outp, namespace, space)
this module contains the output system
Definition: output.F90:117
logical function, public parse_is_defined(namespace, name)
Definition: parser.F90:463
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:680
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
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 transform_rs_densities(hm, mesh, rs_current_density, ff_density, sign)
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, st, tr, time, dt, counter)
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 mxll_propagate_expgauss1(hm, namespace, gr, st, tr, time, dt)
Exponential propagation scheme with Gauss collocation points, s=1.
subroutine, public energy_mxll_calc(gr, st, hm, energy_mxll, rs_field, rs_field_plane_waves)
subroutine, public mxll_propagate_expgauss2(hm, namespace, gr, st, tr, time, dt)
Exponential propagation scheme with Gauss collocation points, s=2.
subroutine, public propagator_mxll_init(gr, namespace, st, hm, tr)
This module implements the basic propagator framework.
Definition: propagator.F90:119
character(len=algo_label_len), parameter, public store_current_status
Definition: propagator.F90:176
This module defines the quantity_t class and the IDs for quantities, which can be exposed by a system...
Definition: quantity.F90:140
Implementation details for regridding.
Definition: regridding.F90:172
logical function, public clean_stop(comm)
returns true if a file named stop exists
Definition: restart.F90:337
integer, parameter, public restart_type_dump
Definition: restart.F90:184
integer, parameter, public restart_td
Definition: restart.F90:156
integer, parameter, public restart_type_load
Definition: restart.F90:184
This module is intended to contain "only mathematical" functions and procedures.
Definition: sort.F90:119
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, 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:120
subroutine, public system_end(this)
Definition: system.F90:1151
integer, parameter, public out_maxwell_max
Definition: td_write.F90:298
subroutine, public td_write_mxll_end(writ)
Definition: td_write.F90:3921
subroutine, public td_write_mxll_init(writ, namespace, iter, dt, grp)
Definition: td_write.F90:3730
subroutine, public td_write_mxll_free_data(namespace, space, gr, st, hm, helmholtz, outp, iter, time)
Definition: td_write.F90:4381
subroutine, public td_write_mxll_iter(writ, space, gr, st, hm, helmholtz, dt, iter, namespace)
Definition: td_write.F90:3937
type(type_t), parameter, public type_integer
Definition: types.F90:137
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:134
This module defines the unit system, used for input and output.
Descriptor of one algorithmic operation.
Definition: algorithm.F90:165
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:196
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:144
Systems (system_t) can expose quantities that can be used to calculate interactions with other system...
Definition: quantity.F90:173
Abstract class for systems.
Definition: system.F90:174
int true(void)