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 math_oct_m
42 use mesh_oct_m
43 use mpi_oct_m
48 use parser_oct_m
55 use space_oct_m
56 use system_oct_m
58
59 implicit none
60
61 private
62 public :: &
67
68 integer, parameter :: &
69 MEDIUM_PARALLELEPIPED = 1, &
71
74 type, extends(system_t) :: linear_medium_t
75 type(space_t) :: space
76 real(real64) :: ep_factor
77 real(real64) :: mu_factor
78 real(real64) :: sigma_e_factor
79 real(real64) :: sigma_m_factor
80 integer :: edge_profile
81 type(output_t) :: outp
82 type(grid_t) :: gr
83 type(multicomm_t) :: mc
84 type(single_medium_box_t) :: medium_box
85
86 contains
87 procedure :: init_interaction => linear_medium_init_interaction
88 procedure :: init_interaction_as_partner => linear_medium_init_interaction_as_partner
89 procedure :: initialize => linear_medium_initialize
90 procedure :: do_algorithmic_operation => linear_medium_do_algorithmic_operation
91 procedure :: is_tolerance_reached => linear_medium_is_tolerance_reached
92 procedure :: copy_quantities_to_interaction => linear_medium_copy_quantities_to_interaction
93 procedure :: restart_write_data => linear_medium_restart_write_data
94 procedure :: restart_read_data => linear_medium_restart_read_data
95 procedure :: update_kinetic_energy => linear_medium_update_kinetic_energy
97 end type linear_medium_t
98
99 interface linear_medium_t
100 procedure linear_medium_constructor
101 end interface linear_medium_t
102
103contains
104
105 ! ---------------------------------------------------------
110 function linear_medium_constructor(namespace) result(sys)
111 class(linear_medium_t), pointer :: sys
112 type(namespace_t), intent(in) :: namespace
113
115
116 allocate(sys)
117
118 call linear_medium_init(sys, namespace)
119
121 end function linear_medium_constructor
122
123 ! ---------------------------------------------------------
128 ! ---------------------------------------------------------
129 subroutine linear_medium_init(this, namespace)
130 class(linear_medium_t), target, intent(inout) :: this
131 type(namespace_t), intent(in) :: namespace
132
133 integer :: nlines, ncols,n_points, n_global_points
134 integer(int64) :: index_range(4)
135 integer, allocatable :: tmp(:)
136 type(block_t) :: blk
137
138 push_sub(linear_medium_init)
139
140 call profiling_in('MEDIUM_BOX_INIT')
141
142 this%namespace = namespace
143 this%space = space_t(this%namespace)
144 call grid_init_stage_1(this%gr, this%namespace, this%space, latt=lattice_vectors_t(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)
158
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
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
350
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 type(restart_basic_t) :: restart
369 integer :: restart_file_unit, ierr
370
372
373 call restart%init(this%namespace, restart_td, restart_type_dump, ierr)
374
375 message(1) = "Linear medium system "//trim(this%namespace%get())//" will only write a dummy restart file."
376 call messages_warning(1, namespace=this%namespace)
377
378 if (ierr==0) then
379 restart_file_unit = restart%open('restart_linear_medium')
380 if (restart%do_i_write()) write(restart_file_unit, *) 'Linear medium restart file'
381 call restart%close(restart_file_unit)
382 message(1) = "Successfully wrote restart data for system "//trim(this%namespace%get())
383 call messages_info(1, namespace=this%namespace)
384 else
385 message(1) = "Could not write restart data for system "//trim(this%namespace%get())
386 call messages_warning(1)
387 endif
388
389 call restart%end()
390
393
394 ! ---------------------------------------------------------
395 ! this function returns true if restart data could be read
396 logical function linear_medium_restart_read_data(this)
397 class(linear_medium_t), intent(inout) :: this
398
399 type(restart_basic_t) :: restart
400 integer :: restart_file_unit, ierr
401
403
404 call restart%init(this%namespace, restart_td, restart_type_load, ierr)
405 restart_file_unit = restart%open('restart_linear_medium')
406 if (ierr == 0) then
407 call restart%close(restart_file_unit)
409 else
410 ! could not open file
412 end if
413 call restart%end()
416 message(1) = "Successfully read restart data for system "//trim(this%namespace%get())
417 call messages_info(1, namespace=this%namespace)
418 end if
419
422
424 class(linear_medium_t), intent(inout) :: this
425
427
428 this%kinetic_energy = m_zero
429
431
433
434 ! ---------------------------------------------------------
435 subroutine linear_medium_finalize(this)
436 type(linear_medium_t), intent(inout) :: this
437
438 push_sub(linear_medium_finalize)
439
440 call single_medium_box_end(this%medium_box)
441 call system_end(this)
442 call multicomm_end(this%mc)
443 call grid_end(this%gr)
444
446 end subroutine linear_medium_finalize
447
448
449 ! Specific routines for this system:
450
451 ! ---------------------------------------------------------
452 subroutine get_medium_box_points_map(medium_box, gr)
453 type(single_medium_box_t), intent(inout) :: medium_box
454 type(grid_t), intent(in) :: gr
455
456 integer :: ip, ip_in, ip_bd, idim
457 integer, allocatable :: tmp_points_map(:), tmp_bdry_map(:)
458 real(real64) :: bounds(2,gr%box%dim), xx(gr%box%dim)
459 logical :: inside
462
463 call profiling_in('GET_MEDIUM_BOX_POINTS_MAP')
464
465 safe_allocate(tmp_points_map(1:gr%np))
466 safe_allocate(tmp_bdry_map(1:gr%np))
467 tmp_points_map = 0
468 tmp_bdry_map = 0
469
470 if (medium_box%box_shape == medium_box_file) then
471
472 call get_points_map_from_file(medium_box%filename, gr, medium_box%points_number,&
473 medium_box%global_points_number, tmp_points_map)
474
475 else
476
477 do idim = 1, 3
478 bounds(1,idim) = medium_box%center(idim) - medium_box%lsize(idim)/m_two
479 bounds(2,idim) = medium_box%center(idim) + medium_box%lsize(idim)/m_two
480 end do
481 ip_in = 0
482 ip_bd = 0
483 do ip = 1, gr%np
484 xx(1:3) = gr%x(ip, 1:3)
485 inside = check_point_in_bounds(xx, bounds(:,:))
486 if (check_point_in_bounds(xx, bounds(:,:))) then
487 ip_in = ip_in + 1
488 tmp_points_map(ip_in) = ip
489 end if
490 if (check_point_on_bounds(xx, bounds(:,:))) then
491 ip_bd = ip_bd + 1
492 tmp_bdry_map(ip_bd) = ip
493 end if
494 end do
495
496 medium_box%points_number = ip_in
497 medium_box%bdry_number = ip_bd
498
499 safe_allocate(medium_box%bdry_map(1:medium_box%bdry_number))
500 medium_box%bdry_map = 0
501 medium_box%bdry_map = tmp_bdry_map(1:medium_box%bdry_number)
502
503 end if
504 call single_medium_box_allocate(medium_box, medium_box%points_number)
505 medium_box%points_map(:) = tmp_points_map(1:medium_box%points_number)
506
507 ! TODO: add some check that medium boxes do not overlap (now they are different systems)
508 safe_deallocate_a(tmp_points_map)
509 safe_deallocate_a(tmp_bdry_map)
510
511 call profiling_out('GET_MEDIUM_BOX_POINTS_MAP')
512
514 contains
515
516 logical pure function check_point_in_bounds(xx, bounds) result (check)
517 real(real64), intent(in) :: xx(:)
518 real(real64), intent(in) :: bounds(:,:)
519
520 check = .false.
521 if ((xx(1) >= bounds(1,1)) .and. (xx(1) <= bounds(2,1)) .and. &
522 (xx(2) >= bounds(1,2)) .and. (xx(2) <= bounds(2,2)) .and. &
523 (xx(3) >= bounds(1,3)) .and. (xx(3) <= bounds(2,3)) ) then
524 check = .true.
525 end if
526
527 end function check_point_in_bounds
528
529 logical pure function check_point_on_bounds(xx, bounds) result (check)
530 real(real64), intent(in) :: xx(:)
531 real(real64), intent(in) :: bounds(:,:)
532
533 check = .false.
534 if (is_close(xx(1), bounds(1,1)) .and. (xx(2) >= bounds(1,2) .and. xx(3) >= bounds(1,3)) &
535 .and. (xx(2) <= bounds(2,2) .and. xx(3) <= bounds(2,3)) .or. &
536 is_close(xx(2), bounds(1,2)) .and. (xx(1) >= bounds(1,1) .and. xx(3) >= bounds(1,3)) &
537 .and. (xx(1) <= bounds(2,1) .and. xx(3) <= bounds(2,3)) .or. &
538 is_close(xx(3), bounds(1,3)) .and. (xx(1) >= bounds(1,1) .and. xx(2) >= bounds(1,2)) &
539 .and. (xx(1) <= bounds(2,1) .and. xx(2) <= bounds(2,2)) .or. &
540 is_close(xx(1), bounds(2,1)) .and. (xx(2) >= bounds(1,2) .and. xx(3) >= bounds(1,3)) &
541 .and. (xx(2) <= bounds(2,2) .and. xx(3) <= bounds(2,3)) .or. &
542 is_close(xx(2), bounds(2,2)) .and. (xx(1) >= bounds(1,1) .and. xx(3) >= bounds(1,3)) &
543 .and. (xx(1) <= bounds(2,1) .and. xx(3) <= bounds(2,3)) .or. &
544 is_close(xx(3), bounds(2,3)) .and. (xx(1) >= bounds(1,1) .and. xx(2) >= bounds(1,2)) &
545 .and. (xx(1) <= bounds(2,1) .and. xx(2) <= bounds(2,2)) ) then
546 check = .true.
547 end if
548
549 end function check_point_on_bounds
550
551 end subroutine get_medium_box_points_map
552
553 ! ----------------------------------------------------------
555 subroutine get_linear_medium_em_properties(this, medium_box, gr)
556 type(linear_medium_t), intent(in) :: this
557 type(single_medium_box_t), intent(inout) :: medium_box
558 type(grid_t), intent(in) :: gr
559
560 integer :: ip, ip_in, ip_bd, ipp
561 real(real64) :: xx(gr%box%dim), xxp(gr%box%dim), dd, dd_max, dd_min
562 real(real64), allocatable :: tmp(:), tmp_grad(:,:)
563
565
566 call profiling_in('GET_LIN_MEDIUM_EM_PROP')
567
568 safe_allocate(tmp(1:gr%np_part))
569 safe_allocate(tmp_grad(1:gr%np_part, 1:gr%box%dim))
570 dd_max = max(2*gr%spacing(1), 2*gr%spacing(2), 2*gr%spacing(3))
571
572 do ip_in = 1, medium_box%points_number
573 ip = medium_box%points_map(ip_in)
574 if (this%edge_profile == option__linearmediumedgeprofile__smooth) then
575 assert(allocated(medium_box%bdry_map))
576
577 xx(1:3) = gr%x(ip,1:3)
578 dd_min = m_huge
579
580 do ip_bd = 1, medium_box%bdry_number
581 ipp = medium_box%bdry_map(ip_bd)
582 xxp(1:3) = gr%x(ipp,1:3)
583 dd = norm2(xx(1:3) - xxp(1:3))
584 if (dd < dd_min) dd_min = dd
585 end do
586
587 medium_box%ep(ip_in) = p_ep + ((p_ep * this%ep_factor - p_ep) &
588 * m_one/(m_one + exp(-m_five/dd_max * (dd_min - m_two*dd_max))))
589 medium_box%mu(ip_in) = p_mu + ((p_mu * this%mu_factor - p_mu) &
590 * m_one/(m_one + exp(-m_five/dd_max * (dd_min - m_two*dd_max))))
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 * m_one/(m_one + exp(-m_five/dd_max * (dd_min - m_two*dd_max)))
594 medium_box%sigma_m(ip_in) = this%sigma_m_factor &
595 * m_one/(m_one + exp(-m_five/dd_max * (dd_min - m_two*dd_max)))
596
597 else if (this%edge_profile == option__linearmediumedgeprofile__edged) then
598
599 medium_box%ep(ip_in) = p_ep * this%ep_factor
600 medium_box%mu(ip_in) = p_mu * this%mu_factor
601 medium_box%c(ip_in) = m_one/sqrt(medium_box%ep(ip_in)*medium_box%mu(ip_in))
602 medium_box%sigma_e(ip_in) = this%sigma_e_factor
603 medium_box%sigma_m(ip_in) = this%sigma_m_factor
604
605 end if
606 end do
607
608 tmp(:) = p_ep
609 do ip_in = 1, medium_box%points_number
610 ip = medium_box%points_map(ip_in)
611 tmp(ip)= medium_box%ep(ip_in)
612 end do
613 call dderivatives_grad(gr%der, tmp, tmp_grad, set_bc = .false.)
614 do ip_in = 1, medium_box%points_number
615 ip = medium_box%points_map(ip_in)
616 medium_box%aux_ep(ip_in, :) = tmp_grad(ip, :)/(m_four * medium_box%ep(ip_in))
617 end do
618
619 tmp(:) = p_mu
620 do ip_in = 1, medium_box%points_number
621 ip = medium_box%points_map(ip_in)
622 tmp(ip) = medium_box%mu(ip_in)
623 end do
624 call dderivatives_grad(gr%der, tmp, tmp_grad, set_bc = .false.)
625 do ip_in = 1, medium_box%points_number
626 ip = medium_box%points_map(ip_in)
627 medium_box%aux_mu(ip_in, :) = tmp_grad(ip, :)/(m_four * medium_box%mu(ip_in))
628 end do
629
630 !TODO: add print information about the medium box
631
632 safe_deallocate_a(tmp)
633 safe_deallocate_a(tmp_grad)
634
635 call profiling_out('GET_LIN_MEDIUM_EM_PROP')
636
639
640 ! ----------------------------------------------------------
642 subroutine get_points_map_from_file(filename, mesh, n_points, global_points_number, tmp_map, scale_factor)
643 character(len=256), intent(in) :: filename
644 class(mesh_t), intent(in) :: mesh
645 integer, intent(out) :: n_points
646 integer, intent(out) :: global_points_number
647 integer, intent(inout) :: tmp_map(:)
648 real(real64), optional, intent(in) :: scale_factor
649
650 integer :: ip_in, ip
651 real(real64) :: xx(3)
652 type(cgal_polyhedra_t) :: cgal_poly
653
655
656 call profiling_in('GET_POINTS_MAP_FROM_FILE')
657
658 call cgal_polyhedron_init(cgal_poly, trim(filename), verbose = .false.)
659
660 ip_in = 0
661 do ip = 1, mesh%np
662 if (present(scale_factor)) then
663 xx(1:3) = scale_factor * mesh%x(ip, 1:3)
664 else
665 xx(1:3) = mesh%x(ip, 1:3)
666 end if
667 if (cgal_polyhedron_point_inside(cgal_poly, xx(1), xx(2), xx(3))) then
668 ip_in = ip_in + 1
669 tmp_map(ip_in) = ip
670 end if
671 end do
672 n_points = ip_in
673 call cgal_polyhedron_end(cgal_poly)
674
675 call mpi_world%allreduce(ip_in, global_points_number, 1, mpi_integer, mpi_sum)
676
677 call profiling_out('GET_POINTS_MAP_FROM_FILE')
678
680
681 end subroutine get_points_map_from_file
682
683 ! ----------------------------------------------------------
685 subroutine medium_box_init(medium_box, namespace)
686 type(single_medium_box_t), intent(inout) :: medium_box
687 type(namespace_t), intent(in) :: namespace
688
689 integer :: nlines, ncols, idim
690 type(block_t) :: blk
691
692 push_sub(medium_box_init)
693
694 !%Variable LinearMediumBoxShape
695 !%Type integer
696 !%Section Maxwell::Medium
697 !%Description
698 !% This variable defines the shape of the linear medium box.
699 !% The default is <tt>medium_parallelepiped</tt>.
700 !%Option medium_parallelepiped 1
701 !% The medium box will be a parallelepiped whose center and dimensions are taken from
702 !% the variable <tt>LinearMediumBoxSize</tt>.
703 !%Option medium_box_file 2
704 !% The simulation box will be read from an external file in OFF format, defined by the variable <tt>LinearMediumBoxFile</tt>.
705 !%End
706 call parse_variable(namespace, 'LinearMediumBoxShape', medium_parallelepiped, medium_box%box_shape)
707
708 if (.not. varinfo_valid_option('LinearMediumBoxShape', medium_box%box_shape)) then
709 call messages_input_error(namespace, 'LinearMediumBoxShape')
710 end if
711
712 select case (medium_box%box_shape)
714
715 !%Variable LinearMediumBoxSize
716 !%Type block
717 !%Section Maxwell::Medium
718 !%Description
719 !% Defines center and size of a parallelepiped linear medium box.
720 !%
721 !% Example:
722 !%
723 !% <tt>%LinearMediumBoxSize
724 !% <br>&nbsp;&nbsp; center_x | center_y | center_z | x_length | y_length | z_length
725 !% <br>%</tt>
726 !%End
727
728 if (parse_block(namespace, 'LinearMediumBoxSize', blk) == 0) then
729 call messages_print_with_emphasis(msg=trim('Linear Medium box center and size:'), namespace=namespace)
730
731 nlines = parse_block_n(blk)
732 if (nlines /= 1) then
733 call messages_input_error(namespace, 'LinearMediumBoxSize', 'should consist of one line')
734 end if
735 ncols = parse_block_cols(blk, 0)
736 if (ncols /= 6) then
737 call messages_input_error(namespace, 'LinearMediumBoxSize', 'should consist of six columns')
738 end if
739 do idim = 1, 3
740 call parse_block_float(blk, 0, idim-1, medium_box%center(idim))
741 call parse_block_float(blk, 0, idim+2, medium_box%lsize(idim))
742 end do
743 write(message(1),'(a,es9.2,a,es9.2,a,es9.2)') 'Box center: ', medium_box%center(1), ' | ',&
744 medium_box%center(2), ' | ', medium_box%center(3)
745 write(message(2),'(a,es9.2,a,es9.2,a,es9.2)') 'Box size : ', medium_box%lsize(1), ' | ', &
746 medium_box%lsize(2), ' | ', medium_box%lsize(3)
747 write(message(3),'(a)') ""
748 call messages_info(3)
749 call parse_block_end(blk)
750
751 call messages_print_with_emphasis(namespace=namespace)
752 else
753 message(1) = "For parallelepiped box shapes, you must provide a LinearMediumBoxSize block."
754 call messages_fatal(1, namespace=namespace)
755 end if
756
757 case (medium_box_file)
758 !%Variable LinearMediumBoxFile
759 !%Type string
760 !%Section Maxwell::Medium
761 !%Description
762 !% File in OFF format with the shape of the linear medium.
763 !%End
764 if (parse_is_defined(namespace, 'LinearMediumBoxFile')) then
765 call parse_variable(namespace, 'LinearMediumBoxFile', 'mediumboxfile', medium_box%filename)
766 else
767 message(1) = "When using box_file as the box shape, you must provide a filename through the LinearMediumBoxFile variable."
768 call messages_fatal(1, namespace=namespace)
769 end if
770
771 !%Variable CheckPointsMediumFromFile
772 !%Type logical
773 !%Default no
774 !%Section Maxwell::Medium
775 !%Description
776 !% Whether to re-calculate the points map by artificially shrinking the coordinate system by a factor of
777 !% 0.99 to check if the points inside the medium surface are properly detected. This works for only one
778 !% medium surface which is centered in the origin of the coordinate system.
779 !%End
780 call parse_variable(namespace, 'CheckPointsMediumFromFile', .false., medium_box%check_medium_points)
781
782 end select
783
784 pop_sub(medium_box_init)
785
786 end subroutine medium_box_init
787
788end module linear_medium_oct_m
789
790!! Local Variables:
791!! mode: f90
792!! coding: utf-8
793!! 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:143
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:192
real(real64), parameter, public m_zero
Definition: global.F90:190
This module implements the underlying real-space grid.
Definition: grid.F90:119
subroutine, public grid_init_stage_1(gr, namespace, space, symm, latt, n_sites, site_position)
First stage of the grid initialization.
Definition: grid.F90:197
subroutine, public grid_init_stage_2(gr, namespace, space, mc, qvector)
Second stage of the grid initialization.
Definition: grid.F90:473
subroutine, public grid_end(gr)
finalize a grid object
Definition: grid.F90:501
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:116
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:117
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:120
subroutine, public messages_print_with_emphasis(msg, iunit, namespace)
Definition: messages.F90:904
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:531
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:162
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:416
subroutine, public messages_input_error(namespace, var, details, row, column)
Definition: messages.F90:697
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:600
type(mpi_grp_t), public mpi_world
Definition: mpi.F90:272
This module handles the communicators for the various parallelization strategies.
Definition: multicomm.F90:147
subroutine, public multicomm_end(mc)
Definition: multicomm.F90:697
subroutine, public multicomm_init(mc, namespace, base_grp, mode_para, n_node, index_range, min_range)
create index and domain communicators
Definition: multicomm.F90:266
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:117
integer function, public parse_block(namespace, name, blk, check_varinfo_)
Definition: parser.F90:621
subroutine, public profiling_out(label)
Increment out counter and sum up difference between entry and exit time.
Definition: profiling.F90:625
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
Definition: profiling.F90:554
This module implements the basic propagator framework.
Definition: propagator.F90:119
This module defines the quantity_t class and the IDs for quantities, which can be exposed by a system...
Definition: quantity.F90:140
Implementation details for regridding.
Definition: regridding.F90:172
integer, parameter, public restart_type_dump
Definition: restart.F90:183
integer, parameter, public restart_td
Definition: restart.F90:156
integer, parameter, public restart_type_load
Definition: restart.F90:183
This module implements the abstract system type.
Definition: system.F90:120
subroutine, public system_end(this)
Definition: system.F90:1149
Descriptor of one algorithmic operation.
Definition: algorithm.F90:165
Description of the grid, containing information on derivatives, stencil, and symmetries.
Definition: grid.F90:171
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:173
contains the information of the meshes and provides the transfer functions
Definition: regridding.F90:198
restart_basic_t stores the basic information about a restart object.
Definition: restart.F90:216
Abstract class for systems.
Definition: system.F90:174
int true(void)