Neko 1.99.1
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
point.f90
Go to the documentation of this file.
1! Copyright (c) 2019-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!
34!
35module point
36 use num_types, only : dp, rp
37 use math, only : abscmp
38 use entity, only : entity_t
39 implicit none
40 private
41
43 type, extends(entity_t), public :: point_t
44 real(kind=dp), dimension(3) :: x
45 contains
46 generic :: init => point_init_array, point_init_xyz
47 procedure, pass(this) :: point_init_array
48 procedure, pass(this) :: point_init_xyz
49 procedure :: point_eq
50 procedure :: point_ne
51 procedure :: point_lt
52 procedure :: point_gt
53 procedure :: point_assign
54 procedure :: point_add
55 procedure :: point_subtract
56 procedure :: point_scalar_mult
57 procedure, pass(p1) :: dist => point_euclid_dist
58 procedure, pass(x) :: point_mat_mult
59 generic :: operator(.eq.) => point_eq
60 generic :: operator(.ne.) => point_ne
61 generic :: operator(.lt.) => point_lt
62 generic :: operator(.gt.) => point_gt
63 generic :: assignment(=) => point_assign
64 generic :: operator(+) => point_add
65 generic :: operator(-) => point_subtract
66 generic :: operator(*) => point_scalar_mult, point_mat_mult
67 end type point_t
68
70 type, public :: point_ptr
71 type(point_t), pointer :: p
72 end type point_ptr
73
74contains
75
77 subroutine point_init_array(this, x, id)
78 class(point_t), intent(inout) :: this
79 real(kind=dp), dimension(3), intent(in) :: x
80 integer, optional, intent(in) :: id
81
82 if (present(id)) then
83 call this%set_id(id)
84 else
85 call this%set_id(-1)
86 end if
87
88 this%x = x
89
90 end subroutine point_init_array
91
93 subroutine point_init_xyz(this, x, y, z, id)
94 class(point_t), intent(inout) :: this
95 real(kind=dp), intent(in) :: x
96 real(kind=dp), intent(in) :: y
97 real(kind=dp), intent(in) :: z
98 integer, optional, intent(in) :: id
99
100 if (present(id)) then
101 call this%set_id(id)
102 else
103 call this%set_id(-1)
104 end if
105
106 this%x(1) = x
107 this%x(2) = y
108 this%x(3) = z
109
110 end subroutine point_init_xyz
111
113 subroutine point_assign(this, x)
114 class(point_t), intent(inout) :: this
115 real(kind=dp), dimension(3), intent(in) :: x
116
117 this%x = x
118
119 end subroutine point_assign
120
123 pure function point_eq(p1, p2) result(res)
124 class(point_t), intent(in) :: p1
125 class(point_t), intent(in) :: p2
126 logical :: res
127
128 if (abscmp(p1%x(1), p2%x(1)) .and. &
129 abscmp(p1%x(2), p2%x(2)) .and. &
130 abscmp(p1%x(3), p2%x(3))) then
131 res = .true.
132 else
133 res = .false.
134 end if
135
136 end function point_eq
137
140 pure function point_ne(p1, p2) result(res)
141 class(point_t), intent(in) :: p1
142 class(point_t), intent(in) :: p2
143 logical :: res
144
145 if (.not. abscmp(p1%x(1), p2%x(1)) .or. &
146 .not. abscmp(p1%x(2), p2%x(2)) .or. &
147 .not. abscmp(p1%x(3), p2%x(3))) then
148 res = .true.
149 else
150 res = .false.
151 end if
152
153 end function point_ne
154
157 pure function point_lt(p1, p2) result(res)
158 class(point_t), intent(in) :: p1
159 class(point_t), intent(in) :: p2
160 logical :: res
161
162 if (p1%x(1) .lt. p2%x(1) .or. &
163 (abscmp(p1%x(1), p2%x(1)) .and. &
164 (p1%x(2) .lt. p2%x(2) .or. &
165 (abscmp(p1%x(2), p2%x(2)) .and. p1%x(3) .lt. p2%x(3))))) then
166 res = .true.
167 else
168 res = .false.
169 end if
170
171 end function point_lt
172
175 pure function point_gt(p1, p2) result(res)
176 class(point_t), intent(in) :: p1
177 class(point_t), intent(in) :: p2
178 logical :: res
179
180 if (point_lt(p1, p2)) then
181 res = .false.
182 else
183 res = .true.
184 end if
185
186 end function point_gt
187
189 function point_add(p1, p2) result(res)
190 class(point_t), intent(in) :: p1
191 class(point_t), intent(in) :: p2
192 type(point_t) :: res
193
194 res%x(1) = p1%x(1) + p2%x(1)
195 res%x(2) = p1%x(2) + p2%x(2)
196 res%x(3) = p1%x(3) + p2%x(3)
197
198 end function point_add
199
201 function point_subtract(p1, p2) result(res)
202 class(point_t), intent(in) :: p1
203 class(point_t), intent(in) :: p2
204 type(point_t) :: res
205
206 res%x(1) = p1%x(1) - p2%x(1)
207 res%x(2) = p1%x(2) - p2%x(2)
208 res%x(3) = p1%x(3) - p2%x(3)
209
210 end function point_subtract
211
213 function point_scalar_mult(p, a) result(res)
214 class(point_t), intent(in) :: p
215 real(kind=rp), intent(in) :: a
216 type(point_t) :: res
217
218 res%x(1) = a * p%x(1)
219 res%x(2) = a * p%x(2)
220 res%x(3) = a * p%x(3)
221
222 end function point_scalar_mult
223
225 pure function point_euclid_dist(p1, p2) result(res)
226 class(point_t), intent(in) :: p1
227 type(point_t), intent(in) :: p2
228 real(kind=rp) :: res
229
230 res = sqrt( (p1%x(1) - p2%x(1))**2 &
231 + (p1%x(2) - p2%x(2))**2 &
232 + (p1%x(3) - p2%x(3))**2 )
233 end function point_euclid_dist
234
236 function point_mat_mult(A,x) result(b)
237 class(point_t), intent(in) :: x
238 real(kind=rp), intent(in) :: a(3,3)
239 type(point_t) :: b
240 integer :: i,j
241
242 b%x = 0.0_rp
243
244 do i = 1, 3
245 do j = 1, 3
246 b%x(i) = b%x(i) + a(i,j) * x%x(j)
247 end do
248 end do
249
250 end function point_mat_mult
251
252end module point
Definition math.f90:60
integer, parameter, public dp
Definition num_types.f90:9
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Implements a point.
Definition point.f90:35
subroutine point_init_array(this, x, id)
Initialize a point from an array x of coordinates.
Definition point.f90:78
pure logical function point_gt(p1, p2)
Check if .
Definition point.f90:176
pure logical function point_eq(p1, p2)
Check if .
Definition point.f90:124
subroutine point_init_xyz(this, x, y, z, id)
Initialize a point from coordinates.
Definition point.f90:94
type(point_t) function point_subtract(p1, p2)
Returns the subtraction of 2 points .
Definition point.f90:202
type(point_t) function point_mat_mult(a, x)
Computes matrix-vector product in : .
Definition point.f90:237
subroutine point_assign(this, x)
Assigns coordinates x to a point.
Definition point.f90:114
type(point_t) function point_scalar_mult(p, a)
Returns the multiplication of a point by a scalar .
Definition point.f90:214
pure logical function point_lt(p1, p2)
Check if .
Definition point.f90:158
pure real(kind=rp) function point_euclid_dist(p1, p2)
Returns the Euclidean distance between two points .
Definition point.f90:226
type(point_t) function point_add(p1, p2)
Returns the addition of 2 points .
Definition point.f90:190
pure logical function point_ne(p1, p2)
Check if .
Definition point.f90:141
Base type for an entity.
Definition entity.f90:38
Defines a pointer to a point type.
Definition point.f90:70
A point in with coordinates .
Definition point.f90:43