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
39 use math_oct_m
41 use mesh_oct_m
42 use mpi_oct_m
47 use parser_oct_m
53 use space_oct_m
54 use system_oct_m
56
57 implicit none
58
59 private
60 public :: &
65
66 integer, parameter :: &
67 MEDIUM_PARALLELEPIPED = 1, &
69
72 type, extends(system_t) :: linear_medium_t
73 type(space_t) :: space
74 real(real64) :: ep_factor
75 real(real64) :: mu_factor
76 real(real64) :: sigma_e_factor
77 real(real64) :: sigma_m_factor
78 integer :: edge_profile
79 type(output_t) :: outp
80 type(grid_t) :: gr
81 type(multicomm_t) :: mc
82 type(single_medium_box_t) :: medium_box
83
84 contains
85 procedure :: init_interaction => linear_medium_init_interaction
86 procedure :: init_interaction_as_partner => linear_medium_init_interaction_as_partner
87 procedure :: initial_conditions => linear_medium_initial_conditions
88 procedure :: do_algorithmic_operation => linear_medium_do_algorithmic_operation
89 procedure :: is_tolerance_reached => linear_medium_is_tolerance_reached
90 procedure :: copy_quantities_to_interaction => linear_medium_copy_quantities_to_interaction
91 procedure :: restart_write_data => linear_medium_restart_write_data
92 procedure :: restart_read_data => linear_medium_restart_read_data
93 procedure :: update_kinetic_energy => linear_medium_update_kinetic_energy
95 end type linear_medium_t
96
97 interface linear_medium_t
98 procedure linear_medium_constructor
99 end interface linear_medium_t
100
101contains
102
103 ! ---------------------------------------------------------
108 function linear_medium_constructor(namespace) result(sys)
109 class(linear_medium_t), pointer :: sys
110 type(namespace_t), intent(in) :: namespace
111
113
114 allocate(sys)
115
116 call linear_medium_init(sys, namespace)
117
119 end function linear_medium_constructor
120
121 ! ---------------------------------------------------------
126 ! ---------------------------------------------------------
127 subroutine linear_medium_init(this, namespace)
128 class(linear_medium_t), target, intent(inout) :: this
129 type(namespace_t), intent(in) :: namespace
130
131 integer :: nlines, ncols,n_points, n_global_points
132 integer(int64) :: index_range(4)
133 integer, allocatable :: tmp(:)
134 type(block_t) :: blk
135
136 push_sub(linear_medium_init)
137
138 call profiling_in('MEDIUM_BOX_INIT')
139
140 this%namespace = namespace
141 this%space = space_t(this%namespace)
142 if (this%space%is_periodic()) then
143 call messages_not_implemented('Linear medium for periodic systems', namespace=namespace)
144 end if
145 call grid_init_stage_1(this%gr, this%namespace, this%space)
146 ! store the ranges for these two indices (serves as initial guess
147 ! for parallelization strategy)
148 index_range(1) = this%gr%np_global ! Number of points in mesh
149 index_range(2) = 1 ! Number of states
150 index_range(3) = 1 ! Number of k-points
151 index_range(4) = 100000 ! Some large number
152
153 ! create index and domain communicators
154 call multicomm_init(this%mc, this%namespace, mpi_world, calc_mode_par, mpi_world%size, &
155 index_range, (/ 5000, 1, 1, 1 /))
156 call grid_init_stage_2(this%gr, this%namespace, this%space, this%mc)
157
158 call medium_box_init(this%medium_box, this%namespace)
160 !%Variable LinearMediumProperties
161 !%Type block
162 !%Section Maxwell::Medium
163 !%Description
164 !% Defines electromagnetic parameters for a linear medium box.
165 !%
166 !% Example:
167 !%
168 !% <tt>%LinearMediumProperties
169 !% <br>&nbsp;&nbsp; epsilon_factor | mu_factor | sigma_e | sigma_m
170 !% <br>%</tt>
171 !%
172 !% Permittivity factor, permeability factor, electric conductivity and magnetic conductivity of the medium box.
173 !%End
175 if (parse_block(namespace, 'LinearMediumProperties', blk) == 0) then
176
177 nlines = parse_block_n(blk)
178 if (nlines /= 1) then
179 call messages_input_error(namespace, 'LinearMediumProperties', 'should consist of one line', nlines)
180 end if
181 ncols = parse_block_cols(blk, 0)
182 if (ncols /= 4) then
183 call messages_input_error(namespace, 'LinearMediumProperties', 'should consist of four columns')
184 end if
185 call parse_block_float(blk, 0, 0, this%ep_factor)
186 call parse_block_float(blk, 0, 1, this%mu_factor)
187 call parse_block_float(blk, 0, 2, this%sigma_e_factor)
188 call parse_block_float(blk, 0, 3, this%sigma_m_factor)
189 write(message(1),'(a,es9.2)') 'Box epsilon factor: ', this%ep_factor
190 write(message(2),'(a,es9.2)') 'Box mu factor: ', this%mu_factor
191 write(message(3),'(a,es9.2)') 'Box electric sigma: ', this%sigma_e_factor
192 write(message(4),'(a,es9.2)') 'Box magnetic sigma: ', this%sigma_m_factor
193 write(message(5),'(a)') ""
194 call messages_info(5, namespace=namespace)
195 call parse_block_end(blk)
196
197 call messages_print_with_emphasis(namespace=namespace)
198 else
199 message(1) = 'You must specify the properties of your linear medium through the LinearMediumProperties block.'
200 call messages_fatal(1, namespace=namespace)
201 end if
202
203 !%Variable LinearMediumEdgeProfile
204 !%Type integer
205 !%Section Maxwell::Medium
206 !%Description
207 !% Defines the type of numerical approximation used for the derivatives at the edges of the medium box.
208 !% Default is edged. When the box shape is read from file, only the edged profile is supported.
209 !%
210 !%Option edged 1
211 !% Medium box edges are considered steep for derivatives.
212 !%Option smooth 2
213 !% Medium box edged and softened for derivatives.
214 !%End
215 call parse_variable(namespace, 'LinearMediumEdgeProfile', option__linearmediumedgeprofile__edged, this%edge_profile)
216 if (this%edge_profile == option__linearmediumedgeprofile__edged) then
217 write(message(1),'(a,a)') 'Box shape: ', 'edged'
218 else if (this%edge_profile == option__linearmediumedgeprofile__smooth) then
219 write(message(1),'(a,a)') 'Box shape: ', 'smooth'
220 end if
221 call messages_info(1, namespace=namespace)
222
223 allocate(this%supported_interactions(0))
224 this%supported_interactions_as_partner = [linear_medium_to_em_field]
225 this%quantities(permittivity)%required = .true.
226 this%quantities(permittivity)%updated_on_demand = .false.
227 this%quantities(permittivity)%always_available = .true.
228
229 this%quantities(permeability)%required = .true.
230 this%quantities(permeability)%updated_on_demand = .false.
231 this%quantities(permeability)%always_available = .true.
232
233 this%quantities(e_conductivity)%required = .true.
234 this%quantities(e_conductivity)%updated_on_demand = .false.
235 this%quantities(e_conductivity)%always_available = .true.
236
237 this%quantities(m_conductivity)%required = .true.
238 this%quantities(m_conductivity)%updated_on_demand = .false.
239 this%quantities(m_conductivity)%always_available = .true.
240
241 call output_linear_medium_init(this%outp, this%namespace, this%space)
242
243 call get_medium_box_points_map(this%medium_box, this%gr)
244 call get_linear_medium_em_properties(this, this%medium_box, this%gr)
245 call output_linear_medium(this%outp, this%namespace, this%space, this%gr, &
246 this%outp%iter_dir, this%medium_box%points_map, this%medium_box%ep, &
247 this%medium_box%mu, this%medium_box%c)
248
249 if (this%medium_box%check_medium_points) then
250 safe_allocate(tmp(1:this%gr%np))
251 n_global_points = 0
252 write(message(1),'(a, a, a)') 'Check of points inside surface of medium ', trim(this%medium_box%filename), ":"
253 call messages_info(1, namespace=this%namespace)
254 call get_points_map_from_file(this%medium_box%filename, this%gr, n_points, n_global_points, tmp, 0.99_real64)
255 write(message(1),'(a, I8)')'Number of points inside medium (normal coordinates):', &
256 this%medium_box%global_points_number
257 write(message(2),'(a, I8)')'Number of points inside medium (rescaled coordinates):', n_global_points
258 write(message(3), '(a)') ""
259 call messages_info(3, namespace=this%namespace)
260 safe_deallocate_a(tmp)
261 end if
262
263 call profiling_out('MEDIUM_BOX_INIT')
264
265 pop_sub(linear_medium_init)
266 end subroutine linear_medium_init
267
268 ! ---------------------------------------------------------
269 subroutine linear_medium_init_interaction(this, interaction)
270 class(linear_medium_t), target, intent(inout) :: this
271 class(interaction_t), intent(inout) :: interaction
272
274
275 select type (interaction)
276 class default
277 message(1) = "Trying to initialize an unsupported interaction by a linear medium."
278 call messages_fatal(1, namespace=this%namespace)
279 end select
280
282 end subroutine linear_medium_init_interaction
283
284 ! ---------------------------------------------------------
285 subroutine linear_medium_init_interaction_as_partner(partner, interaction)
286 class(linear_medium_t), intent(in) :: partner
287 class(interaction_surrogate_t), intent(inout) :: interaction
288
289 type(regridding_t), pointer :: regridding
290 type(single_medium_box_t), pointer :: medium_box_grid
291
293
294 select type (interaction)
296 regridding => regridding_t(interaction%system_gr, &
297 partner%gr, partner%space, partner%namespace)
298
299 ! create a medium box with the size of the partner grid and map all quantities to it
300 medium_box_grid => partner%medium_box%to_grid(partner%gr)
301
302 ! now do the grid transfer to the medium box stored in the interaction which has the size of the system grid
303 call regridding%do_transfer(interaction%medium_box%aux_ep, medium_box_grid%aux_ep)
304 call regridding%do_transfer(interaction%medium_box%aux_mu, medium_box_grid%aux_mu)
305 call regridding%do_transfer(interaction%medium_box%ep, medium_box_grid%ep)
306 call regridding%do_transfer(interaction%medium_box%mu, medium_box_grid%mu)
307 call regridding%do_transfer(interaction%medium_box%c, medium_box_grid%c)
308 call regridding%do_transfer(interaction%medium_box%sigma_e, medium_box_grid%sigma_e)
309 call regridding%do_transfer(interaction%medium_box%sigma_m, medium_box_grid%sigma_m)
310
311 safe_deallocate_p(regridding)
312 safe_deallocate_p(medium_box_grid)
313 class default
314 message(1) = "Trying to initialize an unsupported interaction by a linear medium."
315 call messages_fatal(1, namespace=partner%namespace)
316 end select
317
320
321 ! ---------------------------------------------------------
322 subroutine linear_medium_initial_conditions(this)
323 class(linear_medium_t), intent(inout) :: this
324
326
329
330 ! ---------------------------------------------------------
331 logical function linear_medium_do_algorithmic_operation(this, operation, updated_quantities) result(done)
332 class(linear_medium_t), intent(inout) :: this
333 class(algorithmic_operation_t), intent(in) :: operation
334 integer, allocatable, intent(out) :: updated_quantities(:)
335
337
338 ! Currently there are no linear medium specific algorithmic operations
339 done = .false.
340
343
344 ! ---------------------------------------------------------
345 logical function linear_medium_is_tolerance_reached(this, tol) result(converged)
346 class(linear_medium_t), intent(in) :: this
347 real(real64), intent(in) :: tol
348
350
351 ! this routine is never called at present, no reason to be here
352 assert(.false.)
353 converged = .false.
354
357
358 ! ---------------------------------------------------------
359 subroutine linear_medium_copy_quantities_to_interaction(partner, interaction)
360 class(linear_medium_t), intent(inout) :: partner
361 class(interaction_surrogate_t), intent(inout) :: interaction
364
365 select type (interaction)
367 ! No need to copy anything
368 class default
369 message(1) = "Unsupported interaction."
370 call messages_fatal(1, namespace=partner%namespace)
371 end select
372
375
376 ! ---------------------------------------------------------
377 subroutine linear_medium_restart_write_data(this)
378 class(linear_medium_t), intent(inout) :: this
379
380 integer :: restart_file_unit
381
383
384 message(1) = "Linear medium system "//trim(this%namespace%get())//" will only write a dummy restart file."
385 call messages_warning(1, namespace=this%namespace)
386
387 call io_mkdir('restart/'//td_dir, this%namespace, parents=.true.)
388 restart_file_unit = io_open('restart/'//td_dir// 'restart_linear_medium', this%namespace, action='write')
389 write(restart_file_unit, *) 'Linear medium restart file'
390 call io_close(restart_file_unit)
391
392 message(1) = "Successfully wrote restart data for system "//trim(this%namespace%get())
393 call messages_info(1, namespace=this%namespace)
394
397
398 ! ---------------------------------------------------------
399 ! this function returns true if restart data could be read
400 logical function linear_medium_restart_read_data(this)
401 class(linear_medium_t), intent(inout) :: this
402
403 integer :: restart_file_unit
404
406
407 call io_mkdir('restart/'//td_dir, this%namespace, parents=.true.)
408 restart_file_unit = io_open('restart/'//td_dir// 'restart_linear_medium', this%namespace, action='read', die=.false.)
409 if (restart_file_unit > 0) then
410 call io_close(restart_file_unit)
412 else
413 ! could not open file
415 end if
416
418 message(1) = "Successfully read restart data for system "//trim(this%namespace%get())
419 call messages_info(1, namespace=this%namespace)
420 end if
421
426 class(linear_medium_t), intent(inout) :: this
427
429
430 this%kinetic_energy = m_zero
431
433
435
436 ! ---------------------------------------------------------
437 subroutine linear_medium_finalize(this)
438 type(linear_medium_t), intent(inout) :: this
439
440 push_sub(linear_medium_finalize)
441
442 call single_medium_box_end(this%medium_box)
443 call system_end(this)
444 call multicomm_end(this%mc)
445 call grid_end(this%gr)
446
448 end subroutine linear_medium_finalize
449
450
451 ! Specific routines for this system:
453 ! ---------------------------------------------------------
454 subroutine get_medium_box_points_map(medium_box, gr)
455 type(single_medium_box_t), intent(inout) :: medium_box
456 type(grid_t), intent(in) :: gr
457
458 integer :: ip, ip_in, ip_bd, idim
459 integer, allocatable :: tmp_points_map(:), tmp_bdry_map(:)
460 real(real64) :: bounds(2,gr%box%dim), xx(gr%box%dim)
461 logical :: inside
462
464
465 call profiling_in('GET_MEDIUM_BOX_POINTS_MAP')
466
467 safe_allocate(tmp_points_map(1:gr%np))
468 safe_allocate(tmp_bdry_map(1:gr%np))
469 tmp_points_map = 0
470 tmp_bdry_map = 0
471
472 if (medium_box%box_shape == medium_box_file) then
473
474 call get_points_map_from_file(medium_box%filename, gr, medium_box%points_number,&
475 medium_box%global_points_number, tmp_points_map)
476
477 else
478
479 do idim = 1, 3
480 bounds(1,idim) = medium_box%center(idim) - medium_box%lsize(idim)/m_two
481 bounds(2,idim) = medium_box%center(idim) + medium_box%lsize(idim)/m_two
482 end do
483 ip_in = 0
484 ip_bd = 0
485 do ip = 1, gr%np
486 xx(1:3) = gr%x(ip, 1:3)
487 inside = check_point_in_bounds(xx, bounds(:,:))
488 if (check_point_in_bounds(xx, bounds(:,:))) then
489 ip_in = ip_in + 1
490 tmp_points_map(ip_in) = ip
491 end if
492 if (check_point_on_bounds(xx, bounds(:,:))) then
493 ip_bd = ip_bd + 1
494 tmp_bdry_map(ip_bd) = ip
495 end if
496 end do
497
498 medium_box%points_number = ip_in
499 medium_box%bdry_number = ip_bd
500
501 safe_allocate(medium_box%bdry_map(1:medium_box%bdry_number))
502 medium_box%bdry_map = 0
503 medium_box%bdry_map = tmp_bdry_map(1:medium_box%bdry_number)
504
505 end if
506 call single_medium_box_allocate(medium_box, medium_box%points_number)
507 medium_box%points_map(:) = tmp_points_map(1:medium_box%points_number)
508
509 ! TODO: add some check that medium boxes do not overlap (now they are different systems)
510 safe_deallocate_a(tmp_points_map)
511 safe_deallocate_a(tmp_bdry_map)
512
513 call profiling_out('GET_MEDIUM_BOX_POINTS_MAP')
514
516 contains
517
518 logical pure function check_point_in_bounds(xx, bounds) result (check)
519 real(real64), intent(in) :: xx(:)
520 real(real64), intent(in) :: bounds(:,:)
521
522 check = .false.
523 if ((xx(1) >= bounds(1,1)) .and. (xx(1) <= bounds(2,1)) .and. &
524 (xx(2) >= bounds(1,2)) .and. (xx(2) <= bounds(2,2)) .and. &
525 (xx(3) >= bounds(1,3)) .and. (xx(3) <= bounds(2,3)) ) then
526 check = .true.
527 end if
528
529 end function check_point_in_bounds
531 logical pure function check_point_on_bounds(xx, bounds) result (check)
532 real(real64), intent(in) :: xx(:)
533 real(real64), intent(in) :: bounds(:,:)
534
535 check = .false.
536 if (is_close(xx(1), bounds(1,1)) .and. (xx(2) >= bounds(1,2) .and. xx(3) >= bounds(1,3)) &
537 .and. (xx(2) <= bounds(2,2) .and. xx(3) <= bounds(2,3)) .or. &
538 is_close(xx(2), bounds(1,2)) .and. (xx(1) >= bounds(1,1) .and. xx(3) >= bounds(1,3)) &
539 .and. (xx(1) <= bounds(2,1) .and. xx(3) <= bounds(2,3)) .or. &
540 is_close(xx(3), bounds(1,3)) .and. (xx(1) >= bounds(1,1) .and. xx(2) >= bounds(1,2)) &
541 .and. (xx(1) <= bounds(2,1) .and. xx(2) <= bounds(2,2)) .or. &
542 is_close(xx(1), bounds(2,1)) .and. (xx(2) >= bounds(1,2) .and. xx(3) >= bounds(1,3)) &
543 .and. (xx(2) <= bounds(2,2) .and. xx(3) <= bounds(2,3)) .or. &
544 is_close(xx(2), bounds(2,2)) .and. (xx(1) >= bounds(1,1) .and. xx(3) >= bounds(1,3)) &
545 .and. (xx(1) <= bounds(2,1) .and. xx(3) <= bounds(2,3)) .or. &
546 is_close(xx(3), bounds(2,3)) .and. (xx(1) >= bounds(1,1) .and. xx(2) >= bounds(1,2)) &
547 .and. (xx(1) <= bounds(2,1) .and. xx(2) <= bounds(2,2)) ) then
548 check = .true.
549 end if
550
551 end function check_point_on_bounds
552
553 end subroutine get_medium_box_points_map
554
555 ! ----------------------------------------------------------
557 subroutine get_linear_medium_em_properties(this, medium_box, gr)
558 type(linear_medium_t), intent(in) :: this
559 type(single_medium_box_t), intent(inout) :: medium_box
560 type(grid_t), intent(in) :: gr
561
562 integer :: ip, ip_in, ip_bd, ipp
563 real(real64) :: xx(gr%box%dim), xxp(gr%box%dim), dd, dd_max, dd_min
564 real(real64), allocatable :: tmp(:), tmp_grad(:,:)
565
567
568 call profiling_in('GET_LIN_MEDIUM_EM_PROP')
569
570 safe_allocate(tmp(1:gr%np_part))
571 safe_allocate(tmp_grad(1:gr%np_part, 1:gr%box%dim))
572 dd_max = max(2*gr%spacing(1), 2*gr%spacing(2), 2*gr%spacing(3))
573
574 do ip_in = 1, medium_box%points_number
575 ip = medium_box%points_map(ip_in)
576 if (this%edge_profile == option__linearmediumedgeprofile__smooth) then
577 assert(allocated(medium_box%bdry_map))
578
579 xx(1:3) = gr%x(ip,1:3)
580 dd_min = m_huge
581
582 do ip_bd = 1, medium_box%bdry_number
583 ipp = medium_box%bdry_map(ip_bd)
584 xxp(1:3) = gr%x(ipp,1:3)
585 dd = norm2(xx(1:3) - xxp(1:3))
586 if (dd < dd_min) dd_min = dd
587 end do
588
589 medium_box%ep(ip_in) = p_ep + ((p_ep * this%ep_factor - p_ep) &
590 * m_one/(m_one + exp(-m_five/dd_max * (dd_min - m_two*dd_max))))
591 medium_box%mu(ip_in) = p_mu + ((p_mu * this%mu_factor - p_mu) &
592 * m_one/(m_one + exp(-m_five/dd_max * (dd_min - m_two*dd_max))))
593 medium_box%c(ip_in) = m_one/sqrt(medium_box%ep(ip_in)*medium_box%mu(ip_in))
594 medium_box%sigma_e(ip_in) = this%sigma_e_factor &
595 * m_one/(m_one + exp(-m_five/dd_max * (dd_min - m_two*dd_max)))
596 medium_box%sigma_m(ip_in) = this%sigma_m_factor &
597 * m_one/(m_one + exp(-m_five/dd_max * (dd_min - m_two*dd_max)))
598
599 else if (this%edge_profile == option__linearmediumedgeprofile__edged) then
600
601 medium_box%ep(ip_in) = p_ep * this%ep_factor
602 medium_box%mu(ip_in) = p_mu * this%mu_factor
603 medium_box%c(ip_in) = m_one/sqrt(medium_box%ep(ip_in)*medium_box%mu(ip_in))
604 medium_box%sigma_e(ip_in) = this%sigma_e_factor
605 medium_box%sigma_m(ip_in) = this%sigma_m_factor
606
607 end if
608 end do
609
610 tmp(:) = p_ep
611 do ip_in = 1, medium_box%points_number
612 ip = medium_box%points_map(ip_in)
613 tmp(ip)= medium_box%ep(ip_in)
614 end do
615 call dderivatives_grad(gr%der, tmp, tmp_grad, set_bc = .false.)
616 do ip_in = 1, medium_box%points_number
617 ip = medium_box%points_map(ip_in)
618 medium_box%aux_ep(ip_in, :) = tmp_grad(ip, :)/(m_four * medium_box%ep(ip_in))
619 end do
620
621 tmp(:) = p_mu
622 do ip_in = 1, medium_box%points_number
623 ip = medium_box%points_map(ip_in)
624 tmp(ip) = medium_box%mu(ip_in)
625 end do
626 call dderivatives_grad(gr%der, tmp, tmp_grad, set_bc = .false.)
627 do ip_in = 1, medium_box%points_number
628 ip = medium_box%points_map(ip_in)
629 medium_box%aux_mu(ip_in, :) = tmp_grad(ip, :)/(m_four * medium_box%mu(ip_in))
630 end do
631
632 !TODO: add print information about the medium box
633
634 safe_deallocate_a(tmp)
635 safe_deallocate_a(tmp_grad)
636
637 call profiling_out('GET_LIN_MEDIUM_EM_PROP')
638
641
642 ! ----------------------------------------------------------
644 subroutine get_points_map_from_file(filename, mesh, n_points, global_points_number, tmp_map, scale_factor)
645 character(len=256), intent(in) :: filename
646 class(mesh_t), intent(in) :: mesh
647 integer, intent(out) :: n_points
648 integer, intent(out) :: global_points_number
649 integer, intent(inout) :: tmp_map(:)
650 real(real64), optional, intent(in) :: scale_factor
651
652 integer :: ip_in, ip
653 real(real64) :: xx(3)
654 type(cgal_polyhedra_t) :: cgal_poly
655
657
658 call profiling_in('GET_POINTS_MAP_FROM_FILE')
659
660 call cgal_polyhedron_init(cgal_poly, trim(filename), verbose = .false.)
661
662 ip_in = 0
663 do ip = 1, mesh%np
664 if (present(scale_factor)) then
665 xx(1:3) = scale_factor * mesh%x(ip, 1:3)
666 else
667 xx(1:3) = mesh%x(ip, 1:3)
668 end if
669 if (cgal_polyhedron_point_inside(cgal_poly, xx(1), xx(2), xx(3))) then
670 ip_in = ip_in + 1
671 tmp_map(ip_in) = ip
672 end if
673 end do
674 n_points = ip_in
675 call cgal_polyhedron_end(cgal_poly)
676
677 call mpi_world%allreduce(ip_in, global_points_number, 1, mpi_integer, mpi_sum)
678
679 call profiling_out('GET_POINTS_MAP_FROM_FILE')
680
682
683 end subroutine get_points_map_from_file
684
685 ! ----------------------------------------------------------
687 subroutine medium_box_init(medium_box, namespace)
688 type(single_medium_box_t), intent(inout) :: medium_box
689 type(namespace_t), intent(in) :: namespace
690
691 integer :: nlines, ncols, idim
692 type(block_t) :: blk
693
694 push_sub(medium_box_init)
695
696 !%Variable LinearMediumBoxShape
697 !%Type integer
698 !%Section Maxwell::Medium
699 !%Description
700 !% This variable defines the shape of the linear medium box.
701 !% The default is <tt>medium_parallelepiped</tt>.
702 !%Option medium_parallelepiped 1
703 !% The medium box will be a parallelepiped whose center and dimensions are taken from
704 !% the variable <tt>LinearMediumBoxSize</tt>.
705 !%Option medium_box_file 2
706 !% The simulation box will be read from an external file in OFF format, defined by the variable <tt>LinearMediumBoxFile</tt>.
707 !%End
708 call parse_variable(namespace, 'LinearMediumBoxShape', medium_parallelepiped, medium_box%box_shape)
709
710 if (.not. varinfo_valid_option('LinearMediumBoxShape', medium_box%box_shape)) then
711 call messages_input_error(namespace, 'LinearMediumBoxShape')
712 end if
713
714 select case (medium_box%box_shape)
716
717 !%Variable LinearMediumBoxSize
718 !%Type block
719 !%Section Maxwell::Medium
720 !%Description
721 !% Defines center and size of a parallelepiped linear medium box.
722 !%
723 !% Example:
724 !%
725 !% <tt>%LinearMediumBoxSize
726 !% <br>&nbsp;&nbsp; center_x | center_y | center_z | x_length | y_length | z_length
727 !% <br>%</tt>
728 !%End
729
730 if (parse_block(namespace, 'LinearMediumBoxSize', blk) == 0) then
731 call messages_print_with_emphasis(msg=trim('Linear Medium box center and size:'), namespace=namespace)
732
733 nlines = parse_block_n(blk)
734 if (nlines /= 1) then
735 call messages_input_error(namespace, 'LinearMediumBoxSize', 'should consist of one line')
736 end if
737 ncols = parse_block_cols(blk, 0)
738 if (ncols /= 6) then
739 call messages_input_error(namespace, 'LinearMediumBoxSize', 'should consist of six columns')
740 end if
741 do idim = 1, 3
742 call parse_block_float(blk, 0, idim-1, medium_box%center(idim))
743 call parse_block_float(blk, 0, idim+2, medium_box%lsize(idim))
744 end do
745 write(message(1),'(a,es9.2,a,es9.2,a,es9.2)') 'Box center: ', medium_box%center(1), ' | ',&
746 medium_box%center(2), ' | ', medium_box%center(3)
747 write(message(2),'(a,es9.2,a,es9.2,a,es9.2)') 'Box size : ', medium_box%lsize(1), ' | ', &
748 medium_box%lsize(2), ' | ', medium_box%lsize(3)
749 write(message(3),'(a)') ""
750 call messages_info(3)
751 call parse_block_end(blk)
752
753 call messages_print_with_emphasis(namespace=namespace)
754 else
755 message(1) = "For parallelepiped box shapes, you must provide a LinearMediumBoxSize block."
756 call messages_fatal(1, namespace=namespace)
757 end if
758
759 case (medium_box_file)
760 !%Variable LinearMediumBoxFile
761 !%Type string
762 !%Section Maxwell::Medium
763 !%Description
764 !% File in OFF format with the shape of the linear medium.
765 !%End
766 if (parse_is_defined(namespace, 'LinearMediumBoxFile')) then
767 call parse_variable(namespace, 'LinearMediumBoxFile', 'mediumboxfile', medium_box%filename)
768 else
769 message(1) = "When using box_file as the box shape, you must provide a filename through the LinearMediumBoxFile variable."
770 call messages_fatal(1, namespace=namespace)
771 end if
772
773 !%Variable CheckPointsMediumFromFile
774 !%Type logical
775 !%Default no
776 !%Section Maxwell::Medium
777 !%Description
778 !% Whether to re-calculate the points map by artificially shrinking the coordinate system by a factor of
779 !% 0.99 to check if the points inside the medium surface are properly detected. This works for only one
780 !% medium surface which is centered in the origin of the coordinate system.
781 !%End
782 call parse_variable(namespace, 'CheckPointsMediumFromFile', .false., medium_box%check_medium_points)
783
784 end select
785
786 pop_sub(medium_box_init)
787
788 end subroutine medium_box_init
789
790end module linear_medium_oct_m
791
792!! Local Variables:
793!! mode: f90
794!! coding: utf-8
795!! 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:189
real(real64), parameter, public m_zero
Definition: global.F90:187
character(len= *), parameter, public td_dir
Definition: global.F90:250
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:468
subroutine, public io_mkdir(fname, namespace, parents)
Definition: io.F90:354
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
Definition: io.F90:395
This module defines a linear medium for use in classical electrodynamics calculations.
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_initial_conditions(this)
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:930
subroutine, public messages_not_implemented(feature, namespace)
Definition: messages.F90:1125
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:543
subroutine, public messages_info(no_lines, iunit, verbose_limit, stress, all_nodes, namespace)
Definition: messages.F90:624
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:160
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:420
subroutine, public messages_input_error(namespace, var, details, row, column)
Definition: messages.F90:723
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:891
subroutine, public multicomm_init(mc, namespace, base_grp, mode_para, n_node, index_range, min_range)
create index and domain communicators
Definition: multicomm.F90:264
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 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:137
integer, parameter, public e_conductivity
Definition: quantity.F90:146
integer, parameter, public m_conductivity
Definition: quantity.F90:146
integer, parameter, public permittivity
Definition: quantity.F90:146
integer, parameter, public permeability
Definition: quantity.F90:146
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:1116
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
contains the information of the meshes and provides the transfer functions
Definition: regridding.F90:197
Abstract class for systems.
Definition: system.F90:172
int true(void)