forked from dnowacki-usgs/swanmod
-
Notifications
You must be signed in to change notification settings - Fork 0
/
SwanGridFace.f90
320 lines (319 loc) · 10.7 KB
/
SwanGridFace.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
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
subroutine SwanGridFace ( nfaces, ncells, nverts, xcugrd, ycugrd, kvertf )
!
! --|-----------------------------------------------------------|--
! | 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
!
! Fills face-based data structure
!
! Method
!
! Based on unstructured grid
!
! Modules used
!
use ocpcomm4
use SwanGridobjects
!ADC USE SIZES, ONLY: MYPROC
!
implicit none
!
! Argument variables
!
integer, intent(in) :: ncells ! number of cells in grid
integer, intent(in) :: nfaces ! number of faces in grid
integer, intent(in) :: nverts ! number of vertices in grid
!
integer, dimension(2, nfaces), intent(in) :: kvertf ! vertices of the face
! (must be filled by a gridgenerator!)
!
real, dimension(nverts), intent(in) :: xcugrd ! the x-coordinates of the grid vertices
real, dimension(nverts), intent(in) :: ycugrd ! the y-coordinates of the grid vertices
!
! Local variables
!
integer :: icell ! loop counter over cells
integer :: icell1 ! sequence number of cell 1 adjacent to present face
integer :: icell2 ! sequence number of cell 2 adjacent to present face
integer :: iface ! loop counter over faces
integer, save :: ient = 0 ! number of entries in this subroutine
integer :: ivert ! loop counter over vertices
integer :: j ! loop counter
integer :: jf ! loop counter
integer :: k ! counter
integer :: v1 ! first vertex of present face
integer :: vl1 ! first vertex of local face for given cell
integer :: v2 ! second vertex of present face
integer :: vl2 ! second vertex of local face for given 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
!
real :: carea1 ! area of cell 1 adjacent to present face
real :: carea2 ! area of cell 2 adjacent to present face
real :: lengthf ! length of present face
real :: xdiff ! difference in x-coordinate between vertex 2 and vertex 1
real :: ydiff ! difference in y-coordinate between vertex 2 and vertex 1
!
logical :: facefound ! true if face is found
!
type(celltype), dimension(:), pointer :: cell ! datastructure for cells with their attributes
type(facetype), dimension(:), pointer :: face ! datastructure for faces with their attributes
type(verttype), dimension(:), pointer :: vert ! datastructure for vertices with their attributes
!
! Structure
!
! Description of the pseudo code
!
! Source text
!
if (ltrace) call strace (ient,'SwanGridFace')
!
! point to vertex, cell and face objects
!
vert => gridobject%vert_grid
cell => gridobject%cell_grid
face => gridobject%face_grid
!
! loop over all faces
!
do iface = 1, nfaces
!
! identification number
!
face(iface)%atti(FACEID) = iface
!
! Fill vertices
!
v1 = kvertf(1,iface)
v2 = kvertf(2,iface)
face(iface)%atti(FACEV1) = v1
face(iface)%atti(FACEV2) = v2
!
! compute length of face
!
xdiff = xcugrd(v2) - xcugrd(v1)
ydiff = ycugrd(v2) - ycugrd(v1)
lengthf = sqrt( xdiff*xdiff + ydiff*ydiff )
!
! fill length of face and normal to face
!
face(iface)%attr(FACELEN ) = lengthf
face(iface)%attr(FACENORMX) = ydiff/lengthf
face(iface)%attr(FACENORMY) = -xdiff/lengthf
!
! compute coordinates of midpoint of face
!
face(iface)%attr(FACEMX) = 0.5*(xcugrd(v1) + xcugrd(v2))
face(iface)%attr(FACEMY) = 0.5*(ycugrd(v1) + ycugrd(v2))
!
face(iface)%atti(FACEC1 ) = 0
face(iface)%atti(FACEC2 ) = 0
face(iface)%atti(FMARKER) = 0
!
enddo
!
allocate(cntv1 (nverts ))
allocate(cntv2 (nverts ))
allocate(iflist1(nverts,10))
allocate(iflist2(nverts,10))
!
cntv1 = 0
cntv2 = 0
iflist1 = -1
iflist2 = -2
!
do iface = 1, nfaces
!
v1 = face(iface)%atti(FACEV1)
v2 = face(iface)%atti(FACEV2)
!
k = cntv1(v1) +1
if ( k > 10 ) then
!ADC PRINT *, "SWAN does not like local vertex ",v1," on core ",MYPROC
call msgerr ( 4, 'SwanGridFace: more than 10 faces around vertex ' )
return
endif
cntv1 (v1 ) = k
iflist1(v1,k) = iface
!
k = cntv2(v2) +1
if ( k > 10 ) then
!ADC PRINT *, "SWAN does not like local vertex ",v2," on core ",MYPROC
call msgerr ( 4, 'SwanGridFace: more than 10 faces around vertex ' )
return
endif
cntv2 (v2 ) = k
iflist2(v2,k) = iface
!
enddo
!
! loop over all cells
!
do icell = 1, ncells
!
! loop over all local faces of the cell
!
do jf = 1, cell(icell)%nof
!
! determine vertices of the local face
!
vl1 = cell(icell)%face(jf)%atti(FACEV1)
vl2 = cell(icell)%face(jf)%atti(FACEV2)
!
! search for identification number of that face
!
facefound = .false.
!
kloop: do k = 1, 10
!
iface = iflist1(vl1,k)
!
do j = 1, 10
if ( iflist2(vl2,j) == iface ) then
facefound = .true.
exit kloop
endif
enddo
!
enddo kloop
!
if ( .not.facefound ) then
!
jloop: do j = 1, 10
!
iface = iflist2(vl1,j)
!
do k = 1, 10
if ( iflist1(vl2,k) == iface ) then
facefound = .true.
exit jloop
endif
enddo
!
enddo jloop
!
endif
!
if ( facefound ) then
cell(icell)%face(jf)%atti(FACEID) = iface
else
call msgerr ( 4, 'inconsistency found in SwanGridFace: no face found ' )
return
endif
!
v1 = face(iface)%atti(FACEV1)
v2 = face(iface)%atti(FACEV2)
!
! Requirement: face(iface)%atti(FACEC1) < face(iface)%atti(FACEC2)
!
if ( v1 == vl1 ) then
if ( face(iface)%atti(FACEC1) == 0 ) then
face(iface)%atti(FACEC1) = icell
else
!ADC PRINT *, "SWAN does not like local element ",icell," on core ",MYPROC
call msgerr ( 4, 'SwanGridFace: not all cells have counterclockwise order of vertices ' )
return
endif
else
if ( face(iface)%atti(FACEC2) == 0 ) then
face(iface)%atti(FACEC2) = icell
else
!ADC PRINT *, "SWAN does not like local element ",icell," on core ",MYPROC
call msgerr ( 4, 'SwanGridFace: not all cells have counterclockwise order of vertices ' )
return
endif
endif
!
enddo
!
enddo
!
deallocate(cntv1,cntv2,iflist1,iflist2)
!
! loop over all faces
!
do iface = 1, nfaces
!
! marks boundary face
!
if (face(iface)%atti(FACEC2) == 0) face(iface)%atti(FMARKER) = 1
!
! store interpolation factors
!
if ( face(iface)%atti(FMARKER) == 1 ) then
!
face(iface)%attr(FACELINPF) = 0.
!
else
!
icell1 = face(iface)%atti(FACEC1)
icell2 = face(iface)%atti(FACEC2)
!
carea1 = cell(icell1)%attr(CELLAREA)
carea2 = cell(icell2)%attr(CELLAREA)
!
face(iface)%attr(FACELINPF) = carea1/(carea1+carea2)
!
endif
!
enddo
!
! marks boundary cells
!
do icell = 1, ncells
!
cell(icell)%atti(CMARKER) = 0
!
! loop over all local faces of the cell
!
do jf = 1, cell(icell)%nof
!
iface = cell(icell)%face(jf)%atti(FACEID)
!
if ( face(iface)%atti(FMARKER) == 1 ) then
cell(icell)%atti(CMARKER) = 1
exit
endif
!
enddo
!
enddo
end subroutine SwanGridFace