Neko 0.9.99
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
element.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!
33module element
34 use num_types
35 use entity, only : entity_t
36 use tuple, only : tuple_t
37 use point, only : point_ptr, point_t
38 implicit none
39 private
40
44 type, public, extends(entity_t), abstract :: element_t
45 integer, private :: gdim_
46 integer, private :: npts_
47 type(point_ptr), allocatable :: pts(:)
48 contains
49 procedure, pass(this) :: element => element_init
50 procedure, pass(this) :: free => element_free
51 procedure, pass(this) :: gdim => element_gdim
52 procedure, pass(this) :: npts => element_npts
53 procedure, pass(this) :: p => element_point
54 procedure, pass(this) :: n_points => element_npts
55 procedure, pass(this), non_overridable :: element_point
56 procedure(element_equal), pass(this), deferred :: equal
57 procedure(element_diameter), pass(this), deferred :: diameter
58 procedure(element_centroid), pass(this), deferred :: centroid
59 procedure(element_facet_id), pass(this), deferred :: facet_id
60 procedure(element_facet_order), pass(this), deferred :: facet_order
61 end type element_t
62
63 abstract interface
64 function element_diameter(this) result(res)
65 import :: element_t
66 import :: dp
67 class(element_t), intent(in) :: this
68 real(kind=dp) :: res
69 end function element_diameter
70 end interface
71
72 abstract interface
73 function element_centroid(this) result(res)
74 import :: element_t
75 import :: point_t
76 class(element_t), intent(in) :: this
77 type(point_t) :: res
78 end function element_centroid
79 end interface
80
81 abstract interface
82 pure function element_equal(this, other) result(res)
83 import :: element_t
84 class(element_t), intent(in) :: this
85 class(element_t), intent(in) :: other
86 logical :: res
87 end function element_equal
88 end interface
89
90 abstract interface
91 subroutine element_facet_id(this, t, side)
92 import :: element_t
93 import :: tuple_t
94 class(element_t), intent(in) :: this
95 class(tuple_t), intent(inout) :: t
96 integer, intent(in) :: side
97 end subroutine element_facet_id
98 end interface
99
100 abstract interface
101 subroutine element_facet_order(this, t, side)
102 import :: element_t
103 import :: tuple_t
104 class(element_t), intent(in) :: this
105 class(tuple_t), intent(inout) :: t
106 integer, intent(in) :: side
107 end subroutine element_facet_order
108 end interface
109contains
110
112 subroutine element_init(this, id, gdim, npts)
113 class(element_t), intent(inout) :: this
114 integer, intent(inout) :: id
115 integer, intent(in) :: gdim
116 integer, intent(in) :: npts
117
118 call this%free()
119
120 call this%set_id(id)
121
122 this%gdim_ = gdim
123 this%npts_ = npts
124
125 allocate(this%pts(this%npts_))
126
127 end subroutine element_init
128
130 subroutine element_free(this)
131 class(element_t), intent(inout) :: this
132
133 if (allocated(this%pts)) then
134 deallocate(this%pts)
135 end if
136
137 end subroutine element_free
138
140 pure function element_gdim(this) result(gdim)
141 class(element_t), intent(in) :: this
142 integer :: gdim
143 gdim = this%gdim_
144 end function element_gdim
145
147 pure function element_npts(this) result(npts)
148 class(element_t), intent(in) :: this
149 integer :: npts
150 npts = this%npts_
151 end function element_npts
152
154 function element_point(this, i) result(pt)
155 class(element_t), intent(in) :: this
156 integer, intent(in) :: i
157 type(point_t), pointer :: pt
158 pt => this%pts(i)%p
159 end function element_point
160
161end module element
subroutine element_free(this)
Deallocate an element.
Definition element.f90:131
subroutine element_init(this, id, gdim, npts)
Create an element with npts.
Definition element.f90:113
type(point_t) function, pointer element_point(this, i)
Return a pointer to point i of the element.
Definition element.f90:155
pure integer function element_npts(this)
Get the number of points in an element.
Definition element.f90:148
pure integer function element_gdim(this)
Get the geometric dimension of an element.
Definition element.f90:141
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
Base type for an entity.
Definition entity.f90:38
Defines a pointer to a point type.
Definition point.f90:67
A point in with coordinates .
Definition point.f90:43
Base type for an n-tuple.
Definition tuple.f90:41