Loading [MathJax]/extensions/tex2jax.js
Neko 0.9.99
A portable framework for high-order spectral element flow simulations
All Classes Namespaces Files Functions Variables Typedefs Enumerator Macros Pages
tri.f90
Go to the documentation of this file.
1! Copyright (c) 2021-2023, 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!
34module tri
35 use num_types, only : dp
36 use element, only : element_t
37 use tuple, only : tuple_t, tuple_i4_t
38 use point, only : point_t
39 implicit none
40 private
41
42 integer, public, parameter :: neko_tri_npts = 3
43 integer, public, parameter :: neko_tri_neds = 3
44 integer, public, parameter :: neko_tri_gdim = 2
45 integer, public, parameter :: neko_tri_sdim = 3
46
61 type, public, extends(element_t) :: tri_t
62 contains
63 procedure, pass(this) :: init => tri_init
64 procedure, pass(this) :: facet_id => tri_facet_id
65 procedure, pass(this) :: facet_order => tri_facet_order
66 procedure, pass(this) :: diameter => tri_diameter
67 procedure, pass(this) :: centroid => tri_centroid
68 procedure, pass(this) :: equal => tri_equal
69 generic :: operator(.eq.) => equal
70 end type tri_t
71
86 integer, parameter, dimension(2, 3) :: edge_nodes = reshape((/1,3,&
87 2,3,&
88 1,2 /),&
89 (/2,3/))
90
91contains
92
94 subroutine tri_init(this, id, p1, p2, p3)
95 class(tri_t), intent(inout) :: this
96 integer, intent(inout) :: id
97 type(point_t), target, intent(in) :: p1, p2, p3
98
99 call this%element(id, neko_tri_gdim, neko_tri_npts)
100
101 this%pts(1)%p => p1
102 this%pts(2)%p => p2
103 this%pts(3)%p => p3
104
105 end subroutine tri_init
106
109 subroutine tri_facet_id(this, t, side)
110 class(tri_t), intent(in) :: this
111 class(tuple_t), intent(inout) :: t
112 integer, intent(in) :: side
113 type(point_t), pointer :: p1, p2
114
115 p1 => this%p(edge_nodes(1, side))
116 p2 => this%p(edge_nodes(2, side))
117
118 select type (t)
119 type is (tuple_i4_t)
120 if (p1 .lt. p2) then
121 t%x = (/ p1%id(), p2%id() /)
122 else
123 t%x = (/ p2%id(), p1%id() /)
124 end if
125 end select
126
127 end subroutine tri_facet_id
128
130 subroutine tri_facet_order(this, t, side)
131 class(tri_t), intent(in) :: this
132 class(tuple_t), intent(inout) :: t
133 integer, intent(in) :: side
134 type(point_t), pointer :: p1, p2
135
136 p1 => this%p(edge_nodes(1, side))
137 p2 => this%p(edge_nodes(2, side))
138
139 select type (t)
140 type is (tuple_i4_t)
141 t%x = (/ p1%id(), p2%id() /)
142 end select
143
144 end subroutine tri_facet_order
145
147 function tri_diameter(this) result(res)
148 class(tri_t), intent(in) :: this
149 real(kind=dp) :: d1, d2, d3, res
150 type(point_t), pointer :: p1, p2, p3
151 integer :: i
152
153 d1 = 0d0
154 d2 = 0d0
155 d3 = 0d0
156
157 p1 => this%p(1)
158 p2 => this%p(2)
159 p3 => this%p(3)
160
161 do i = 1, neko_tri_sdim
162 d1 = d1 + (p2%x(i) - p1%x(i))**2
163 d2 = d2 + (p3%x(i) - p2%x(i))**2
164 d3 = d3 + (p1%x(i) - p3%x(i))**2
165 end do
166
167 res = sqrt(max(max(d1, d2), d3))
168
169 end function tri_diameter
170
172 function tri_centroid(this) result(res)
173 class(tri_t), intent(in) :: this
174 type(point_t) :: res
175 type(point_t), pointer :: p1, p2, p3
176 integer :: i
177
178 p1 => this%p(1)
179 p2 => this%p(2)
180 p3 => this%p(3)
181 res%x = 0d0
182
183 do i = 1, neko_tri_sdim
184 res%x(i) = 1d0/3d0 * (p1%x(i) + p2%x(i) + p3%x(i))
185 end do
186 end function tri_centroid
187
190 pure function tri_equal(this, other) result(res)
191 class(tri_t), intent(in) :: this
192 class(element_t), intent(in) :: other
193 integer :: i
194 logical :: res
195
196 res = .false.
197 select type (other)
198 class is (tri_t)
199 if ((this%gdim() .eq. other%gdim()) .and. &
200 (this%npts() .eq. other%npts())) then
201 do i = 1, neko_tri_npts
202 if (this%pts(i)%p .ne. other%pts(i)%p) then
203 return
204 end if
205 end do
206 res = .true.
207 end if
208 end select
209
210 end function tri_equal
211
212end module tri
integer, parameter, public dp
Definition num_types.f90:9
Implements a point.
Definition point.f90:35
Defines a triangular element.
Definition tri.f90:34
integer, parameter, public neko_tri_npts
Number of points.
Definition tri.f90:42
integer, parameter, public neko_tri_sdim
Dimension of space.
Definition tri.f90:45
integer, parameter, public neko_tri_gdim
Geometric dimension.
Definition tri.f90:44
real(kind=dp) function tri_diameter(this)
Compute the diameter of a triangular element.
Definition tri.f90:148
subroutine tri_facet_order(this, t, side)
Return the ordered edge for face i as a 2-tuple t.
Definition tri.f90:131
integer, parameter, public neko_tri_neds
Number of edges.
Definition tri.f90:43
subroutine tri_facet_id(this, t, side)
Return the edge id for face i as a 2-tuple t.
Definition tri.f90:110
type(point_t) function tri_centroid(this)
Compute the centroid of a triangular element.
Definition tri.f90:173
subroutine tri_init(this, id, p1, p2, p3)
Create a trinagular element based upon three points.
Definition tri.f90:95
integer, dimension(2, 3), parameter edge_nodes
Edge node ids.
Definition tri.f90:86
pure logical function tri_equal(this, other)
Check if two triangle elements are equal.
Definition tri.f90:191
Implements a n-tuple.
Definition tuple.f90:34
Base type for an element.
Definition element.f90:44
A point in with coordinates .
Definition point.f90:43
Triangular element.
Definition tri.f90:61
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