Octopus
linear_medium.F90
Go to the documentation of this file.
1!! Copyright (C) 2020 F. Bonafé, H. Appel, R. Jestädt
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
27 use comm_oct_m
28 use debug_oct_m
30 use global_oct_m
31 use grid_oct_m
35 use iso_c_binding
36 use io_oct_m
37 use math_oct_m
40 use mesh_oct_m
41 use mpi_oct_m
46 use parser_oct_m
52 use space_oct_m
53 use system_oct_m
55
56 implicit none
57
58 private
59 public :: &
64
65 integer, parameter :: &
66 MEDIUM_PARALLELEPIPED = 1, &
68
71 type, extends(system_t) :: linear_medium_t
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
79 type(grid_t) :: gr
80 type(multicomm_t) :: mc
81 type(single_medium_box_t) :: medium_box
82
83 contains
84 procedure :: init_interaction => linear_medium_init_interaction
85 procedure :: init_interaction_as_partner => linear_medium_init_interaction_as_partner
86 procedure :: initialize => linear_medium_initialize
87 procedure :: do_algorithmic_operation => linear_medium_do_algorithmic_operation
88 procedure :: is_tolerance_reached => linear_medium_is_tolerance_reached
89 procedure :: copy_quantities_to_interaction => linear_medium_copy_quantities_to_interaction
90 procedure :: restart_write_data => linear_medium_restart_write_data
91 procedure :: restart_read_data => linear_medium_restart_read_data
92 procedure :: update_kinetic_energy => linear_medium_update_kinetic_energy
94 end type linear_medium_t
95
96 interface linear_medium_t
97 procedure linear_medium_constructor
98 end interface linear_medium_t
99
100contains
101
102 ! ---------------------------------------------------------
107 function linear_medium_constructor(namespace) result(sys)
108 class(linear_medium_t), pointer :: sys
109 type(namespace_t), intent(in) :: namespace
110
112
113 allocate(sys)
114
115 call linear_medium_init(sys, namespace)
118 end function linear_medium_constructor
119
120 ! ---------------------------------------------------------
125 ! ---------------------------------------------------------
126 subroutine linear_medium_init(this, namespace)
127 class(linear_medium_t), target, intent(inout) :: this
128 type(namespace_t), intent(in) :: namespace
129
130 integer :: nlines, ncols,n_points, n_global_points
131 integer(int64) :: index_range(4)
132 integer, allocatable :: tmp(:)
133 type(block_t) :: blk
134
135 push_sub(linear_medium_init)
136
137 call profiling_in('MEDIUM_BOX_INIT')
138
139 this%namespace = namespace
140 this%space = space_t(this%namespace)
141 if (this%space%is_periodic()) then
142 call messages_not_implemented('Linear medium for periodic systems', namespace=namespace)
143 end if
144 call grid_init_stage_1(this%gr, this%namespace, this%space)
145 ! store the ranges for these two indices (serves as initial guess
146 ! for parallelization strategy)
147 index_range(1) = this%gr%np_global ! Number of points in mesh
148 index_range(2) = 1 ! Number of states
149 index_range(3) = 1 ! Number of k-points
150 index_range(4) = 100000 ! Some large number
151
152 ! create index and domain communicators
153 call multicomm_init(this%mc, this%namespace, mpi_world, calc_mode_par, mpi_world%size, &
154 index_range, (/ 5000, 1, 1, 1 /))
155 call grid_init_stage_2(this%gr, this%namespace, this%space, this%mc)
156
157 call medium_box_init(this%medium_box, this%namespace)
159 !%Variable LinearMediumProperties
160 !%Type block
161 !%Section Maxwell::Medium
162 !%Description
163 !% Defines electromagnetic parameters for a linear medium box.
164 !%
165 !% Example:
166 !%
167 !% <tt>%LinearMediumProperties
168 !% <br>&nbsp;&nbsp; epsilon_factor | mu_factor | sigma_e | sigma_m
169 !% <br>%</tt>
170 !%
171 !% Permittivity factor, permeability factor, electric conductivity and magnetic conductivity of the medium box.
172 !%End
174 if (parse_block(namespace, 'LinearMediumProperties', blk) == 0) then
175
176 nlines = parse_block_n(blk)
177 if (nlines /= 1) then
178 call messages_input_error(namespace, 'LinearMediumProperties', 'should consist of one line', nlines)
179 end if
180 ncols = parse_block_cols(blk, 0)
181 if (ncols /= 4) then
182 call messages_input_error(namespace, 'LinearMediumProperties', 'should consist of four columns')
183 end if
184 call parse_block_float(blk, 0, 0, this%ep_factor)
185 call parse_block_float(blk, 0, 1, this%mu_factor)
186 call parse_block_float(blk, 0, 2, this%sigma_e_factor)
187 call parse_block_float(blk, 0, 3, this%sigma_m_factor)
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
192 write(message(5),'(a)') ""
193 call messages_info(5, namespace=namespace)
194 call parse_block_end(blk)
195
196 call messages_print_with_emphasis(namespace=namespace)
197 else
198 message(1) = 'You must specify the properties of your linear medium through the LinearMediumProperties block.'
199 call messages_fatal(1, namespace=namespace)
200 end if
201
202 !%Variable LinearMediumEdgeProfile
203 !%Type integer
204 !%Section Maxwell::Medium
205 !%Description
206 !% Defines the type of numerical approximation used for the derivatives at the edges of the medium box.
207 !% Default is edged. When the box shape is read from file, only the edged profile is supported.
208 !%
209 !%Option edged 1
210 !% Medium box edges are considered steep for derivatives.
211 !%Option smooth 2
212 !% Medium box edged and softened for derivatives.
213 !%End
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'
219 end if
220 call messages_info(1, namespace=namespace)
221
222 allocate(this%supported_interactions(0))
223 this%supported_interactions_as_partner = [linear_medium_to_em_field]
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.))
228
229 call output_linear_medium_init(this%outp, this%namespace, this%space)
230
231 call get_medium_box_points_map(this%medium_box, this%gr)
232 call get_linear_medium_em_properties(this, this%medium_box, this%gr)
233 call output_linear_medium(this%outp, this%namespace, this%space, this%gr, &
234 this%outp%iter_dir, this%medium_box%points_map, this%medium_box%ep, &
235 this%medium_box%mu, this%medium_box%c)
236
237 if (this%medium_box%check_medium_points) then
238 safe_allocate(tmp(1:this%gr%np))
239 n_global_points = 0
240 write(message(1),'(a, a, a)') 'Check of points inside surface of medium ', trim(this%medium_box%filename), ":"
241 call messages_info(1, namespace=this%namespace)
242 call get_points_map_from_file(this%medium_box%filename, this%gr, n_points, n_global_points, tmp, 0.99_real64)
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
246 write(message(3), '(a)') ""
247 call messages_info(3, namespace=this%namespace)
248 safe_deallocate_a(tmp)
249 end if
250
251 call profiling_out('MEDIUM_BOX_INIT')
252
253 pop_sub(linear_medium_init)
254 end subroutine linear_medium_init
255
256 ! ---------------------------------------------------------
257 subroutine linear_medium_init_interaction(this, interaction)
258 class(linear_medium_t), target, intent(inout) :: this
259 class(interaction_t), intent(inout) :: interaction
260
262
263 select type (interaction)
264 class default
265 message(1) = "Trying to initialize an unsupported interaction by a linear medium."
266 call messages_fatal(1, namespace=this%namespace)
267 end select
268
270 end subroutine linear_medium_init_interaction
271
272 ! ---------------------------------------------------------
273 subroutine linear_medium_init_interaction_as_partner(partner, interaction)
274 class(linear_medium_t), intent(in) :: partner
275 class(interaction_surrogate_t), intent(inout) :: interaction
276
277 type(regridding_t), pointer :: regridding
278 type(single_medium_box_t), pointer :: medium_box_grid
279
281
282 select type (interaction)
284 regridding => regridding_t(interaction%system_gr, &
285 partner%gr, partner%space, partner%namespace)
286
287 ! create a medium box with the size of the partner grid and map all quantities to it
288 medium_box_grid => partner%medium_box%to_grid(partner%gr)
289
290 ! now do the grid transfer to the medium box stored in the interaction which has the size of the system grid
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)
298
299 safe_deallocate_p(regridding)
300 safe_deallocate_p(medium_box_grid)
301 class default
302 message(1) = "Trying to initialize an unsupported interaction by a linear medium."
303 call messages_fatal(1, namespace=partner%namespace)
304 end select
305
308
309 ! ---------------------------------------------------------
310 subroutine linear_medium_initialize(this)
311 class(linear_medium_t), intent(inout) :: this
312
314
316 end subroutine linear_medium_initialize
317
318 ! ---------------------------------------------------------
319 logical function linear_medium_do_algorithmic_operation(this, operation, updated_quantities) result(done)
320 class(linear_medium_t), intent(inout) :: this
321 class(algorithmic_operation_t), intent(in) :: operation
322 character(len=:), allocatable, intent(out) :: updated_quantities(:)
323
325
326 ! Currently there are no linear medium specific algorithmic operations
327 done = .false.
328
331
332 ! ---------------------------------------------------------
333 logical function linear_medium_is_tolerance_reached(this, tol) result(converged)
334 class(linear_medium_t), intent(in) :: this
335 real(real64), intent(in) :: tol
336
338
339 ! this routine is never called at present, no reason to be here
340 assert(.false.)
341 converged = .false.
342
345
346 ! ---------------------------------------------------------
347 subroutine linear_medium_copy_quantities_to_interaction(partner, interaction)
348 class(linear_medium_t), intent(inout) :: partner
349 class(interaction_surrogate_t), intent(inout) :: interaction
352
353 select type (interaction)
355 ! No need to copy anything
356 class default
357 message(1) = "Unsupported interaction."
358 call messages_fatal(1, namespace=partner%namespace)
359 end select
360
363
364 ! ---------------------------------------------------------
365 subroutine linear_medium_restart_write_data(this)
366 class(linear_medium_t), intent(inout) :: this
367
368 integer :: restart_file_unit
369
371
372 message(1) = "Linear medium system "//trim(this%namespace%get())//" will only write a dummy restart file."
373 call messages_warning(1, namespace=this%namespace)
374
375 call io_mkdir('restart/'//td_dir, this%namespace, parents=.true.)
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'
378 call io_close(restart_file_unit)
379
380 message(1) = "Successfully wrote restart data for system "//trim(this%namespace%get())
381 call messages_info(1, namespace=this%namespace)
382
385
386 ! ---------------------------------------------------------
387 ! this function returns true if restart data could be read
388 logical function linear_medium_restart_read_data(this)
389 class(linear_medium_t), intent(inout) :: this
390
391 integer :: restart_file_unit
392
394
395 call io_mkdir('restart/'//td_dir, this%namespace, parents=.true.)
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
398 call io_close(restart_file_unit)
400 else
401 ! could not open file
403 end if
404
406 message(1) = "Successfully read restart data for system "//trim(this%namespace%get())
407 call messages_info(1, namespace=this%namespace)
408 end if
409
414 class(linear_medium_t), intent(inout) :: this
415
417
418 this%kinetic_energy = m_zero
419
421
423
424 ! ---------------------------------------------------------
425 subroutine linear_medium_finalize(this)
426 type(linear_medium_t), intent(inout) :: this
427
428 push_sub(linear_medium_finalize)
429
430 call single_medium_box_end(this%medium_box)
431 call system_end(this)
432 call multicomm_end(this%mc)
433 call grid_end(this%gr)
434
436 end subroutine linear_medium_finalize
437
438
439 ! Specific routines for this system:
441 ! ---------------------------------------------------------
442 subroutine get_medium_box_points_map(medium_box, gr)
443 type(single_medium_box_t), intent(inout) :: medium_box
444 type(grid_t), intent(in) :: gr
445
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)
449 logical :: inside
450
452
453 call profiling_in('GET_MEDIUM_BOX_POINTS_MAP')
454
455 safe_allocate(tmp_points_map(1:gr%np))
456 safe_allocate(tmp_bdry_map(1:gr%np))
457 tmp_points_map = 0
458 tmp_bdry_map = 0
459
460 if (medium_box%box_shape == medium_box_file) then
461
462 call get_points_map_from_file(medium_box%filename, gr, medium_box%points_number,&
463 medium_box%global_points_number, tmp_points_map)
464
465 else
466
467 do idim = 1, 3
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
470 end do
471 ip_in = 0
472 ip_bd = 0
473 do ip = 1, gr%np
474 xx(1:3) = gr%x(ip, 1:3)
475 inside = check_point_in_bounds(xx, bounds(:,:))
476 if (check_point_in_bounds(xx, bounds(:,:))) then
477 ip_in = ip_in + 1
478 tmp_points_map(ip_in) = ip
479 end if
480 if (check_point_on_bounds(xx, bounds(:,:))) then
481 ip_bd = ip_bd + 1
482 tmp_bdry_map(ip_bd) = ip
483 end if
484 end do
485
486 medium_box%points_number = ip_in
487 medium_box%bdry_number = ip_bd
488
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)
492
493 end if
494 call single_medium_box_allocate(medium_box, medium_box%points_number)
495 medium_box%points_map(:) = tmp_points_map(1:medium_box%points_number)
496
497 ! TODO: add some check that medium boxes do not overlap (now they are different systems)
498 safe_deallocate_a(tmp_points_map)
499 safe_deallocate_a(tmp_bdry_map)
500
501 call profiling_out('GET_MEDIUM_BOX_POINTS_MAP')
502
504 contains
505
506 logical pure function check_point_in_bounds(xx, bounds) result (check)
507 real(real64), intent(in) :: xx(:)
508 real(real64), intent(in) :: bounds(:,:)
509
510 check = .false.
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
514 check = .true.
515 end if
516
517 end function check_point_in_bounds
519 logical pure function check_point_on_bounds(xx, bounds) result (check)
520 real(real64), intent(in) :: xx(:)
521 real(real64), intent(in) :: bounds(:,:)
522
523 check = .false.
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
536 check = .true.
537 end if
538
539 end function check_point_on_bounds
540
541 end subroutine get_medium_box_points_map
542
543 ! ----------------------------------------------------------
545 subroutine get_linear_medium_em_properties(this, medium_box, gr)
546 type(linear_medium_t), intent(in) :: this
547 type(single_medium_box_t), intent(inout) :: medium_box
548 type(grid_t), intent(in) :: gr
549
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(:,:)
553
555
556 call profiling_in('GET_LIN_MEDIUM_EM_PROP')
557
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))
561
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))
566
567 xx(1:3) = gr%x(ip,1:3)
568 dd_min = m_huge
569
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
575 end do
576
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)))
586
587 else if (this%edge_profile == option__linearmediumedgeprofile__edged) then
588
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
594
595 end if
596 end do
597
598 tmp(:) = p_ep
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)
602 end do
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))
607 end do
608
609 tmp(:) = p_mu
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)
613 end do
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))
618 end do
619
620 !TODO: add print information about the medium box
621
622 safe_deallocate_a(tmp)
623 safe_deallocate_a(tmp_grad)
624
625 call profiling_out('GET_LIN_MEDIUM_EM_PROP')
626
629
630 ! ----------------------------------------------------------
632 subroutine get_points_map_from_file(filename, mesh, n_points, global_points_number, tmp_map, scale_factor)
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
639
640 integer :: ip_in, ip
641 real(real64) :: xx(3)
642 type(cgal_polyhedra_t) :: cgal_poly
643
645
646 call profiling_in('GET_POINTS_MAP_FROM_FILE')
647
648 call cgal_polyhedron_init(cgal_poly, trim(filename), verbose = .false.)
649
650 ip_in = 0
651 do ip = 1, mesh%np
652 if (present(scale_factor)) then
653 xx(1:3) = scale_factor * mesh%x(ip, 1:3)
654 else
655 xx(1:3) = mesh%x(ip, 1:3)
656 end if
657 if (cgal_polyhedron_point_inside(cgal_poly, xx(1), xx(2), xx(3))) then
658 ip_in = ip_in + 1
659 tmp_map(ip_in) = ip
660 end if
661 end do
662 n_points = ip_in
663 call cgal_polyhedron_end(cgal_poly)
664
665 call mpi_world%allreduce(ip_in, global_points_number, 1, mpi_integer, mpi_sum)
666
667 call profiling_out('GET_POINTS_MAP_FROM_FILE')
668
670
671 end subroutine get_points_map_from_file
672
673 ! ----------------------------------------------------------
675 subroutine medium_box_init(medium_box, namespace)
676 type(single_medium_box_t), intent(inout) :: medium_box
677 type(namespace_t), intent(in) :: namespace
678
679 integer :: nlines, ncols, idim
680 type(block_t) :: blk
681
682 push_sub(medium_box_init)
683
684 !%Variable LinearMediumBoxShape
685 !%Type integer
686 !%Section Maxwell::Medium
687 !%Description
688 !% This variable defines the shape of the linear medium box.
689 !% The default is <tt>medium_parallelepiped</tt>.
690 !%Option medium_parallelepiped 1
691 !% The medium box will be a parallelepiped whose center and dimensions are taken from
692 !% the variable <tt>LinearMediumBoxSize</tt>.
693 !%Option medium_box_file 2
694 !% The simulation box will be read from an external file in OFF format, defined by the variable <tt>LinearMediumBoxFile</tt>.
695 !%End
696 call parse_variable(namespace, 'LinearMediumBoxShape', medium_parallelepiped, medium_box%box_shape)
697
698 if (.not. varinfo_valid_option('LinearMediumBoxShape', medium_box%box_shape)) then
699 call messages_input_error(namespace, 'LinearMediumBoxShape')
700 end if
701
702 select case (medium_box%box_shape)
704
705 !%Variable LinearMediumBoxSize
706 !%Type block
707 !%Section Maxwell::Medium
708 !%Description
709 !% Defines center and size of a parallelepiped linear medium box.
710 !%
711 !% Example:
712 !%
713 !% <tt>%LinearMediumBoxSize
714 !% <br>&nbsp;&nbsp; center_x | center_y | center_z | x_length | y_length | z_length
715 !% <br>%</tt>
716 !%End
717
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)
720
721 nlines = parse_block_n(blk)
722 if (nlines /= 1) then
723 call messages_input_error(namespace, 'LinearMediumBoxSize', 'should consist of one line')
724 end if
725 ncols = parse_block_cols(blk, 0)
726 if (ncols /= 6) then
727 call messages_input_error(namespace, 'LinearMediumBoxSize', 'should consist of six columns')
728 end if
729 do idim = 1, 3
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))
732 end do
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)
740
741 call messages_print_with_emphasis(namespace=namespace)
742 else
743 message(1) = "For parallelepiped box shapes, you must provide a LinearMediumBoxSize block."
744 call messages_fatal(1, namespace=namespace)
745 end if
746
747 case (medium_box_file)
748 !%Variable LinearMediumBoxFile
749 !%Type string
750 !%Section Maxwell::Medium
751 !%Description
752 !% File in OFF format with the shape of the linear medium.
753 !%End
754 if (parse_is_defined(namespace, 'LinearMediumBoxFile')) then
755 call parse_variable(namespace, 'LinearMediumBoxFile', 'mediumboxfile', medium_box%filename)
756 else
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)
759 end if
760
761 !%Variable CheckPointsMediumFromFile
762 !%Type logical
763 !%Default no
764 !%Section Maxwell::Medium
765 !%Description
766 !% Whether to re-calculate the points map by artificially shrinking the coordinate system by a factor of
767 !% 0.99 to check if the points inside the medium surface are properly detected. This works for only one
768 !% medium surface which is centered in the origin of the coordinate system.
769 !%End
770 call parse_variable(namespace, 'CheckPointsMediumFromFile', .false., medium_box%check_medium_points)
771
772 end select
773
774 pop_sub(medium_box_init)
775
776 end subroutine medium_box_init
777
778end module linear_medium_oct_m
779
780!! Local Variables:
781!! mode: f90
782!! coding: utf-8
783!! End:
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.
Definition: algorithm.F90:141
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
Definition: global.F90:190
real(real64), parameter, public m_zero
Definition: global.F90:188
character(len= *), parameter, public td_dir
Definition: global.F90:251
This module implements the underlying real-space grid.
Definition: grid.F90:117
subroutine, public grid_init_stage_1(gr, namespace, space, symm, latt, n_sites, site_position)
First stage of the grid initialization.
Definition: grid.F90:194
subroutine, public grid_init_stage_2(gr, namespace, space, mc, qvector)
Second stage of the grid initialization.
Definition: grid.F90:465
subroutine, public grid_end(gr)
finalize a grid object
Definition: grid.F90:493
integer, parameter, public linear_medium_to_em_field
This module defines the abstract interaction_t class, and some auxiliary classes for interactions.
Definition: io.F90:114
subroutine, public io_close(iunit, grp)
Definition: io.F90:418
subroutine, public io_mkdir(fname, namespace, parents)
Definition: io.F90:311
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
Definition: io.F90:352
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.
Definition: math.F90:115
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
subroutine, public messages_print_with_emphasis(msg, iunit, namespace)
Definition: messages.F90:921
subroutine, public messages_not_implemented(feature, namespace)
Definition: messages.F90:1114
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:538
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:161
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:415
subroutine, public messages_input_error(namespace, var, details, row, column)
Definition: messages.F90:714
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:617
type(mpi_grp_t), public mpi_world
Definition: mpi.F90:266
This module handles the communicators for the various parallelization strategies.
Definition: multicomm.F90:145
subroutine, public multicomm_end(mc)
Definition: multicomm.F90:897
subroutine, public multicomm_init(mc, namespace, base_grp, mode_para, n_node, index_range, min_range)
create index and domain communicators
Definition: multicomm.F90:268
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
Definition: output_low.F90:115
integer function, public parse_block(namespace, name, blk, check_varinfo_)
Definition: parser.F90:618
subroutine, public profiling_out(label)
Increment out counter and sum up difference between entry and exit time.
Definition: profiling.F90:623
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
Definition: profiling.F90:552
This module implements the basic propagator framework.
Definition: propagator.F90:117
This module defines the quantity_t class and the IDs for quantities, which can be exposed by a system...
Definition: quantity.F90:138
Implementation details for regridding.
Definition: regridding.F90:170
This module implements the abstract system type.
Definition: system.F90:118
subroutine, public system_end(this)
Definition: system.F90:1145
Descriptor of one algorithmic operation.
Definition: algorithm.F90:163
Description of the grid, containing information on derivatives, stencil, and symmetries.
Definition: grid.F90:168
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...
Definition: quantity.F90:171
contains the information of the meshes and provides the transfer functions
Definition: regridding.F90:196
Abstract class for systems.
Definition: system.F90:172
int true(void)