65 integer,
parameter :: &
66 MEDIUM_PARALLELEPIPED = 1, &
72 type(space_t) :: space
73 real(real64) :: ep_factor
74 real(real64) :: mu_factor
75 real(real64) :: sigma_e_factor
76 real(real64) :: sigma_m_factor
77 integer :: edge_profile
78 type(output_t) :: outp
80 type(multicomm_t) :: mc
81 type(single_medium_box_t) :: medium_box
97 procedure linear_medium_constructor
108 class(linear_medium_t),
pointer :: sys
109 type(namespace_t),
intent(in) :: namespace
127 class(linear_medium_t),
target,
intent(inout) :: this
128 type(namespace_t),
intent(in) :: namespace
130 integer :: nlines, ncols,n_points, n_global_points
131 integer(int64) :: index_range(4)
132 integer,
allocatable :: tmp(:)
139 this%namespace = namespace
140 this%space =
space_t(this%namespace)
141 if (this%space%is_periodic())
then
147 index_range(1) = this%gr%np_global
150 index_range(4) = 100000
154 index_range, (/ 5000, 1, 1, 1 /))
174 if (
parse_block(namespace,
'LinearMediumProperties', blk) == 0)
then
177 if (nlines /= 1)
then
188 write(
message(1),
'(a,es9.2)')
'Box epsilon factor: ', this%ep_factor
189 write(
message(2),
'(a,es9.2)')
'Box mu factor: ', this%mu_factor
190 write(
message(3),
'(a,es9.2)')
'Box electric sigma: ', this%sigma_e_factor
191 write(
message(4),
'(a,es9.2)')
'Box magnetic sigma: ', this%sigma_m_factor
198 message(1) =
'You must specify the properties of your linear medium through the LinearMediumProperties block.'
214 call parse_variable(namespace,
'LinearMediumEdgeProfile', option__linearmediumedgeprofile__edged, this%edge_profile)
215 if (this%edge_profile == option__linearmediumedgeprofile__edged)
then
216 write(
message(1),
'(a,a)')
'Box shape: ',
'edged'
217 else if (this%edge_profile == option__linearmediumedgeprofile__smooth)
then
218 write(
message(1),
'(a,a)')
'Box shape: ',
'smooth'
222 allocate(this%supported_interactions(0))
224 call this%quantities%add(
quantity_t(
"permitivity", updated_on_demand = .false., always_available = .
true.))
225 call this%quantities%add(
quantity_t(
"permeability", updated_on_demand = .false., always_available = .
true.))
226 call this%quantities%add(
quantity_t(
"E conductivity", updated_on_demand = .false., always_available = .
true.))
227 call this%quantities%add(
quantity_t(
"M conductivity", updated_on_demand = .false., always_available = .
true.))
234 this%outp%iter_dir, this%medium_box%points_map, this%medium_box%ep, &
235 this%medium_box%mu, this%medium_box%c)
237 if (this%medium_box%check_medium_points)
then
238 safe_allocate(tmp(1:this%gr%np))
240 write(
message(1),
'(a, a, a)')
'Check of points inside surface of medium ', trim(this%medium_box%filename),
":"
243 write(
message(1),
'(a, I8)')
'Number of points inside medium (normal coordinates):', &
244 this%medium_box%global_points_number
245 write(
message(2),
'(a, I8)')
'Number of points inside medium (rescaled coordinates):', n_global_points
248 safe_deallocate_a(tmp)
263 select type (interaction)
265 message(1) =
"Trying to initialize an unsupported interaction by a linear medium."
282 select type (interaction)
285 partner%gr, partner%space, partner%namespace)
288 medium_box_grid => partner%medium_box%to_grid(partner%gr)
291 call regridding%do_transfer(interaction%medium_box%aux_ep, medium_box_grid%aux_ep)
292 call regridding%do_transfer(interaction%medium_box%aux_mu, medium_box_grid%aux_mu)
293 call regridding%do_transfer(interaction%medium_box%ep, medium_box_grid%ep)
294 call regridding%do_transfer(interaction%medium_box%mu, medium_box_grid%mu)
295 call regridding%do_transfer(interaction%medium_box%c, medium_box_grid%c)
296 call regridding%do_transfer(interaction%medium_box%sigma_e, medium_box_grid%sigma_e)
297 call regridding%do_transfer(interaction%medium_box%sigma_m, medium_box_grid%sigma_m)
299 safe_deallocate_p(regridding)
300 safe_deallocate_p(medium_box_grid)
302 message(1) =
"Trying to initialize an unsupported interaction by a linear medium."
322 character(len=:),
allocatable,
intent(out) :: updated_quantities(:)
335 real(real64),
intent(in) :: tol
353 select type (interaction)
357 message(1) =
"Unsupported interaction."
368 integer :: restart_file_unit
372 message(1) =
"Linear medium system "//trim(this%namespace%get())//
" will only write a dummy restart file."
376 restart_file_unit =
io_open(
'restart/'//
td_dir//
'restart_linear_medium', this%namespace, action=
'write')
377 write(restart_file_unit, *)
'Linear medium restart file'
380 message(1) =
"Successfully wrote restart data for system "//trim(this%namespace%get())
391 integer :: restart_file_unit
396 restart_file_unit =
io_open(
'restart/'//
td_dir//
'restart_linear_medium', this%namespace, action=
'read', die=.false.)
397 if (restart_file_unit /= -1)
then
406 message(1) =
"Successfully read restart data for system "//trim(this%namespace%get())
418 this%kinetic_energy =
m_zero
444 type(
grid_t),
intent(in) :: gr
446 integer :: ip, ip_in, ip_bd, idim
447 integer,
allocatable :: tmp_points_map(:), tmp_bdry_map(:)
448 real(real64) :: bounds(2,gr%box%dim), xx(gr%box%dim)
455 safe_allocate(tmp_points_map(1:gr%np))
456 safe_allocate(tmp_bdry_map(1:gr%np))
463 medium_box%global_points_number, tmp_points_map)
468 bounds(1,idim) = medium_box%center(idim) - medium_box%lsize(idim)/
m_two
469 bounds(2,idim) = medium_box%center(idim) + medium_box%lsize(idim)/
m_two
474 xx(1:3) = gr%x(ip, 1:3)
478 tmp_points_map(ip_in) = ip
482 tmp_bdry_map(ip_bd) = ip
486 medium_box%points_number = ip_in
487 medium_box%bdry_number = ip_bd
489 safe_allocate(medium_box%bdry_map(1:medium_box%bdry_number))
490 medium_box%bdry_map = 0
491 medium_box%bdry_map = tmp_bdry_map(1:medium_box%bdry_number)
495 medium_box%points_map(:) = tmp_points_map(1:medium_box%points_number)
498 safe_deallocate_a(tmp_points_map)
499 safe_deallocate_a(tmp_bdry_map)
506 logical pure function check_point_in_bounds(xx, bounds) result (check)
507 real(real64),
intent(in) :: xx(:)
508 real(real64),
intent(in) :: bounds(:,:)
511 if ((xx(1) >= bounds(1,1)) .and. (xx(1) <= bounds(2,1)) .and. &
512 (xx(2) >= bounds(1,2)) .and. (xx(2) <= bounds(2,2)) .and. &
513 (xx(3) >= bounds(1,3)) .and. (xx(3) <= bounds(2,3)) )
then
519 logical pure function check_point_on_bounds(xx, bounds) result (check)
520 real(real64),
intent(in) :: xx(:)
521 real(real64),
intent(in) :: bounds(:,:)
524 if (
is_close(xx(1), bounds(1,1)) .and. (xx(2) >= bounds(1,2) .and. xx(3) >= bounds(1,3)) &
525 .and. (xx(2) <= bounds(2,2) .and. xx(3) <= bounds(2,3)) .or. &
526 is_close(xx(2), bounds(1,2)) .and. (xx(1) >= bounds(1,1) .and. xx(3) >= bounds(1,3)) &
527 .and. (xx(1) <= bounds(2,1) .and. xx(3) <= bounds(2,3)) .or. &
528 is_close(xx(3), bounds(1,3)) .and. (xx(1) >= bounds(1,1) .and. xx(2) >= bounds(1,2)) &
529 .and. (xx(1) <= bounds(2,1) .and. xx(2) <= bounds(2,2)) .or. &
530 is_close(xx(1), bounds(2,1)) .and. (xx(2) >= bounds(1,2) .and. xx(3) >= bounds(1,3)) &
531 .and. (xx(2) <= bounds(2,2) .and. xx(3) <= bounds(2,3)) .or. &
532 is_close(xx(2), bounds(2,2)) .and. (xx(1) >= bounds(1,1) .and. xx(3) >= bounds(1,3)) &
533 .and. (xx(1) <= bounds(2,1) .and. xx(3) <= bounds(2,3)) .or. &
534 is_close(xx(3), bounds(2,3)) .and. (xx(1) >= bounds(1,1) .and. xx(2) >= bounds(1,2)) &
535 .and. (xx(1) <= bounds(2,1) .and. xx(2) <= bounds(2,2)) )
then
547 type(single_medium_box_t),
intent(inout) :: medium_box
548 type(grid_t),
intent(in) :: gr
550 integer :: ip, ip_in, ip_bd, ipp
551 real(real64) :: xx(gr%box%dim), xxp(gr%box%dim), dd, dd_max, dd_min
552 real(real64),
allocatable :: tmp(:), tmp_grad(:,:)
556 call profiling_in(
'GET_LIN_MEDIUM_EM_PROP')
558 safe_allocate(tmp(1:gr%np_part))
559 safe_allocate(tmp_grad(1:gr%np_part, 1:gr%box%dim))
560 dd_max = max(2*gr%spacing(1), 2*gr%spacing(2), 2*gr%spacing(3))
562 do ip_in = 1, medium_box%points_number
563 ip = medium_box%points_map(ip_in)
564 if (this%edge_profile == option__linearmediumedgeprofile__smooth)
then
565 assert(
allocated(medium_box%bdry_map))
567 xx(1:3) = gr%x(ip,1:3)
570 do ip_bd = 1, medium_box%bdry_number
571 ipp = medium_box%bdry_map(ip_bd)
572 xxp(1:3) = gr%x(ipp,1:3)
573 dd = norm2(xx(1:3) - xxp(1:3))
574 if (dd < dd_min) dd_min = dd
577 medium_box%ep(ip_in) = p_ep + ((p_ep * this%ep_factor - p_ep) &
578 * m_one/(m_one +
exp(-m_five/dd_max * (dd_min - m_two*dd_max))))
579 medium_box%mu(ip_in) = p_mu + ((p_mu * this%mu_factor - p_mu) &
580 * m_one/(m_one +
exp(-m_five/dd_max * (dd_min - m_two*dd_max))))
581 medium_box%c(ip_in) = m_one/
sqrt(medium_box%ep(ip_in)*medium_box%mu(ip_in))
582 medium_box%sigma_e(ip_in) = this%sigma_e_factor &
583 * m_one/(m_one +
exp(-m_five/dd_max * (dd_min - m_two*dd_max)))
584 medium_box%sigma_m(ip_in) = this%sigma_m_factor &
585 * m_one/(m_one +
exp(-m_five/dd_max * (dd_min - m_two*dd_max)))
587 else if (this%edge_profile == option__linearmediumedgeprofile__edged)
then
589 medium_box%ep(ip_in) = p_ep * this%ep_factor
590 medium_box%mu(ip_in) = p_mu * this%mu_factor
591 medium_box%c(ip_in) = m_one/
sqrt(medium_box%ep(ip_in)*medium_box%mu(ip_in))
592 medium_box%sigma_e(ip_in) = this%sigma_e_factor
593 medium_box%sigma_m(ip_in) = this%sigma_m_factor
599 do ip_in = 1, medium_box%points_number
600 ip = medium_box%points_map(ip_in)
601 tmp(ip)= medium_box%ep(ip_in)
603 call dderivatives_grad(gr%der, tmp, tmp_grad, set_bc = .false.)
604 do ip_in = 1, medium_box%points_number
605 ip = medium_box%points_map(ip_in)
606 medium_box%aux_ep(ip_in, :) = tmp_grad(ip, :)/(m_four * medium_box%ep(ip_in))
610 do ip_in = 1, medium_box%points_number
611 ip = medium_box%points_map(ip_in)
612 tmp(ip) = medium_box%mu(ip_in)
614 call dderivatives_grad(gr%der, tmp, tmp_grad, set_bc = .false.)
615 do ip_in = 1, medium_box%points_number
616 ip = medium_box%points_map(ip_in)
617 medium_box%aux_mu(ip_in, :) = tmp_grad(ip, :)/(m_four * medium_box%mu(ip_in))
622 safe_deallocate_a(tmp)
623 safe_deallocate_a(tmp_grad)
625 call profiling_out(
'GET_LIN_MEDIUM_EM_PROP')
633 character(len=256),
intent(in) :: filename
634 class(mesh_t),
intent(in) :: mesh
635 integer,
intent(out) :: n_points
636 integer,
intent(out) :: global_points_number
637 integer,
intent(inout) :: tmp_map(:)
638 real(real64),
optional,
intent(in) :: scale_factor
641 real(real64) :: xx(3)
642 type(cgal_polyhedra_t) :: cgal_poly
646 call profiling_in(
'GET_POINTS_MAP_FROM_FILE')
648 call cgal_polyhedron_init(cgal_poly, trim(filename), verbose = .false.)
652 if (
present(scale_factor))
then
653 xx(1:3) = scale_factor * mesh%x(ip, 1:3)
655 xx(1:3) = mesh%x(ip, 1:3)
657 if (cgal_polyhedron_point_inside(cgal_poly, xx(1), xx(2), xx(3)))
then
663 call cgal_polyhedron_end(cgal_poly)
665 call mpi_world%allreduce(ip_in, global_points_number, 1, mpi_integer, mpi_sum)
667 call profiling_out(
'GET_POINTS_MAP_FROM_FILE')
676 type(single_medium_box_t),
intent(inout) :: medium_box
677 type(namespace_t),
intent(in) :: namespace
679 integer :: nlines, ncols, idim
698 if (.not. varinfo_valid_option(
'LinearMediumBoxShape', medium_box%box_shape))
then
699 call messages_input_error(namespace,
'LinearMediumBoxShape')
702 select case (medium_box%box_shape)
718 if (parse_block(namespace,
'LinearMediumBoxSize', blk) == 0)
then
719 call messages_print_with_emphasis(msg=trim(
'Linear Medium box center and size:'), namespace=namespace)
721 nlines = parse_block_n(blk)
722 if (nlines /= 1)
then
723 call messages_input_error(namespace,
'LinearMediumBoxSize',
'should consist of one line')
725 ncols = parse_block_cols(blk, 0)
727 call messages_input_error(namespace,
'LinearMediumBoxSize',
'should consist of six columns')
730 call parse_block_float(blk, 0, idim-1, medium_box%center(idim))
731 call parse_block_float(blk, 0, idim+2, medium_box%lsize(idim))
733 write(message(1),
'(a,es9.2,a,es9.2,a,es9.2)')
'Box center: ', medium_box%center(1),
' | ',&
734 medium_box%center(2),
' | ', medium_box%center(3)
735 write(message(2),
'(a,es9.2,a,es9.2,a,es9.2)')
'Box size : ', medium_box%lsize(1),
' | ', &
736 medium_box%lsize(2),
' | ', medium_box%lsize(3)
737 write(message(3),
'(a)')
""
738 call messages_info(3)
739 call parse_block_end(blk)
741 call messages_print_with_emphasis(namespace=namespace)
743 message(1) =
"For parallelepiped box shapes, you must provide a LinearMediumBoxSize block."
744 call messages_fatal(1, namespace=namespace)
754 if (parse_is_defined(namespace,
'LinearMediumBoxFile'))
then
755 call parse_variable(namespace,
'LinearMediumBoxFile',
'mediumboxfile', medium_box%filename)
757 message(1) =
"When using box_file as the box shape, you must provide a filename through the LinearMediumBoxFile variable."
758 call messages_fatal(1, namespace=namespace)
770 call parse_variable(namespace,
'CheckPointsMediumFromFile', .false., medium_box%check_medium_points)
double exp(double __x) __attribute__((__nothrow__
double sqrt(double __x) __attribute__((__nothrow__
logical pure function check_point_on_bounds(xx, bounds)
logical pure function check_point_in_bounds(xx, bounds)
This module implements the basic elements defining algorithms.
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.
real(real64), parameter, public m_two
real(real64), parameter, public m_zero
character(len= *), parameter, public td_dir
This module implements the underlying real-space grid.
subroutine, public grid_init_stage_1(gr, namespace, space, symm, latt, n_sites, site_position)
First stage of the grid initialization.
subroutine, public grid_init_stage_2(gr, namespace, space, mc, qvector)
Second stage of the grid initialization.
subroutine, public grid_end(gr)
finalize a grid object
integer, parameter, public linear_medium_to_em_field
This module defines the abstract interaction_t class, and some auxiliary classes for interactions.
subroutine, public io_close(iunit, grp)
subroutine, public io_mkdir(fname, namespace, parents)
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
This module defines a linear medium for use in classical electrodynamics calculations.
subroutine linear_medium_initialize(this)
logical function linear_medium_is_tolerance_reached(this, tol)
integer, parameter medium_parallelepiped
subroutine linear_medium_copy_quantities_to_interaction(partner, interaction)
subroutine linear_medium_restart_write_data(this)
class(linear_medium_t) function, pointer linear_medium_constructor(namespace)
The factory routine (or constructor) allocates a pointer of the corresponding type and then calls the...
logical function linear_medium_restart_read_data(this)
subroutine, public get_medium_box_points_map(medium_box, gr)
subroutine linear_medium_init_interaction(this, interaction)
subroutine, public medium_box_init(medium_box, namespace)
Parse and store geometry of medium box.
subroutine linear_medium_init_interaction_as_partner(partner, interaction)
subroutine linear_medium_update_kinetic_energy(this)
integer, parameter medium_box_file
logical function linear_medium_do_algorithmic_operation(this, operation, updated_quantities)
subroutine linear_medium_finalize(this)
subroutine get_points_map_from_file(filename, mesh, n_points, global_points_number, tmp_map, scale_factor)
Populate list of point indices for points inside the polyhedron.
subroutine, public linear_medium_init(this, namespace)
The init routine is a module level procedure This has the advantage that different classes can have d...
subroutine get_linear_medium_em_properties(this, medium_box, gr)
Evaluate electromagnetic properties of linear medium.
subroutine, public single_medium_box_allocate(medium_box, n_points)
Allocation of medium_box components.
subroutine, public single_medium_box_end(medium_box)
Deallocation of medium_box components.
This module is intended to contain "only mathematical" functions and procedures.
This module defines the meshes, which are used in Octopus.
subroutine, public messages_print_with_emphasis(msg, iunit, namespace)
subroutine, public messages_not_implemented(feature, namespace)
subroutine, public messages_warning(no_lines, all_nodes, namespace)
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
subroutine, public messages_input_error(namespace, var, details, row, column)
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
type(mpi_grp_t), public mpi_world
This module handles the communicators for the various parallelization strategies.
subroutine, public multicomm_end(mc)
subroutine, public multicomm_init(mc, namespace, base_grp, mode_para, n_node, index_range, min_range)
create index and domain communicators
subroutine, public output_linear_medium(outp, namespace, space, mesh, dir, points_map, ep, mu, cc)
subroutine, public output_linear_medium_init(outp, namespace, space)
this module contains the low-level part of the output system
integer function, public parse_block(namespace, name, blk, check_varinfo_)
subroutine, public profiling_out(label)
Increment out counter and sum up difference between entry and exit time.
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
This module implements the basic propagator framework.
This module defines the quantity_t class and the IDs for quantities, which can be exposed by a system...
Implementation details for regridding.
This module implements the abstract system type.
subroutine, public system_end(this)
Descriptor of one algorithmic operation.
Description of the grid, containing information on derivatives, stencil, and symmetries.
abstract interaction class
surrogate interaction class to avoid circular dependencies between modules.
linear medium for classical electrodynamics
Systems (system_t) can expose quantities that can be used to calculate interactions with other system...
contains the information of the meshes and provides the transfer functions
Abstract class for systems.