forked from dnowacki-usgs/swanmod
-
Notifications
You must be signed in to change notification settings - Fork 0
/
SwanCreateEdges.f90
196 lines (196 loc) · 6.1 KB
/
SwanCreateEdges.f90
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 SwanCreateEdges
!
! --|-----------------------------------------------------------|--
! | 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-2015 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
!
! Updates
!
! 40.80, July 2007: New subroutine
!
! Purpose
!
! Generates edge-based data structure
!
! Method
!
! Faces of triangles are computed from elements and stored in Swan data structure
!
! Modules used
!
use ocpcomm4
use SwanGriddata
!
implicit none
!
! Local variables
!
integer, save :: ient = 0 ! number of entries in this subroutine
integer :: iface ! actual face of present cell
integer :: istat ! indicate status of allocation
integer :: j ! loop counter
integer :: k ! loop counter
integer :: m ! counter
integer :: maxnf ! maximum number of faces
integer :: n ! counter
integer :: v1 ! first vertex of present cell
integer :: v2 ! second vertex of present cell
integer :: v3 ! third vertex of present cell
integer, dimension(2,3) :: vf ! vertices of faces of present cell
!
integer, dimension(: ), allocatable :: cntv1 ! array of vertex counter for for vertex 1
integer, dimension(: ), allocatable :: cntv2 ! array of vertex counter for for vertex 2
integer, dimension(:,:), allocatable :: iflist1 ! list of index faces stored for vertex 1
integer, dimension(:,:), allocatable :: iflist2 ! list of index faces stored for vertex 2
!
logical :: facefound ! true if face is found
!
! Structure
!
! Description of the pseudo code
!
! Source text
!
if (ltrace) call strace (ient,'SwanCreateEdges')
!
! determine and store vertices of faces
!
maxnf = 3*ncells ! maximum number of faces
!
istat = 0
if(.not.allocated(kvertf)) allocate (kvertf(2,maxnf), stat = istat)
if ( istat /= 0 ) then
call msgerr ( 4, 'Allocation problem in SwanCreateEdges: array kvertf ' )
return
endif
kvertf = 0
!
allocate(cntv1 (nverts ))
allocate(cntv2 (nverts ))
allocate(iflist1(nverts,10))
allocate(iflist2(nverts,10))
!
cntv1 = 0
cntv2 = 0
iflist1 = -1
iflist2 = -2
!
nfaces = 0
!
do j = 1, ncells
!
v1 = kvertc(1,j)
v2 = kvertc(2,j)
v3 = kvertc(3,j)
!
vf(1,1) = v2
vf(2,1) = v3
vf(1,2) = v3
vf(2,2) = v1
vf(1,3) = v1
vf(2,3) = v2
!
kloop: do k = 1, 3
!
! get two vertices of a face
!
v1 = vf(1,k)
v2 = vf(2,k)
!
! search for identification number of that face
!
facefound = .false.
!
mloop: do m = 1, 10
!
iface = iflist1(v1,m)
!
do n = 1, 10
if ( iflist2(v2,n) == iface ) then
facefound = .true.
exit mloop
endif
enddo
!
enddo mloop
!
if ( .not.facefound ) then
!
nloop: do n = 1, 10
!
iface = iflist2(v1,n)
!
do m = 1, 10
if ( iflist1(v2,m) == iface ) then
facefound = .true.
exit nloop
endif
enddo
!
enddo nloop
!
endif
!
if ( facefound ) cycle kloop
!
nfaces = nfaces + 1
if ( nfaces > maxnf ) then
call msgerr ( 4, 'inconsistency found in SwanCreateEdges: more than maximum allowable faces found ' )
return
endif
!
m = cntv1(v1) +1
if ( m > 10 ) then
call msgerr ( 4, 'SwanCreateEdges: more than 10 faces around vertex ' )
return
endif
cntv1 (v1 ) = m
iflist1(v1,m) = nfaces
!
m = cntv2(v2) +1
if ( m > 10 ) then
call msgerr ( 4, 'SwanCreateEdges: more than 10 faces around vertex ' )
return
endif
cntv2 (v2 ) = m
iflist2(v2,m) = nfaces
!
kvertf(1,nfaces) = v1
kvertf(2,nfaces) = v2
!
enddo kloop
!
enddo
!
deallocate(cntv1,cntv2,iflist1,iflist2)
!
end subroutine SwanCreateEdges