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