Neko 0.9.99
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
tet_mesh.f90
Go to the documentation of this file.
1! Copyright (c) 2021, 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!
36 use mesh
37 use tet
38 use point
39 use utils
40 implicit none
41 private
42
43 integer, public, parameter :: tet_msh_otpv = 1, tet_msh_fvtc = 2, &
44 tet_msh_svtc = 3
45
46 type, public :: tet_mesh_t
47 type(tet_t), allocatable :: el(:)
48 type(mesh_t), pointer :: msh
49 integer :: nelv
50 contains
51 procedure, pass(this) :: init => tet_mesh_init
52 procedure, pass(this) :: free => tet_mesh_free
53 end type tet_mesh_t
54
55contains
56
58 subroutine tet_mesh_init(this, msh, mthd)
59 class(tet_mesh_t), intent(inout) :: this
60 type(mesh_t), intent(in), target :: msh
61 integer, intent(in), optional :: mthd
62 integer :: bsct_mthd
63
64 call this%free()
65
66 this%msh => msh
67
68 if (present(mthd)) then
69 bsct_mthd = mthd
70 else
71 bsct_mthd = tet_msh_svtc
72 end if
73
74 if (bsct_mthd .eq. tet_msh_otpv) then
75 this%nelv = msh%nelv * 8
76 allocate(this%el(this%nelv))
77 call tet_mesh_bisect_otpv(this)
78 else if (bsct_mthd .eq. tet_msh_fvtc) then
79 this%nelv = msh%nelv * 5
80 allocate(this%el(this%nelv))
81 call tet_mesh_bisect_fvtc(this)
82 else if (bsct_mthd .eq. tet_msh_svtc) then
83 this%nelv = msh%nelv * 6
84 allocate(this%el(this%nelv))
85 call tet_mesh_bisect_svtc(this)
86 else
87 call neko_error('Invalid bisection strategy')
88 end if
89
90 end subroutine tet_mesh_init
91
93 subroutine tet_mesh_free(this)
94 class(tet_mesh_t), intent(inout) :: this
95
96 if (allocated(this%el)) then
97 deallocate(this%el)
98 end if
99
100 nullify(this%msh)
101
102 end subroutine tet_mesh_free
103
104 ! Bisect hexahedral mesh into a tetrahedral mesh using
105 ! one tet per vertex as described in,
106 ! P. D. Bello-Maldonado and P. F. Fischer
107 ! SIAM J. Sci. Comput., 41(5), S2–S18, 2019.
108 subroutine tet_mesh_bisect_otpv(tet_msh)
109 type(tet_mesh_t), intent(inout) :: tet_msh
110 integer :: i, j
111 type(point_t), pointer :: p1, p2, p3, p4
112
113 j = 0
114 do i = 1, tet_msh%msh%nelv
115
116 j = j + 1
117 p1 => tet_msh%msh%elements(i)%e%pts(1)%p
118 p2 => tet_msh%msh%elements(i)%e%pts(2)%p
119 p3 => tet_msh%msh%elements(i)%e%pts(5)%p
120 p4 => tet_msh%msh%elements(i)%e%pts(4)%p
121 call tet_msh%el(j)%init(j, p1, p2, p3, p4)
122
123 j = j + 1
124 p1 => tet_msh%msh%elements(i)%e%pts(1)%p
125 p2 => tet_msh%msh%elements(i)%e%pts(2)%p
126 p3 => tet_msh%msh%elements(i)%e%pts(3)%p
127 p4 => tet_msh%msh%elements(i)%e%pts(6)%p
128 call tet_msh%el(j)%init(j, p1, p2, p3, p4)
129
130 j = j + 1
131 p1 => tet_msh%msh%elements(i)%e%pts(4)%p
132 p2 => tet_msh%msh%elements(i)%e%pts(3)%p
133 p3 => tet_msh%msh%elements(i)%e%pts(2)%p
134 p4 => tet_msh%msh%elements(i)%e%pts(7)%p
135 call tet_msh%el(j)%init(j, p1, p2, p3, p4)
136
137 j = j + 1
138 p1 => tet_msh%msh%elements(i)%e%pts(4)%p
139 p2 => tet_msh%msh%elements(i)%e%pts(3)%p
140 p3 => tet_msh%msh%elements(i)%e%pts(8)%p
141 p4 => tet_msh%msh%elements(i)%e%pts(1)%p
142 call tet_msh%el(j)%init(j, p1, p2, p3, p4)
143
144 j = j + 1
145 p1 => tet_msh%msh%elements(i)%e%pts(1)%p
146 p2 => tet_msh%msh%elements(i)%e%pts(5)%p
147 p3 => tet_msh%msh%elements(i)%e%pts(6)%p
148 p4 => tet_msh%msh%elements(i)%e%pts(8)%p
149 call tet_msh%el(j)%init(j, p1, p2, p3, p4)
150
151 j = j + 1
152 p1 => tet_msh%msh%elements(i)%e%pts(2)%p
153 p2 => tet_msh%msh%elements(i)%e%pts(5)%p
154 p3 => tet_msh%msh%elements(i)%e%pts(6)%p
155 p4 => tet_msh%msh%elements(i)%e%pts(7)%p
156 call tet_msh%el(j)%init(j, p1, p2, p3, p4)
157
158 j = j + 1
159 p1 => tet_msh%msh%elements(i)%e%pts(3)%p
160 p2 => tet_msh%msh%elements(i)%e%pts(8)%p
161 p3 => tet_msh%msh%elements(i)%e%pts(7)%p
162 p4 => tet_msh%msh%elements(i)%e%pts(6)%p
163 call tet_msh%el(j)%init(j, p1, p2, p3, p4)
164
165 j = j + 1
166 p1 => tet_msh%msh%elements(i)%e%pts(4)%p
167 p2 => tet_msh%msh%elements(i)%e%pts(8)%p
168 p3 => tet_msh%msh%elements(i)%e%pts(7)%p
169 p4 => tet_msh%msh%elements(i)%e%pts(5)%p
170 call tet_msh%el(j)%init(j, p1, p2, p3, p4)
171
172 end do
173
174 end subroutine tet_mesh_bisect_otpv
175
177 subroutine tet_mesh_bisect_fvtc(tet_msh)
178 type(tet_mesh_t), intent(inout) :: tet_msh
179 integer :: i, j
180 type(point_t), pointer :: p1, p2, p3, p4
181
182 j = 0
183 do i = 1, tet_msh%msh%nelv
184
185 j = j + 1
186 p1 => tet_msh%msh%elements(i)%e%pts(1)%p
187 p2 => tet_msh%msh%elements(i)%e%pts(2)%p
188 p3 => tet_msh%msh%elements(i)%e%pts(3)%p
189 p4 => tet_msh%msh%elements(i)%e%pts(6)%p
190 call tet_msh%el(j)%init(j, p1, p2, p3, p4)
191
192 j = j + 1
193 p1 => tet_msh%msh%elements(i)%e%pts(1)%p
194 p2 => tet_msh%msh%elements(i)%e%pts(4)%p
195 p3 => tet_msh%msh%elements(i)%e%pts(3)%p
196 p4 => tet_msh%msh%elements(i)%e%pts(8)%p
197 call tet_msh%el(j)%init(j, p1, p2, p3, p4)
198
199 j = j + 1
200 p1 => tet_msh%msh%elements(i)%e%pts(3)%p
201 p2 => tet_msh%msh%elements(i)%e%pts(8)%p
202 p3 => tet_msh%msh%elements(i)%e%pts(7)%p
203 p4 => tet_msh%msh%elements(i)%e%pts(6)%p
204 call tet_msh%el(j)%init(j, p1, p2, p3, p4)
205
206 j = j + 1
207 p1 => tet_msh%msh%elements(i)%e%pts(1)%p
208 p2 => tet_msh%msh%elements(i)%e%pts(5)%p
209 p3 => tet_msh%msh%elements(i)%e%pts(8)%p
210 p4 => tet_msh%msh%elements(i)%e%pts(6)%p
211 call tet_msh%el(j)%init(j, p1, p2, p3, p4)
212
213 j = j + 1
214 p1 => tet_msh%msh%elements(i)%e%pts(8)%p
215 p2 => tet_msh%msh%elements(i)%e%pts(3)%p
216 p3 => tet_msh%msh%elements(i)%e%pts(1)%p
217 p4 => tet_msh%msh%elements(i)%e%pts(6)%p
218 call tet_msh%el(j)%init(j, p1, p2, p3, p4)
219 end do
220 end subroutine tet_mesh_bisect_fvtc
221
223 subroutine tet_mesh_bisect_svtc(tet_msh)
224 type(tet_mesh_t), intent(inout) :: tet_msh
225 integer :: i, j
226 type(point_t), pointer :: p1, p2, p3, p4
227
228 j = 0
229 do i = 1, tet_msh%msh%nelv
230
231 j = j + 1
232 p1 => tet_msh%msh%elements(i)%e%pts(5)%p
233 p2 => tet_msh%msh%elements(i)%e%pts(8)%p
234 p3 => tet_msh%msh%elements(i)%e%pts(1)%p
235 p4 => tet_msh%msh%elements(i)%e%pts(7)%p
236 call tet_msh%el(j)%init(j, p1, p2, p3, p4)
237
238 j = j + 1
239 p1 => tet_msh%msh%elements(i)%e%pts(5)%p
240 p2 => tet_msh%msh%elements(i)%e%pts(1)%p
241 p3 => tet_msh%msh%elements(i)%e%pts(7)%p
242 p4 => tet_msh%msh%elements(i)%e%pts(6)%p
243 call tet_msh%el(j)%init(j, p1, p2, p3, p4)
244
245 j = j + 1
246 p1 => tet_msh%msh%elements(i)%e%pts(8)%p
247 p2 => tet_msh%msh%elements(i)%e%pts(4)%p
248 p3 => tet_msh%msh%elements(i)%e%pts(7)%p
249 p4 => tet_msh%msh%elements(i)%e%pts(1)%p
250 call tet_msh%el(j)%init(j, p1, p2, p3, p4)
251
252 j = j + 1
253 p1 => tet_msh%msh%elements(i)%e%pts(1)%p
254 p2 => tet_msh%msh%elements(i)%e%pts(7)%p
255 p3 => tet_msh%msh%elements(i)%e%pts(2)%p
256 p4 => tet_msh%msh%elements(i)%e%pts(6)%p
257 call tet_msh%el(j)%init(j, p1, p2, p3, p4)
258
259 j = j + 1
260 p1 => tet_msh%msh%elements(i)%e%pts(1)%p
261 p2 => tet_msh%msh%elements(i)%e%pts(2)%p
262 p3 => tet_msh%msh%elements(i)%e%pts(7)%p
263 p4 => tet_msh%msh%elements(i)%e%pts(3)%p
264 call tet_msh%el(j)%init(j, p1, p2, p3, p4)
265
266 j = j + 1
267 p1 => tet_msh%msh%elements(i)%e%pts(1)%p
268 p2 => tet_msh%msh%elements(i)%e%pts(7)%p
269 p3 => tet_msh%msh%elements(i)%e%pts(3)%p
270 p4 => tet_msh%msh%elements(i)%e%pts(4)%p
271 call tet_msh%el(j)%init(j, p1, p2, p3, p4)
272
273 end do
274
275
276 end subroutine tet_mesh_bisect_svtc
277
278end module tet_mesh
Defines a mesh.
Definition mesh.f90:34
Implements a point.
Definition point.f90:35
Defines a tetrahedral mesh.
Definition tet_mesh.f90:35
subroutine tet_mesh_bisect_svtc(tet_msh)
Bisect each hexahedron into six tetrahedrons.
Definition tet_mesh.f90:224
subroutine tet_mesh_bisect_otpv(tet_msh)
Definition tet_mesh.f90:109
subroutine tet_mesh_bisect_fvtc(tet_msh)
Bisect each hexahedron into five tetrahedrons.
Definition tet_mesh.f90:178
integer, parameter, public tet_msh_otpv
Definition tet_mesh.f90:43
subroutine tet_mesh_init(this, msh, mthd)
Initialise a tetrahedral mesh based on a hexahedral mesh msh.
Definition tet_mesh.f90:59
subroutine tet_mesh_free(this)
Deallocate a tetrahedral mesh.
Definition tet_mesh.f90:94
integer, parameter, public tet_msh_fvtc
Definition tet_mesh.f90:43
integer, parameter, public tet_msh_svtc
Definition tet_mesh.f90:43
Defines a tetrahedral element.
Definition tet.f90:34
Utilities.
Definition utils.f90:35
A point in with coordinates .
Definition point.f90:43
Tetrahedral element.
Definition tet.f90:63