Neko  0.8.1
A portable framework for high-order spectral element flow simulations
hex.f90
Go to the documentation of this file.
1 ! Copyright (c) 2019-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 !
34 module hex
35  use num_types, only : dp
36  use element, only : element_t
37  use tuple, only : tuple_t, tuple_i4_t, tuple4_i4_t
38  use point, only : point_t
39  implicit none
40  private
41 
42  integer, public, parameter :: neko_hex_npts = 8
43  integer, public, parameter :: neko_hex_nfcs = 6
44  integer, public, parameter :: neko_hex_neds = 12
45  integer, public, parameter :: neko_hex_gdim = 3
46 
47 
63  type, public, extends(element_t) :: hex_t
64  contains
65  procedure, pass(this) :: init => hex_init
66  procedure, pass(this) :: facet_id => hex_facet_id
67  procedure, pass(this) :: facet_order => hex_facet_order
68  procedure, pass(this) :: diameter => hex_diameter
69  procedure, pass(this) :: centroid => hex_centroid
70  procedure, pass(this) :: edge_id => hex_edge_id
71  procedure, pass(this) :: equal => hex_equal
72  generic :: operator(.eq.) => equal
73  end type hex_t
74 
93  integer, parameter, dimension(4, 6) :: face_nodes = reshape((/1,5,7,3,&
94  2,6,8,4,&
95  1,2,6,5,&
96  3,4,8,7,&
97  1,2,4,3,&
98  5,6,8,7/),&
99  (/4,6/))
100 
119  integer, parameter, dimension(2, 12) :: edge_nodes = reshape((/1,2,&
120  3,4,&
121  5,6,&
122  7,8,&
123  1,3,&
124  2,4,&
125  5,7,&
126  6,8,&
127  1,5,&
128  2,6,&
129  3,7,&
130  4,8/),&
131  (/2,12/))
132 
133 contains
134 
136  subroutine hex_init(this, id, p1, p2, p3, p4, p5, p6, p7, p8)
137  class(hex_t), intent(inout) :: this
138  integer, intent(inout) :: id
139  type(point_t), target, intent(in) :: p1, p2, p3, p4, p5, p6, p7, p8
140 
141  call this%element(id, neko_hex_gdim, neko_hex_npts)
142 
143  this%pts(1)%p => p1
144  this%pts(2)%p => p2
145  this%pts(3)%p => p3
146  this%pts(4)%p => p4
147  this%pts(5)%p => p5
148  this%pts(6)%p => p6
149  this%pts(7)%p => p7
150  this%pts(8)%p => p8
151 
152  end subroutine hex_init
153 
155  subroutine hex_facet_id(this, t, side)
156  class(hex_t), intent(in) :: this
157  class(tuple_t), intent(inout) :: t
158  integer, intent(in) :: side
159  integer :: i, j, temp
160  type(point_t), pointer :: p1,p2,p3,p4
161 
162  p1 => this%p(face_nodes(1, side))
163  p2 => this%p(face_nodes(2, side))
164  p3 => this%p(face_nodes(3, side))
165  p4 => this%p(face_nodes(4, side))
166 
167  select type(t)
168  type is(tuple4_i4_t)
169  t%x = (/ p1%id(), p2%id(), p3%id(), p4%id() /)
170  do i = 1, 3
171  do j = i+1,4
172  if(t%x(j) .lt. t%x(i)) then
173  temp = t%x(i)
174  t%x(i) = t%x(j)
175  t%x(j) = temp
176  endif
177  enddo
178  enddo
179  end select
180 
181  end subroutine hex_facet_id
182 
184  subroutine hex_facet_order(this, t, side)
185  class(hex_t), intent(in) :: this
186  class(tuple_t), intent(inout) :: t
187  integer, intent(in) :: side
188  type(point_t), pointer :: p1,p2,p3,p4
189 
190  p1 => this%p(face_nodes(1, side))
191  p2 => this%p(face_nodes(2, side))
192  p3 => this%p(face_nodes(3, side))
193  p4 => this%p(face_nodes(4, side))
194 
195  select type(t)
196  type is(tuple4_i4_t)
197  t%x = (/ p1%id(), p2%id(), p3%id(), p4%id() /)
198  end select
199 
200  end subroutine hex_facet_order
201 
202 
204  subroutine hex_edge_id(this, t, side)
205  class(hex_t), intent(in) :: this
206  class(tuple_t), intent(inout) :: t
207  integer, intent(in) :: side
208  type(point_t), pointer :: p1,p2
209 
210  p1 => this%p(edge_nodes(1, side))
211  p2 => this%p(edge_nodes(2, side))
212 
213  select type(t)
214  type is(tuple_i4_t)
215  if (p1%id() .lt. p2%id()) then
216  t%x = (/ p1%id(), p2%id() /)
217  else
218  t%x = (/ p2%id(), p1%id() /)
219  endif
220 
221  end select
222 
223  end subroutine hex_edge_id
224 
226  function hex_diameter(this) result(res)
227  class(hex_t), intent(in) :: this
228  real(kind=dp) :: d1, d2, d3, d4, res
229  type(point_t), pointer :: p1, p2, p3, p4, p5, p6, p7, p8
230  integer :: i
231 
232  d1 = 0d0
233  d2 = 0d0
234  d3 = 0d0
235  d4 = 0d0
236 
237  p1 => this%p(1)
238  p2 => this%p(2)
239  p3 => this%p(3)
240  p4 => this%p(4)
241  p5 => this%p(5)
242  p6 => this%p(6)
243  p7 => this%p(7)
244  p8 => this%p(8)
245 
246  do i = 1, neko_hex_gdim
247  d1 = d1 + (p8%x(i) - p1%x(i))**2
248  d2 = d2 + (p7%x(i) - p2%x(i))**2
249  d3 = d3 + (p5%x(i) - p4%x(i))**2
250  d4 = d4 + (p6%x(i) - p3%x(i))**2
251  end do
252 
253  res = sqrt(max(max(d1, d2), max(d3, d4)))
254 
255  end function hex_diameter
256 
258  function hex_centroid(this) result(res)
259  class(hex_t), intent(in) :: this
260  type(point_t) :: res
261  type(point_t), pointer :: p1, p2, p3, p4, p5, p6, p7, p8
262  integer :: i
263 
264  p1 => this%p(1)
265  p2 => this%p(2)
266  p3 => this%p(3)
267  p4 => this%p(4)
268  p5 => this%p(5)
269  p6 => this%p(6)
270  p7 => this%p(7)
271  p8 => this%p(8)
272  res%x = 0d0
273 
274  do i = 1, this%gdim()
275  res%x(i) = 0.125 * (p1%x(i) + p2%x(i) + p3%x(i) + p4%x(i) + &
276  p5%x(i) + p6%x(i) + p7%x(i) + p8%x(i))
277  end do
278 
279  end function hex_centroid
280 
283  pure function hex_equal(this, other) result(res)
284  class(hex_t), intent(in) :: this
285  class(element_t), intent(in) :: other
286  integer :: i
287  logical :: res
288 
289  res = .false.
290  select type(other)
291  type is (hex_t)
292  if ((this%gdim() .eq. other%gdim()) .and. &
293  (this%npts() .eq. other%npts())) then
294  do i = 1, this%npts()
295  if (this%pts(i)%p .ne. other%pts(i)%p) then
296  return
297  end if
298  end do
299  res = .true.
300  end if
301  end select
302 
303  end function hex_equal
304 
305 end module hex
Defines a hexahedron element.
Definition: hex.f90:34
type(point_t) function hex_centroid(this)
Compute the centroid of a hexahedron element.
Definition: hex.f90:259
subroutine hex_facet_order(this, t, side)
Return the ordered points for face i as a 4-tuple t.
Definition: hex.f90:185
subroutine hex_edge_id(this, t, side)
Return the edge id for an edge i as a 2-tuple t.
Definition: hex.f90:205
integer, dimension(4, 6), parameter face_nodes
Face node ids.
Definition: hex.f90:93
integer, parameter, public neko_hex_gdim
Geometric dimension.
Definition: hex.f90:45
integer, parameter, public neko_hex_npts
Number of points.
Definition: hex.f90:42
integer, parameter, public neko_hex_nfcs
Number of faces.
Definition: hex.f90:43
integer, parameter, public neko_hex_neds
Number of edges.
Definition: hex.f90:44
real(kind=dp) function hex_diameter(this)
Compute the diameter of a hexahedron element.
Definition: hex.f90:227
pure logical function hex_equal(this, other)
Check if two hex elements are equal.
Definition: hex.f90:284
integer, dimension(2, 12), parameter edge_nodes
Edge node ids.
Definition: hex.f90:119
subroutine hex_init(this, id, p1, p2, p3, p4, p5, p6, p7, p8)
Create a hexahedron element based upon eight points.
Definition: hex.f90:137
subroutine hex_facet_id(this, t, side)
Return the facet id for face i as a 4-tuple t.
Definition: hex.f90:156
integer, parameter, public dp
Definition: num_types.f90:9
Implements a point.
Definition: point.f90:35
Implements a n-tuple.
Definition: tuple.f90:34
Base type for an element.
Definition: element.f90:44
Hexahedron element.
Definition: hex.f90:63
A point in with coordinates .
Definition: point.f90:43
Integer based 4-tuple.
Definition: tuple.f90:69
Integer based 2-tuple.
Definition: tuple.f90:51
Base type for an n-tuple.
Definition: tuple.f90:41
#define max(a, b)
Definition: tensor.cu:40