Octopus
box_cgal.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-2006 M. Marques, A. Castro, A. Rubio, G. Bertsch
2!! Copyright (C) 2021 M. Oliveira, K. Lively, A. Obzhirov, I. Albar
3!! Copyright (C) 2022 S. Ohlmann
4!!
5!! This program is free software; you can redistribute it and/or modify
6!! it under the terms of the GNU General Public License as published by
7!! the Free Software Foundation; either version 2, or (at your option)
8!! any later version.
9!!
10!! This program is distributed in the hope that it will be useful,
11!! but WITHOUT ANY WARRANTY; without even the implied warranty of
12!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13!! GNU General Public License for more details.
14!!
15!! You should have received a copy of the GNU General Public License
16!! along with this program; if not, write to the Free Software
17!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
18!! 02110-1301, USA.
19!!
20
21#include "global.h"
22
23module box_cgal_oct_m
24 use box_oct_m
27 use debug_oct_m
28 use global_oct_m
31 use parser_oct_m
33 use string_oct_m
34 use unit_oct_m
35
36 implicit none
37
38 private
39 public :: box_cgal_t
40
42 type, extends(box_shape_t) :: box_cgal_t
43 private
44 type(cgal_polyhedra_t) :: cgal_poly
45 character(len=1024) :: filename
46 contains
47 procedure :: shape_contains_points => box_cgal_shape_contains_points
48 procedure :: write_info => box_cgal_write_info
49 procedure :: short_info => box_cgal_short_info
50 final :: box_cgal_finalize
51 end type box_cgal_t
52
53 interface box_cgal_t
54 procedure box_cgal_constructor
55 end interface box_cgal_t
56
57contains
58
59 !--------------------------------------------------------------
60 function box_cgal_constructor(dim, center, filename, length, namespace) result(box)
61 integer, intent(in) :: dim
62 real(real64), intent(in) :: center(dim)
63 character(len=*), intent(in) :: filename
64 real(real64), intent(in) :: length(dim)
65 type(namespace_t), intent(in) :: namespace
66 class(box_cgal_t), pointer :: box
67
68 push_sub(box_cgal_constructor)
69
70 if (dim /= 3) then
71 message(1) = "CGAL boxes can only be used in 3 dimensions."
72 call messages_fatal(1)
73 end if
74
75 ! Allocate memory
76 safe_allocate(box)
77
78 ! Initialize box
79 call box_shape_init(box, namespace, dim, center, bounding_box_min=-m_half*length, &
80 bounding_box_max=m_half*length)
81 ! initialize cgal part
82 box%filename = trim(filename)
83 call cgal_polyhedron_init(box%cgal_poly, trim(box%filename), verbose = .false.)
84
85 box%bounding_box_l = m_half*length + abs(center)
86
88 end function box_cgal_constructor
89
90 !--------------------------------------------------------------
91 subroutine box_cgal_finalize(this)
92 type(box_cgal_t), intent(inout) :: this
93
94 push_sub(box_cgal_finalize)
95
96 call cgal_polyhedron_end(this%cgal_poly)
97 call box_shape_end(this)
98
99 pop_sub(box_cgal_finalize)
100 end subroutine box_cgal_finalize
101
102 !--------------------------------------------------------------
103 function box_cgal_shape_contains_points(this, nn, xx, tol) result(contained)
104 class(box_cgal_t), intent(in) :: this
105 integer, intent(in) :: nn
106 real(real64), contiguous, intent(in) :: xx(:,:)
107 real(real64), optional, intent(in) :: tol
108 logical :: contained(1:nn)
109
110 integer :: ip
111
112 ! no push_sub/pop_sub, called too often
113
114 ip = 1
115 contained(:) = .false.
117 do ip = 1, nn
118 contained(ip) = cgal_polyhedron_point_inside(this%cgal_poly, xx(ip, 1), xx(ip, 2), xx(ip, 3)) &
119 .neqv. this%is_inside_out()
120 end do
121
123
124 !--------------------------------------------------------------
125 subroutine box_cgal_write_info(this, iunit, namespace)
126 class(box_cgal_t), intent(in) :: this
127 integer, optional, intent(in) :: iunit
128 type(namespace_t), optional, intent(in) :: namespace
129
130 push_sub(box_cgal_write_info)
131
132 write(message(1),'(2x,a)') 'Type = cgal'
133 write(message(2),'(2x,a,1x,a)') 'Loaded shape data from file', trim(this%filename)
134 call messages_info(2, iunit=iunit, namespace=namespace)
136 pop_sub(box_cgal_write_info)
137 end subroutine box_cgal_write_info
139 !--------------------------------------------------------------
140 character(len=BOX_INFO_LEN) function box_cgal_short_info(this, unit_length) result(info)
141 class(box_cgal_t), intent(in) :: this
142 type(unit_t), intent(in) :: unit_length
144 push_sub(box_cgal_short_info)
145
146 write(info,'(3a)') 'BoxShape = cgal; BoxCgalFile = ', trim(this%filename)
147
148 pop_sub(box_cgal_short_info)
149 end function box_cgal_short_info
150
151end module box_cgal_oct_m
152
153!! Local Variables:
154!! mode: f90
155!! coding: utf-8
156!! End:
subroutine info()
Definition: em_resp.F90:1096
logical function, dimension(1:nn) box_cgal_shape_contains_points(this, nn, xx, tol)
Definition: box_cgal.F90:197
class(box_cgal_t) function, pointer box_cgal_constructor(dim, center, filename, length, namespace)
Definition: box_cgal.F90:154
subroutine box_cgal_write_info(this, iunit, namespace)
Definition: box_cgal.F90:219
subroutine box_cgal_finalize(this)
Definition: box_cgal.F90:185
character(len=box_info_len) function box_cgal_short_info(this, unit_length)
Definition: box_cgal.F90:234
subroutine, public box_shape_init(this, namespace, dim, center, bounding_box_min, bounding_box_max, axes)
Definition: box_shape.F90:177
subroutine, public box_shape_end(this)
Definition: box_shape.F90:273
subroutine, public cgal_polyhedron_init(cgal_poly, fname, verbose)
logical function, public cgal_polyhedron_point_inside(cgal_poly, xq, yq, zq)
subroutine, public cgal_polyhedron_end(cgal_poly)
real(real64), parameter, public m_half
Definition: global.F90:194
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:414
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:616
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:132
Class implementing a box defined by a file read using the cgal library.
Definition: box_cgal.F90:135
Base class for more specialized boxes that are defined by a shape and have a center and basis vectors...
Definition: box_shape.F90:134