Octopus
states_mxll.F90
Go to the documentation of this file.
1!! Copyright (C) 2019 R. Jestaedt, F. Bonafe, H. Appel, A. Rubio
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#include "global.h"
19
21 use accel_oct_m
22 use batch_oct_m
25 use comm_oct_m
26 use debug_oct_m
29 use global_oct_m
30 use grid_oct_m
32 use math_oct_m
33 use mesh_oct_m
36 use mpi_oct_m
39#ifdef HAVE_OPENMP
40 use omp_lib
41#endif
42 use parser_oct_m
46 use space_oct_m
51 use types_oct_m
52 use unit_oct_m
55
56 implicit none
57
58 private
59
60 public :: &
80
81 type :: states_mxll_t
82 ! Components are public by default
83 integer :: dim
84 integer :: rs_sign
85 logical :: pack_states
86 logical :: parallel_in_states = .false.
87 integer, public :: nst
88 logical, public :: packed
89
90 complex(real64), allocatable :: rs_state_plane_waves(:,:)
91 complex(real64), allocatable :: rs_state(:,:)
92 complex(real64), allocatable :: rs_state_prev(:,:)
93 complex(real64), allocatable :: rs_state_trans(:,:)
94 complex(real64), allocatable :: rs_state_long(:,:)
95 complex(real64), allocatable :: rs_current_density_t1(:,:)
96 complex(real64), allocatable :: rs_current_density_t2(:,:)
97
98 logical :: rs_current_density_restart = .false.
99 complex(real64), allocatable :: rs_current_density_restart_t1(:,:)
100 complex(real64), allocatable :: rs_current_density_restart_t2(:,:)
101
102 type(batch_t) :: rs_stateb
103 type(batch_t) :: rs_state_prevb
104 type(batch_t) :: inhomogeneousb
105 type(batch_t) :: rs_state_plane_wavesb
106
107 real(real64), allocatable :: ep(:)
108 real(real64), allocatable :: mu(:)
109
110 integer, allocatable :: rs_state_fft_map(:,:,:)
111 integer, allocatable :: rs_state_fft_map_inv(:,:)
112
113 real(real64) :: energy_rate
114 real(real64) :: delta_energy
115 real(real64) :: energy_via_flux_calc
116
117 real(real64) :: trans_energy_rate
118 real(real64) :: trans_delta_energy
119 real(real64) :: trans_energy_via_flux_calc
120
121 real(real64) :: plane_waves_energy_rate
122 real(real64) :: plane_waves_delta_energy
123 real(real64) :: plane_waves_energy_via_flux_calc
124
125 real(real64) :: poynting_vector_box_surface(1:2,1:3,1:3) = m_zero
126 real(real64) :: poynting_vector_box_surface_plane_waves(1:2,1:3,1:3) = m_zero
127 real(real64) :: electric_field_box_surface(1:2,1:3,1:3) = m_zero
128 real(real64) :: electric_field_box_surface_plane_waves(1:2,1:3,1:3) = m_zero
129 real(real64) :: magnetic_field_box_surface(1:2,1:3,1:3) = m_zero
130 real(real64) :: magnetic_field_box_surface_plane_waves(1:2,1:3,1:3) = m_zero
131
132 logical :: rs_state_const_external = .false.
133 complex(real64), allocatable :: rs_state_const(:)
134 complex(real64), allocatable :: rs_state_const_amp(:,:)
135 type(tdf_t), allocatable :: rs_state_const_td_function(:)
136
137 integer :: inner_points_number
138 integer, allocatable :: inner_points_map(:)
139 logical, allocatable :: inner_points_mask(:)
140 integer :: boundary_points_number
141 integer, allocatable :: boundary_points_map(:)
142 logical, allocatable :: boundary_points_mask(:)
143 type(accel_mem_t) :: buff_inner_points_map, buff_boundary_points_map
144
145 integer :: surface_points_number(3)
146 integer, allocatable :: surface_points_map(:,:,:)
147 real(real64) :: surface_element(3)
148
149 integer :: surface_grid_rows_number(3)
150 integer, allocatable :: surface_grid_points_number(:,:,:)
151 integer(int64), allocatable :: surface_grid_points_map(:,:,:,:,:)
152 integer, allocatable :: surface_grid_center(:,:,:,:)
153 real(real64) :: surface_grid_element(3)
154
155 type(mesh_plane_t) :: surface(2,3)
156
157 integer :: selected_points_number
158 real(real64), allocatable :: selected_points_coordinate(:,:)
159 complex(real64), allocatable :: selected_points_rs_state(:,:)
160 complex(real64), allocatable :: selected_points_rs_state_long(:,:)
161 complex(real64), allocatable :: selected_points_rs_state_trans(:,:)
162 integer, allocatable :: selected_points_map(:)
163 type(accel_mem_t) :: buff_selected_points_map
164 real(real64) :: rs_state_trans_var
165
166 real(real64), allocatable :: grid_rho(:,:)
167 complex(real64), allocatable :: kappa_psi(:,:)
168
169 character(len=1024), allocatable :: user_def_e_field(:)
170 character(len=1024), allocatable :: user_def_b_field(:)
171
172 integer :: energy_incident_waves_calc_iter
173 logical :: energy_incident_waves_calc
174
175 ! external current variables
176 integer :: external_current_number
177 integer, allocatable :: external_current_modus(:)
178 character(len=1024), allocatable :: external_current_string(:,:)
179 real(real64), allocatable :: external_current_amplitude(:,:,:)
180 type(tdf_t), allocatable :: external_current_td_function(:)
181 type(tdf_t), allocatable :: external_current_td_phase(:)
182 real(real64), allocatable :: external_current_omega(:)
183 real(real64), allocatable :: external_current_phase(:)
184
186 character(len=1024), allocatable :: user_def_states(:,:,:)
187 logical :: fromscratch = .true.
188 type(mpi_grp_t) :: mpi_grp
189 type(mpi_grp_t) :: dom_st_mpi_grp
190 type(mpi_grp_t) :: system_grp
192#ifdef HAVE_SCALAPACK
193 type(blacs_proc_grid_t) :: dom_st_proc_grid
194#endif
195 type(distributed_t) :: dist
196 logical :: scalapack_compatible
197 integer :: lnst
198 integer :: st_start, st_end
199 integer, allocatable :: node(:)
201 type(poisson_t) :: poisson
202 integer :: transverse_field_mode
204 end type states_mxll_t
206 integer, public, parameter :: &
210contains
211
212 ! ---------------------------------------------------------
213 subroutine states_mxll_init(st, namespace, space)
214 type(states_mxll_t), target, intent(inout) :: st
215 type(namespace_t), intent(in) :: namespace
216 class(space_t), intent(in) :: space
218 type(block_t) :: blk
219 integer :: idim, nlines, ncols, il
220 real(real64), allocatable :: pos(:)
221 integer :: ix_max, iy_max, iz_max
225 call profiling_in('STATES_MXLL_INIT')
226
227 st%fromScratch = .true. ! this will be reset if restart_read is called
229 assert(space%dim == 3)
230 st%dim = space%dim
231 st%nst = 1
233 safe_allocate(st%user_def_e_field(1:st%dim))
234 safe_allocate(st%user_def_b_field(1:st%dim))
236 st%st_start = 1
237 st%st_end = st%nst
238 st%lnst = st%nst
239
240 safe_allocate(st%node(1:st%nst))
241 st%node(1:st%nst) = 0
243 call mpi_grp_init(st%mpi_grp, mpi_comm_undefined)
244 st%parallel_in_states = .false.
245 st%packed = .false.
247 ! The variable StatesPack is documented in states_elec.F90.
248 ! We cannot include the documentation twice.
249 ! TODO: We should think whether these variables could be moved to a higher (abstract) class.
251 call parse_variable(namespace, 'StatesPack', .true., st%pack_states)
253 call messages_print_var_value('StatesPack', st%pack_states, namespace=namespace)
255 !rs_sign is not defined any more by the user, since it does not influence
256 !the results of the simulations.
257 st%rs_sign = 1
259 !%Variable MaxwellFieldsCoordinate
260 !%Type block
261 !%Section Maxwell::Output
262 !%Description
263 !% The Maxwell MaxwellFieldsCoordinate block allows to output Maxwell fields at particular
264 !% points in space. For each point a new line with three columns has to be added to the block,
265 !% where the columns denote the x, y, and z coordinate of the point.
266 !%
267 !% <tt>%MaxwellFieldsCoordinate
268 !% <br>&nbsp;&nbsp; -1.0 | 2.0 | 4.0
269 !% <br>&nbsp;&nbsp; 0.0 | 1.0 | -2.0
270 !% <br>%</tt>
271 !%
272 !%End
274 safe_allocate(pos(1:st%dim))
275 st%selected_points_number = 1
276 if (parse_block(namespace, 'MaxwellFieldsCoordinate', blk) == 0) then
277 nlines = parse_block_n(blk)
278 st%selected_points_number = nlines
279 safe_allocate(st%selected_points_coordinate(1:st%dim,1:nlines))
280 safe_allocate(st%selected_points_rs_state(1:st%dim,1:nlines))
281 safe_allocate(st%selected_points_rs_state_long(1:st%dim,1:nlines))
282 safe_allocate(st%selected_points_rs_state_trans(1:st%dim,1:nlines))
283 safe_allocate(st%selected_points_map(1:nlines))
284 do il = 1, nlines
285 ncols = parse_block_cols(blk,0)
286 if (ncols < 3 .or. ncols > 3) then
287 message(1) = 'MaxwellFieldCoordinate must have 3 columns.'
288 call messages_fatal(1, namespace=namespace)
289 end if
290 do idim = 1, st%dim
291 call parse_block_float(blk, il-1, idim-1, pos(idim), units_inp%length)
292 end do
293 st%selected_points_coordinate(:,il) = pos
294 st%selected_points_rs_state(:,il) = m_z0
295 st%selected_points_rs_state_long(:,il) = m_z0
296 st%selected_points_rs_state_trans(:,il) = m_z0
297 end do
298 call parse_block_end(blk)
299 else
300 safe_allocate(st%selected_points_coordinate(1:st%dim, 1))
301 safe_allocate(st%selected_points_rs_state(1:st%dim, 1))
302 safe_allocate(st%selected_points_rs_state_long(1:st%dim, 1))
303 safe_allocate(st%selected_points_rs_state_trans(1:st%dim, 1))
304 safe_allocate(st%selected_points_map(1))
305 st%selected_points_coordinate(:,:) = m_zero
306 st%selected_points_rs_state(:,:) = m_z0
307 st%selected_points_rs_state_long(:,:) = m_z0
308 st%selected_points_rs_state_trans(:,:) = m_z0
309 st%selected_points_map(:) = -1
310 end if
311
312 safe_deallocate_a(pos)
313
314 st%surface_grid_rows_number(1) = 3
315 ix_max = st%surface_grid_rows_number(1)
316 st%surface_grid_rows_number(2) = 3
317 iy_max = st%surface_grid_rows_number(2)
318 st%surface_grid_rows_number(3) = 3
319 iz_max = st%surface_grid_rows_number(3)
320
321 safe_allocate(st%surface_grid_center(1:2, 1:st%dim, 1:ix_max, 1:iy_max))
322 safe_allocate(st%surface_grid_points_number(1:st%dim, 1:ix_max, 1:iy_max))
323
324 !%Variable TransverseFieldCalculation
325 !%Type integer
326 !%Default no
327 !%Section Maxwell
328 !%Description
329 !% This variable selects the method for the calculation of the transverse field.
330 !%Option helmholtz 1
331 !% Transverse field calculated from Helmholtz decompisition (unreliable at the moment).
332 !%Option total_minus_long 2
333 !% Total field minus longitudinal field.
334 !%End
335 call parse_variable(namespace, 'TransverseFieldCalculation', transverse_from_helmholtz, &
336 st%transverse_field_mode)
337
338 call profiling_out('STATES_MXLL_INIT')
339
340 pop_sub(states_mxll_init)
341
342 end subroutine states_mxll_init
343
344 ! ---------------------------------------------------------
346 subroutine states_mxll_allocate(st, mesh)
347 type(states_mxll_t), intent(inout) :: st
348 class(mesh_t), intent(in) :: mesh
349
350
351 push_sub(states_mxll_allocate)
352
353 call profiling_in('STATES_MXLL_ALLOCATE')
354
355 safe_allocate(st%rs_state(1:mesh%np_part, 1:st%dim))
356 st%rs_state(:,:) = m_z0
357
358 safe_allocate(st%rs_state_prev(1:mesh%np_part, 1:st%dim))
359 st%rs_state_prev(:,:) = m_z0
360
361 safe_allocate(st%rs_state_trans(1:mesh%np_part, 1:st%dim))
362 st%rs_state_trans(:,:) = m_z0
363
364 safe_allocate(st%rs_state_long(1:mesh%np_part, 1:st%dim))
365 st%rs_state_long(:,:) = m_z0
366
367 safe_allocate(st%rs_state_plane_waves(1:mesh%np_part, 1:st%dim))
368 st%rs_state_plane_waves(:,:) = m_z0
369
370 safe_allocate(st%rs_current_density_t1(1:mesh%np, 1:st%dim))
371 st%rs_current_density_t1 = m_z0
372
373 safe_allocate(st%rs_current_density_t2(1:mesh%np, 1:st%dim))
374 st%rs_current_density_t2 = m_z0
375
376 safe_allocate(st%rs_current_density_restart_t1(1:mesh%np_part, 1:st%dim))
377 st%rs_current_density_restart_t1 = m_z0
378
379 safe_allocate(st%rs_current_density_restart_t2(1:mesh%np_part, 1:st%dim))
380 st%rs_current_density_restart_t2 = m_z0
381
382 safe_allocate(st%ep(1:mesh%np_part))
383 safe_allocate(st%mu(1:mesh%np_part))
384 st%ep = p_ep
385 st%mu = p_mu
386
387 call profiling_out('STATES_MXLL_ALLOCATE')
388
389 pop_sub(states_mxll_allocate)
390 end subroutine states_mxll_allocate
391
392 ! ---------------------------------------------------------
393 subroutine states_mxll_end(st)
394 type(states_mxll_t), intent(inout) :: st
395
396
397 push_sub(states_mxll_end)
398
399 call profiling_in('STATES_MXLL_END')
400
401 safe_deallocate_a(st%rs_state)
402 safe_deallocate_a(st%rs_state_prev)
403 safe_deallocate_a(st%rs_state_trans)
404 safe_deallocate_a(st%selected_points_coordinate)
405 safe_deallocate_a(st%selected_points_rs_state)
406 safe_deallocate_a(st%selected_points_rs_state_long)
407 safe_deallocate_a(st%selected_points_rs_state_trans)
408 safe_deallocate_a(st%rs_current_density_t1)
409 safe_deallocate_a(st%rs_current_density_t2)
410 safe_deallocate_a(st%rs_state_long)
411 safe_deallocate_a(st%rs_current_density_restart_t1)
412 safe_deallocate_a(st%rs_current_density_restart_t2)
413 safe_deallocate_a(st%user_def_e_field)
414 safe_deallocate_a(st%user_def_b_field)
415
416 safe_deallocate_a(st%rs_state_const)
417 safe_deallocate_a(st%rs_state_const_td_function)
418 safe_deallocate_a(st%rs_state_const_amp)
419 safe_deallocate_a(st%rs_state_plane_waves)
420
421 safe_deallocate_a(st%surface_grid_center)
422 safe_deallocate_a(st%surface_grid_points_number)
423 safe_deallocate_a(st%surface_grid_points_map)
424 safe_deallocate_a(st%inner_points_map)
425 safe_deallocate_a(st%boundary_points_map)
426 safe_deallocate_a(st%inner_points_mask)
427 safe_deallocate_a(st%boundary_points_mask)
428 safe_deallocate_a(st%ep)
429 safe_deallocate_a(st%mu)
430 if (accel_is_enabled()) then
431 call accel_free_buffer(st%buff_inner_points_map)
432 call accel_free_buffer(st%buff_boundary_points_map)
433 call accel_free_buffer(st%buff_selected_points_map)
434 end if
435#ifdef HAVE_SCALAPACK
436 call blacs_proc_grid_end(st%dom_st_proc_grid)
437#endif
438 safe_deallocate_a(st%external_current_modus)
439 safe_deallocate_a(st%external_current_string)
440 safe_deallocate_a(st%external_current_amplitude)
441 safe_deallocate_a(st%external_current_td_function)
442 safe_deallocate_a(st%external_current_omega)
443 safe_deallocate_a(st%external_current_td_phase)
444
445 call distributed_end(st%dist)
446 safe_deallocate_a(st%node)
447
448 call profiling_out('STATES_MXLL_END')
449
450 pop_sub(states_mxll_end)
451 end subroutine states_mxll_end
452
453
454 !----------------------------------------------------------
455 subroutine build_rs_vector(e_vector, b_vector, rs_sign, rs_vector, ep_element, mu_element)
456 real(real64), intent(in) :: e_vector(:), b_vector(:)
457 complex(real64), intent(inout) :: rs_vector(:)
458 integer, intent(in) :: rs_sign
459 real(real64), optional, intent(in) :: ep_element
460 real(real64), optional, intent(in) :: mu_element
461
462 ! no PUSH_SUB, called too often
463
464 if (present(ep_element) .and. present(mu_element)) then
465 rs_vector = sqrt(ep_element/m_two) * e_vector + m_zi * rs_sign * sqrt(m_one/(m_two*mu_element)) * b_vector
466 else
467 rs_vector = sqrt(p_ep/m_two) * e_vector + m_zi * rs_sign * sqrt(m_one/(m_two*p_mu)) * b_vector
468 end if
469
470 end subroutine build_rs_vector
471
472
473 !----------------------------------------------------------
474 subroutine build_rs_state(e_field, b_field, rs_sign, rs_state, mesh, ep_field, mu_field, np)
475 real(real64), intent(in) :: e_field(:,:), b_field(:,:)
476 complex(real64), intent(inout) :: rs_state(:,:)
477 integer, intent(in) :: rs_sign
478 class(mesh_t), intent(in) :: mesh
479 real(real64), optional, intent(in) :: ep_field(:)
480 real(real64), optional, intent(in) :: mu_field(:)
481 integer, optional, intent(in) :: np
482
483 integer :: ip, np_
484
485 push_sub(build_rs_state)
486
487 call profiling_in('BUILD_RS_STATE')
489 np_ = optional_default(np, mesh%np)
490
491 do ip = 1, np_
492 if (present(ep_field) .and. present(mu_field)) then
493 rs_state(ip, :) = sqrt(ep_field(ip)/m_two) * e_field(ip, :) &
494 + m_zi * rs_sign * sqrt(m_one/(m_two*mu_field(ip))) * b_field(ip, :)
495 else
496 rs_state(ip, :) = sqrt(p_ep/m_two) * e_field(ip, :) &
497 + m_zi * rs_sign * sqrt(m_one/(m_two*p_mu)) * b_field(ip, :)
498 end if
499 end do
500
501 call profiling_out('BUILD_RS_STATE')
502
503 pop_sub(build_rs_state)
504
505 end subroutine build_rs_state
506
507
508 !----------------------------------------------------------
509 subroutine build_rs_current_state(current_state, mesh, rs_current_state, ep_field, np)
510 real(real64), intent(in) :: current_state(:,:)
511 class(mesh_t), intent(in) :: mesh
512 complex(real64), intent(inout) :: rs_current_state(:,:)
513 real(real64), optional, intent(in) :: ep_field(:)
514 integer, optional, intent(in) :: np
515
516 integer :: ip, idim, np_, ff_dim
517
518 ! no PUSH_SUB, called too often
519
520 call profiling_in("BUILD_RS_CURRENT_STATE")
521
522 np_ = optional_default(np, mesh%np)
523 ff_dim = size(current_state, dim=2)
524
525 if (present(ep_field)) then
526 do idim = 1, ff_dim
527 do ip = 1, np_
528 rs_current_state(ip, idim) = m_one/sqrt(m_two*ep_field(ip)) * current_state(ip, idim)
529 end do
530 end do
531 else
532 do idim = 1, ff_dim
533 do ip = 1, np_
534 rs_current_state(ip, idim) = m_one/sqrt(m_two*p_ep) * current_state(ip, idim)
535 end do
536 end do
537 end if
538
539 call profiling_out("BUILD_RS_CURRENT_STATE")
540
541 end subroutine build_rs_current_state
542
543
544 !----------------------------------------------------------
545 subroutine get_electric_field_vector(rs_state_vector, electric_field_vector, ep_element)
546 complex(real64), intent(in) :: rs_state_vector(:)
547 real(real64), intent(out) :: electric_field_vector(:)
548 real(real64), optional, intent(in) :: ep_element
549
550 ! no PUSH_SUB, called too often
551
552 if (present(ep_element)) then
553 electric_field_vector(:) = sqrt(m_two/ep_element) * real(rs_state_vector(:), real64)
554 else
555 electric_field_vector(:) = sqrt(m_two/p_ep) * real(rs_state_vector(:), real64)
556 end if
557
558 end subroutine get_electric_field_vector
559
560
561 !----------------------------------------------------------
562 subroutine get_magnetic_field_vector(rs_state_vector, rs_sign, magnetic_field_vector, mu_element)
563 complex(real64), intent(in) :: rs_state_vector(:)
564 integer, intent(in) :: rs_sign
565 real(real64), intent(out) :: magnetic_field_vector(:)
566 real(real64), optional, intent(in) :: mu_element
567
568 ! no PUSH_SUB, called too often
570 if (present(mu_element)) then
571 magnetic_field_vector(:) = sqrt(m_two*mu_element) * rs_sign * aimag(rs_state_vector(:))
572 else
573 magnetic_field_vector(:) = sqrt(m_two*p_mu) * rs_sign * aimag(rs_state_vector(:))
574 end if
575
576 end subroutine get_magnetic_field_vector
577
578
579 !----------------------------------------------------------
580 subroutine get_electric_field_state(rs_state, mesh, electric_field, ep_field, np)
581 complex(real64), intent(in) :: rs_state(:,:)
582 class(mesh_t), intent(in) :: mesh
583 real(real64), intent(out) :: electric_field(:,:)
584 real(real64), optional, intent(in) :: ep_field(:)
585 integer, optional, intent(in) :: np
586
587 integer :: ip, np_
588
590
591 call profiling_in('GET_ELECTRIC_FIELD_STATE')
592
593 np_ = optional_default(np, mesh%np)
594
595 do ip = 1, np_
596 if (present(ep_field)) then
597 electric_field(ip, :) = sqrt(m_two/ep_field(ip)) * real(rs_state(ip, :), real64)
598 else
599 electric_field(ip,:) = sqrt(m_two/p_ep) * real(rs_state(ip, :), real64)
600 end if
601 end do
602
603 call profiling_out('GET_ELECTRIC_FIELD_STATE')
606
607 end subroutine get_electric_field_state
608
609
610 !----------------------------------------------------------
611 subroutine get_magnetic_field_state(rs_state, mesh, rs_sign, magnetic_field, mu_field, np)
612 complex(real64), intent(in) :: rs_state(:,:)
613 class(mesh_t), intent(in) :: mesh
614 integer, intent(in) :: rs_sign
615 real(real64), intent(out) :: magnetic_field(:,:)
616 real(real64), optional, intent(in) :: mu_field(:)
617 integer, optional, intent(in) :: np
618
619 integer :: ip, np_
620
622
623 call profiling_in('GET_MAGNETIC_FIELD_STATE')
624
625 np_ = optional_default(np, mesh%np)
626
627 if (present(mu_field)) then
628 do ip = 1, np_
629 magnetic_field(ip, :) = sqrt(m_two*mu_field(ip)) * rs_sign * aimag(rs_state(ip, :))
630 end do
631 else
632 do ip = 1, np_
633 magnetic_field(ip, :) = sqrt(m_two*p_mu) * rs_sign * aimag(rs_state(ip, :))
634 end do
635 end if
636
637 call profiling_out('GET_MAGNETIC_FIELD_STATE')
638
641 end subroutine get_magnetic_field_state
642
643 !----------------------------------------------------------
644 subroutine get_rs_state_at_point(rs_state_point, rs_state, pos, st, mesh)
645
646 complex(real64), intent(inout) :: rs_state_point(:,:)
647 complex(real64), intent(in) :: rs_state(:,:)
648 real(real64), intent(in) :: pos(:,:)
649 type(states_mxll_t), intent(in) :: st
650 class(mesh_t), intent(in) :: mesh
651
652 integer :: ip, pos_index, rankmin
653 real(real64) :: dmin
654 complex(real64), allocatable :: ztmp(:)
655
656 push_sub(get_rs_state_at_point)
658 safe_allocate(ztmp(1:size(rs_state, dim=2)))
659
660 do ip = 1, st%selected_points_number
661 pos_index = mesh_nearest_point(mesh, pos(:,ip), dmin, rankmin)
662 if (mesh%mpi_grp%rank == rankmin) then
663 ztmp(:) = rs_state(pos_index, :)
664 end if
665 if (mesh%parallel_in_domains) then
666 call mesh%mpi_grp%bcast(ztmp, st%dim, mpi_double_complex, rankmin)
667 end if
668 rs_state_point(:, ip) = ztmp(:)
669 end do
670
671 safe_deallocate_a(ztmp)
672
673
674 pop_sub(get_rs_state_at_point)
676
677
678 !----------------------------------------------------------
679 subroutine get_rs_state_batch_selected_points(rs_state_point, rs_stateb, st, mesh)
680 complex(real64), contiguous, intent(inout) :: rs_state_point(:,:)
681 type(batch_t), intent(in) :: rs_stateb
682 type(states_mxll_t), intent(in) :: st
683 class(mesh_t), intent(in) :: mesh
684
685 integer :: ip_in, ip
686 complex(real64) :: rs_state_tmp(1:st%dim, 1:st%selected_points_number)
687 type(accel_kernel_t), save :: kernel
688 type(accel_mem_t) :: buff_points
689 integer(int64), dimension(3) :: gsizes, bsizes
690
692
693 rs_state_tmp(:,:) = m_z0
694
695 select case (rs_stateb%status())
696 case (batch_not_packed)
697 do ip_in = 1, st%selected_points_number
698 ip = st%selected_points_map(ip_in)
699 if (ip >= 0) then
700 rs_state_tmp(1:st%dim, ip_in) = rs_stateb%zff_linear(ip, 1:st%dim)
701 end if
702 end do
703 case (batch_packed)
704 do ip_in = 1, st%selected_points_number
705 ip = st%selected_points_map(ip_in)
706 if (ip >= 0) then
707 rs_state_tmp(1:st%dim, ip_in) = rs_stateb%zff_pack(1:st%dim, ip)
708 end if
709 end do
711 call accel_kernel_start_call(kernel, 'get_points.cu', 'get_selected_points')
712
714 st%selected_points_number*st%dim)
715 call accel_set_buffer_to_zero(buff_points, type_integer, st%selected_points_number*st%dim)
716
717 call accel_set_kernel_arg(kernel, 0, st%selected_points_number)
718 call accel_set_kernel_arg(kernel, 1, st%buff_selected_points_map)
719 call accel_set_kernel_arg(kernel, 2, rs_stateb%ff_device)
720 call accel_set_kernel_arg(kernel, 3, log2(int(rs_stateb%pack_size_real(1), int32)))
721 call accel_set_kernel_arg(kernel, 4, buff_points)
722 call accel_set_kernel_arg(kernel, 5, st%dim*2)
723
724 ! Compute the grid (extend to another dimensions if the size of the problem is too big)
725 call accel_grid_size_extend_dim(int(st%selected_points_number, int64), rs_stateb%pack_size_real(1), &
726 gsizes, bsizes, kernel)
727
728 call accel_kernel_run(kernel, gsizes, bsizes)
729
730 call accel_read_buffer(buff_points, st%dim, st%selected_points_number, rs_state_tmp)
731 call accel_free_buffer(buff_points)
732 end select
733
734 call mesh%mpi_grp%allreduce(rs_state_tmp, rs_state_point, st%selected_points_number*st%dim, mpi_double_complex, mpi_sum)
735
738
739 !----------------------------------------------------------
740 subroutine get_divergence_field(gr, field, field_div, charge_density)
741 type(grid_t), intent(in) :: gr
742 real(real64), contiguous, intent(inout) :: field(:,:)
743 real(real64), contiguous, intent(inout) :: field_div(:)
744 logical, intent(in) :: charge_density
745
746 push_sub(get_divergence_field)
747
748 call dderivatives_div(gr%der, field, field_div)
749
750 if (optional_default(charge_density,.false.)) then
751 field_div = p_ep * field_div
752 end if
753
754 pop_sub(get_divergence_field)
755 end subroutine get_divergence_field
756
757
758 ! ---------------------------------------------------------
759 subroutine get_poynting_vector(mesh, st, rs_state, rs_sign, poynting_vector, ep_field, mu_field)
760 class(mesh_t), intent(in) :: mesh
761 type(states_mxll_t), intent(in) :: st
762 complex(real64), intent(in) :: rs_state(:,:)
763 integer, intent(in) :: rs_sign
764 real(real64), intent(out) :: poynting_vector(:,:)
765 real(real64), optional, intent(in) :: ep_field(:)
766 real(real64), optional, intent(in) :: mu_field(:)
767
768 integer :: ip
769
770 push_sub(get_poynting_vector)
771 if (present(ep_field) .and. present(mu_field)) then
772 do ip = 1, mesh%np
773 poynting_vector(ip, 1:3) = m_one/mu_field(ip) * sqrt(m_two/ep_field(ip)) &
774 * sqrt(m_two*mu_field(ip)) &
775 * dcross_product(real(rs_state(ip, 1:3), real64) , &
776 rs_sign*aimag(rs_state(ip,1:3)))
777 end do
778 else
779 do ip = 1, mesh%np
780 poynting_vector(ip, 1:3) = m_one/st%mu(ip) * sqrt(m_two/st%ep(ip)) &
781 * sqrt(m_two*st%mu(ip)) &
782 * dcross_product(real(rs_state(ip, 1:3), real64) , &
783 rs_sign*aimag(rs_state(ip, 1:3)))
784 end do
785 end if
786
787 pop_sub(get_poynting_vector)
788 end subroutine get_poynting_vector
789
790
791 ! ---------------------------------------------------------
792 subroutine get_orbital_angular_momentum(mesh, poynting_vector, orbital_angular_momentum)
793 type(mesh_t), intent(in) :: mesh
794 real(real64), intent(in) :: poynting_vector(:,:)
795 real(real64), intent(out) :: orbital_angular_momentum(:,:)
796
797 integer :: ip
798
800
801 do ip = 1, mesh%np
802 orbital_angular_momentum(ip,1:3) = dcross_product(real(mesh%x(1:3, ip), real64) , &
803 poynting_vector(ip, 1:3))
804 end do
805
807 end subroutine get_orbital_angular_momentum
808
809 ! ---------------------------------------------------------
810 subroutine mxll_set_batch(rs_stateb, rs_state, np, dim, offset)
811 type(batch_t), intent(inout) :: rs_stateb
812 complex(real64), contiguous, intent(in) :: rs_state(:, :)
813 integer, intent(in) :: np
814 integer, intent(in) :: dim
815 integer, optional, intent(in) :: offset
816
817 integer :: offset_, idir
818
819 push_sub(mxll_set_batch)
820
821 offset_ = optional_default(offset, 1)
822
823 do idir = offset_, offset_ + dim - 1
824 call batch_set_state(rs_stateb, idir, np, rs_state(:, idir))
825 end do
826
827 pop_sub(mxll_set_batch)
828 end subroutine mxll_set_batch
829
830 ! ---------------------------------------------------------
831 subroutine mxll_get_batch(rs_stateb, rs_state, np, dim, offset)
832 type(batch_t), intent(in) :: rs_stateb
833 complex(real64), contiguous, intent(out) :: rs_state(:, :)
834 integer, intent(in) :: np
835 integer, intent(in) :: dim
836 integer, optional, intent(in) :: offset
837
838 integer :: offset_, idir
839
840 push_sub(mxll_get_batch)
841
842 offset_ = optional_default(offset, 1)
843
844 do idir = offset_, offset_ + dim - 1
845 call batch_get_state(rs_stateb, idir, np, rs_state(:, idir))
846 end do
847
848 pop_sub(mxll_get_batch)
849 end subroutine mxll_get_batch
850
851 !----------------------------------------------------------
852 subroutine get_transverse_rs_state(helmholtz, st, namespace)
853 type(helmholtz_decomposition_t), intent(inout) :: helmholtz
854 type(states_mxll_t), intent(inout) :: st
855 type(namespace_t), intent(in) :: namespace
856
857
859
860 call profiling_in('GET_TRANSVERSE_RS_STATE')
861
862 select case (st%transverse_field_mode)
864 call helmholtz%get_trans_field(namespace, st%rs_state_trans, total_field=st%rs_state)
866 call helmholtz%get_long_field(namespace, st%rs_state_long, total_field=st%rs_state)
867 st%rs_state_trans = st%rs_state - st%rs_state_long
868 case default
869 message(1) = 'Unknown transverse field calculation mode.'
870 call messages_fatal(1, namespace=namespace)
871 end select
872
873 call profiling_out('GET_TRANSVERSE_RS_STATE')
874
876
877 end subroutine get_transverse_rs_state
878
879
880end module states_mxll_oct_m
881
882
883
884!! Local Variables:
885!! mode: f90
886!! coding: utf-8
887!! End:
There are several ways how to call batch_set_state and batch_get_state:
Definition: batch_ops.F90:203
Prints out to iunit a message in the form: ["InputVariable" = value] where "InputVariable" is given b...
Definition: messages.F90:181
subroutine, public accel_free_buffer(this, async)
Definition: accel.F90:1005
subroutine, public accel_kernel_start_call(this, file_name, kernel_name, flags)
Definition: accel.F90:1413
integer, parameter, public accel_mem_read_write
Definition: accel.F90:185
pure logical function, public accel_is_enabled()
Definition: accel.F90:402
This module implements batches of mesh functions.
Definition: batch.F90:135
integer, parameter, public batch_not_packed
functions are stored in CPU memory, unpacked order
Definition: batch.F90:286
integer, parameter, public batch_device_packed
functions are stored in device memory in packed order
Definition: batch.F90:286
integer, parameter, public batch_packed
functions are stored in CPU memory, in transposed (packed) order
Definition: batch.F90:286
This module implements common operations on batches of mesh functions.
Definition: batch_ops.F90:118
This module provides the BLACS processor grid.
subroutine, public blacs_proc_grid_end(this)
This module calculates the derivatives (gradients, Laplacians, etc.) of a function.
subroutine, public dderivatives_div(der, ff, op_ff, ghost_update, set_bc, to_cartesian)
apply the divergence operator to a vector of mesh functions
subroutine, public distributed_end(this)
real(real64), parameter, public m_two
Definition: global.F90:193
real(real64), parameter, public m_zero
Definition: global.F90:191
real(real64), parameter, public p_mu
Definition: global.F90:238
real(real64), parameter, public p_ep
Definition: global.F90:237
complex(real64), parameter, public m_z0
Definition: global.F90:201
complex(real64), parameter, public m_zi
Definition: global.F90:205
real(real64), parameter, public m_one
Definition: global.F90:192
This module implements the underlying real-space grid.
Definition: grid.F90:119
The Helmholtz decomposition is intended to contain "only mathematical" functions and procedures to co...
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:117
pure real(real64) function, dimension(1:3), public dcross_product(a, b)
Definition: math.F90:1886
This module defines various routines, operating on mesh functions.
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:120
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
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:409
type(mpi_comm), parameter, public mpi_comm_undefined
used to indicate a communicator has not been initialized
Definition: mpi.F90:138
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
integer function, public parse_block(namespace, name, blk, check_varinfo_)
Definition: parser.F90:615
subroutine, public profiling_out(label)
Increment out counter and sum up difference between entry and exit time.
Definition: profiling.F90:631
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
Definition: profiling.F90:554
This module handles spin dimensions of the states and the k-point distribution.
This module handles groups of electronic batches and their parallel distribution.
subroutine, public mxll_set_batch(rs_stateb, rs_state, np, dim, offset)
subroutine, public get_orbital_angular_momentum(mesh, poynting_vector, orbital_angular_momentum)
integer, parameter, public transverse_from_helmholtz
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 build_rs_vector(e_vector, b_vector, rs_sign, rs_vector, ep_element, mu_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 build_rs_current_state(current_state, mesh, rs_current_state, ep_field, np)
integer, parameter, public transverse_as_total_minus_long
subroutine, public build_rs_state(e_field, b_field, rs_sign, rs_state, mesh, ep_field, mu_field, np)
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_poynting_vector(mesh, st, rs_state, rs_sign, poynting_vector, ep_field, mu_field)
subroutine, public get_divergence_field(gr, field, field_div, charge_density)
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)
type(type_t), public type_cmplx
Definition: types.F90:136
type(type_t), 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.
type(unit_system_t), public units_inp
the units systems for reading and writing
Class defining batches of mesh functions.
Definition: batch.F90:161
Distribution of N instances over mpi_grpsize processes, for the local rank mpi_grprank....
Description of the grid, containing information on derivatives, stencil, and symmetries.
Definition: grid.F90:171
Describes mesh distribution to nodes.
Definition: mesh.F90:187
This is defined even when running serial.
Definition: mpi.F90:144
int true(void)