-
Notifications
You must be signed in to change notification settings - Fork 0
/
SwanVertlist.ftn90
196 lines (196 loc) · 6.33 KB
/
SwanVertlist.ftn90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
subroutine SwanVertlist
!
! --|-----------------------------------------------------------|--
! | Delft University of Technology |
! | Faculty of Civil Engineering and Geosciences |
! | Environmental Fluid Mechanics Section |
! | P.O. Box 5048, 2600 GA Delft, The Netherlands |
! | |
! | Programmer: Marcel Zijlema |
! --|-----------------------------------------------------------|--
!
!
! SWAN (Simulating WAves Nearshore); a third generation wave model
! Copyright (C) 1993-2016 Delft University of Technology
!
! This program is free software; you can redistribute it and/or
! modify it under the terms of the GNU General Public License as
! published by the Free Software Foundation; either version 2 of
! the License, or (at your option) any later version.
!
! This program is distributed in the hope that it will be useful,
! but WITHOUT ANY WARRANTY; without even the implied warranty of
! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
! GNU General Public License for more details.
!
! A copy of the GNU General Public License is available at
! http://www.gnu.org/copyleft/gpl.html#SEC3
! or by writing to the Free Software Foundation, Inc.,
! 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
!
!
! Authors
!
! 40.80: Marcel Zijlema
! 41.07: Casey Dietrich
! 41.48: Marcel Zijlema
!
! Updates
!
! 40.80, July 2007: New subroutine
! 41.07, July 2009: small fix (assign ref.point to deepest point in case of no b.c.)
! 41.48, March 2013: including order along a user-given direction
!
! Purpose
!
! Makes vertex list with specific order
!
! Method
!
! Sorting based on distance with respect to a reference point
! This reference point can be either a vertex with boundary condition or a deepest point
!
! Modules used
!
use ocpcomm4
use m_genarr, only: DEPTH
use SwanGriddata
use SwanGridobjects
use SwanCompdata
!
implicit none
!
! Local variables
!
integer, save :: ient = 0 ! number of entries in this subroutine
integer :: istat ! indicate status of allocation
integer :: itmp ! temporary stored integer for swapping
integer :: j ! loop counter over vertices
integer :: k ! counter
integer, dimension(1) :: kd ! location of minimum value in array dist
integer, dimension(1) :: kx ! location of minimum value in array of x-coordinates of boundary vertices
integer, dimension(1) :: ky ! location of minimum value in array of y-coordinates of boundary vertices
!
real :: d1 ! distance of a point to origin
real :: d2 ! distance of another point to origin
real :: depmax ! maximum depth found
real :: rtmp ! temporary stored real for swapping
real :: x0 ! x-coordinate of reference point
real :: y0 ! y-coordinate of reference point
!
real, dimension(:), allocatable :: dist ! distance of each point with respect to reference point
!
type(verttype), dimension(:), pointer :: vert ! datastructure for vertices with their attributes
!
! Structure
!
! Description of the pseudo code
!
! Source text
!
if (ltrace) call strace (ient,'SwanVertlist')
!
! point to vertex object
!
vert => gridobject%vert_grid
!
istat = 0
if(.not.allocated(vlist)) allocate (vlist(nverts), stat = istat)
if ( istat /= 0 ) then
call msgerr ( 4, 'Allocation problem in SwanVertlist: array vlist ' )
return
endif
!
allocate (dist(nverts))
!
if ( asort /= -999. ) then
!
! order direction asort
!
do j = 1, nverts
dist(j) = vert(j)%attr(VERTX) * cos(asort) + vert(j)%attr(VERTY) * sin(asort)
enddo
!
else
!
! determine reference point
!
if ( all(mask=vert(:)%atti(VBC)==0) ) then
!
! if no boundary condition is given then find the vertex with the maximum depth
!
depmax = -999.
!
do j = 1, nverts
!
if ( DEPTH(j) > depmax ) then
!
depmax = DEPTH(j)
kx(1) = j
ky(1) = j
!
endif
!
enddo
!
else
!
! reference point is one of the vertices nearest to the origin where boundary condition is given
!
kx = minloc(vert(:)%attr(VERTX), vert(:)%atti(VBC)/=0)
ky = minloc(vert(:)%attr(VERTY), vert(:)%atti(VBC)/=0)
!
endif
!
if ( kx(1) == ky(1) ) then
x0 = vert(kx(1))%attr(VERTX)
y0 = vert(ky(1))%attr(VERTY)
else
!
d1 = sqrt((vert(kx(1))%attr(VERTX))**2+(vert(kx(1))%attr(VERTY))**2)
d2 = sqrt((vert(ky(1))%attr(VERTX))**2+(vert(ky(1))%attr(VERTY))**2)
!
if ( d1 < d2 ) then
x0 = vert(kx(1))%attr(VERTX)
y0 = vert(kx(1))%attr(VERTY)
else
x0 = vert(ky(1))%attr(VERTX)
y0 = vert(ky(1))%attr(VERTY)
endif
!
endif
!
! calculate distance of each point with respect to reference point
!
do j = 1, nverts
dist(j) = sqrt((vert(j)%attr(VERTX)-x0)**2 + (vert(j)%attr(VERTY)-y0)**2)
enddo
!
endif
!
! sort vertex list in order of increasing distance
!
vlist=(/ (j, j=1, nverts) /)
!
do j = 1, nverts-1
!
kd = minloc(dist(j:nverts))
k = kd(1) + j-1
!
if ( k /= j ) then
!
rtmp = dist(j)
dist(j) = dist(k)
dist(k) = rtmp
!
itmp = vlist(j)
vlist(j) = vlist(k)
vlist(k) = itmp
!
endif
!
enddo
!
deallocate(dist)
!
end subroutine SwanVertlist