Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
stl_file.f90
Go to the documentation of this file.
1! Copyright (c) 2022-2025, The Neko Authors
2! All rights reserved.
3!
4! Redistribution and use in source and binary forms, with or without
5! modification, are permitted provided that the following conditions
6! are met:
7!
8! * Redistributions of source code must retain the above copyright
9! notice, this list of conditions and the following disclaimer.
10!
11! * Redistributions in binary form must reproduce the above
12! copyright notice, this list of conditions and the following
13! disclaimer in the documentation and/or other materials provided
14! with the distribution.
15!
16! * Neither the name of the authors nor the names of its
17! contributors may be used to endorse or promote products derived
18! from this software without specific prior written permission.
19!
20! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
21! "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
22! LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
23! FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE
24! COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT,
25! INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING,
26! BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
27! LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
28! CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
29! LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN
30! ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
31! POSSIBILITY OF SUCH DAMAGE.
32!
35 use num_types, only : rp
37 use tri_mesh, only : tri_mesh_t
38 use logger, only : neko_log
39 use point, only : point_t
41 use mpi_f08, only : mpi_mode_rdonly, mpi_info_null, &
42 mpi_file_open, mpi_file_read_all, mpi_file, mpi_status, &
43 mpi_file_close
44 use comm, only : neko_comm
45 use stl, only : stl_hdr_t, stl_triangle_t
46 implicit none
47 private
48
50 type, public, extends(generic_file_t) :: stl_file_t
51 contains
52 procedure :: read => stl_file_read
53 procedure :: write => stl_file_write
54 end type stl_file_t
55
56contains
57
58 subroutine stl_file_write(this, data, t)
59 class(stl_file_t), intent(inout) :: this
60 class(*), target, intent(in) :: data
61 real(kind=rp), intent(in), optional :: t
62 call neko_log%error('Not implemented')
63 end subroutine stl_file_write
64
65 subroutine stl_file_read(this, data)
66 class(stl_file_t) :: this
67 class(*), target, intent(inout) :: data
68 type(tri_mesh_t), pointer :: tri_msh => null()
69 type(mpi_status) :: status
70 type(mpi_file) :: fh
71 type(point_t), target :: p1, p2, p3
72 type(stl_hdr_t) :: stl_hdr
73 type(stl_triangle_t), allocatable :: stl_tri(:)
74 integer :: i, p_idx, ierr
75
76 call this%check_exists()
77
78 select type (data)
79 type is (tri_mesh_t)
80 tri_msh => data
81 class default
82 call neko_log%error('Invalid data')
83 end select
84
85 call mpi_file_open(neko_comm, trim(this%get_fname()), &
86 mpi_mode_rdonly, mpi_info_null, fh, ierr)
87 call mpi_file_read_all(fh, stl_hdr, 1, mpi_stl_header, status, ierr)
88
89 if (stl_hdr%hdr(1:6) .eq. 'solid') then
90 call neko_log%error('Invalid STL file (ASCII)')
91 end if
92
93 call tri_msh%init(stl_hdr%ntri)
94 allocate(stl_tri(stl_hdr%ntri))
95
96 call mpi_file_read_all(fh, stl_tri, stl_hdr%ntri, &
97 mpi_stl_triangle, status, ierr)
98
99 p_idx = 0
100 do i = 1, stl_hdr%ntri
101 p_idx = p_idx + 1
102 call p1%init(dble(stl_tri(i)%v1), p_idx)
103 p_idx = p_idx + 1
104 call p2%init(dble(stl_tri(i)%v2), p_idx)
105 p_idx = p_idx + 1
106 call p3%init(dble(stl_tri(i)%v3), p_idx)
107 call tri_msh%add_element(p1, p2, p3)
108 end do
109
110 deallocate(stl_tri)
111
112 call mpi_file_close(fh, ierr)
113
114 end subroutine stl_file_read
115
116end module stl_file
Definition comm.F90:1
type(mpi_comm), public neko_comm
MPI communicator.
Definition comm.F90:43
Logging routines.
Definition log.f90:34
type(log_t), public neko_log
Global log stream.
Definition log.f90:79
MPI derived types.
Definition mpi_types.f90:34
type(mpi_datatype), public mpi_stl_header
MPI Derived type for a STL header.
Definition mpi_types.f90:62
type(mpi_datatype), public mpi_stl_triangle
MPI derived type for a STL triangle.
Definition mpi_types.f90:63
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Implements a point.
Definition point.f90:35
Stereolithography (STL) file.
Definition stl_file.f90:34
subroutine stl_file_write(this, data, t)
Definition stl_file.f90:59
subroutine stl_file_read(this, data)
Definition stl_file.f90:66
Stereolithography format.
Definition stl.f90:2
Defines a triangular surface mesh.
Definition tri_mesh.f90:35
A generic file handler.
A point in with coordinates .
Definition point.f90:43
Defines a STL hdr.
Definition stl.f90:8
Defines a STL triangle.
Definition stl.f90:14
Interface for STL files.
Definition stl_file.f90:50