-
Notifications
You must be signed in to change notification settings - Fork 0
/
SwanReadADCGrid.ftn90
211 lines (211 loc) · 6.96 KB
/
SwanReadADCGrid.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
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
subroutine SwanReadADCGrid
!
! --|-----------------------------------------------------------|--
! | 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
! 40.95: Marcel Zijlema
! 41.07: Casey Dietrich
!
! Updates
!
! 40.80, December 2007: New subroutine
! 40.95, June 2008: parallelization of unSWAN using MESSENGER of ADCIRC
! 41.07, August 2009: use ADCIRC boundary info to mark all boundary vertices
!
! Purpose
!
! Reads ADCIRC grid described in fort.14
!
! Method
!
! Grid coordinates of vertices are read from file fort.14 and stored in Swan data structure
! Vertices of triangles are read from file fort.14 and stored in Swan data structure
!
! Bottom topography from file fort.14 will also be stored
!
! Modules used
!
!PUN use ocpcomm2
use ocpcomm4
use m_genarr
use SwanGriddata
!PUN use SIZES
!PUN use MESSENGER
!
implicit none
!
! Local variables
!
character(80) :: grdfil ! name of grid file including path
integer, save :: ient = 0 ! number of entries in this subroutine
integer :: idum ! dummy integer
integer :: ii ! auxiliary integer
integer :: iostat ! I/O status in call FOR
integer :: istat ! indicate status of allocation
integer :: itype ! ADCIRC boundary type
integer :: ivert ! vertex index
integer :: ivert1 ! another vertex index
integer :: j ! loop counter
integer :: k ! loop counter
integer :: n1 ! auxiliary integer
integer :: n2 ! another auxiliary integer
integer :: ndsd ! unit reference number of file
integer :: nopbc ! number of open boundaries in ADCIRC
integer :: vm ! boundary marker
character(80) :: line ! auxiliary textline
logical :: stpnow ! indicate whether program must be terminated or not
!
! Structure
!
! Description of the pseudo code
!
! Source text
!
if (ltrace) call strace (ient,'SwanReadADCGrid')
!
! open file fort.14
!
ndsd = 0
iostat = 0
grdfil = 'fort.14'
!PUN grdfil = trim(INPUTDIR)//DIRCH2//trim(grdfil)
call for (ndsd, grdfil, 'OF', iostat)
if (stpnow()) goto 900
!
! skip first line
!
read(ndsd,'(a80)', end=950, err=910) line
!
! read number of elements and number of vertices
!
read(ndsd, *, end=950, err=910) ncells, nverts
istat = 0
if(.not.allocated(xcugrd)) allocate (xcugrd(nverts), stat = istat)
if ( istat == 0 ) then
if(.not.allocated(ycugrd)) allocate (ycugrd(nverts), stat = istat)
endif
if ( istat == 0 ) then
if(.not.allocated(DEPTH)) allocate (DEPTH(nverts), stat = istat)
endif
if ( istat /= 0 ) then
call msgerr ( 4, 'Allocation problem in SwanReadADCGrid: array xcugrd, ycugrd or depth ' )
goto 900
endif
!
! read coordinates of vertices and bottom topography
!
do j = 1, nverts
read(ndsd, *, end=950, err=910) ii, xcugrd(ii), ycugrd(ii), DEPTH(ii)
if ( ii/=j ) call msgerr ( 1, 'numbering of vertices is not sequential in grid file fort.14 ' )
enddo
!
if(.not.allocated(kvertc)) allocate (kvertc(3,ncells), stat = istat)
if ( istat /= 0 ) then
call msgerr ( 4, 'Allocation problem in SwanReadADCGrid: array kvertc ' )
goto 900
endif
!
! read vertices of triangles
!
do j = 1, ncells
read(ndsd, *, end=950, err=910) ii, idum, kvertc(1,ii), kvertc(2,ii), kvertc(3,ii)
if ( ii/=j ) call msgerr ( 1, 'numbering of triangles is not sequential in grid file fort.14 ' )
enddo
!
if(.not.allocated(vmark)) allocate (vmark(nverts), stat = istat)
if ( istat /= 0 ) then
call msgerr ( 4, 'Allocation problem in SwanReadADCGrid: array vmark ' )
goto 900
endif
vmark = 0
!
! read ADCIRC boundary information and store boundary markers
!
read(ndsd, *, end=950, err=910) nopbc
read(ndsd, *, end=950, err=910) idum
do j = 1, nopbc
!PUN if ( MNPROC==1 ) then
vm = j
read(ndsd, *, end=950, err=910) n2
!PUN else
!PUN read(ndsd, *, end=950, err=910) n2, vm
!PUN endif
do k = 1, n2
read(ndsd, *, end=950, err=910) ivert
vmark(ivert) = vm
enddo
enddo
!
read(ndsd, *, end=950, err=910) n1
read(ndsd, *, end=950, err=910) idum
do j = 1, n1
!PUN if ( MNPROC==1 ) then
vm = nopbc + j
read(ndsd, *, end=950, err=910) n2, itype
!PUN else
!PUN read(ndsd, *, end=950, err=910) n2, itype, vm
!PUN endif
if ( itype /= 4 .and. itype /= 24 ) then
do k = 1, n2
read(ndsd, *, end=950, err=910) ivert
vmark(ivert) = vm
enddo
else
do k = 1, n2
read(ndsd, *, end=950, err=910) ivert, ivert1
vmark(ivert ) = vm
vmark(ivert1) = vm
enddo
endif
enddo
!
! close file fort.14
!
close(ndsd)
!
!PUN ! ghost vertices are marked with +999
!PUN !
!PUN do j = 1, NEIGHPROC
!PUN do k = 1, NNODRECV(j)
!PUN ivert = IRECVLOC(k,j)
!PUN vmark(ivert) = 999
!PUN enddo
!PUN enddo
!PUN !
900 return
!
910 call msgerr (4, 'error reading data from grid file fort.14' )
goto 900
950 call msgerr (4, 'unexpected end of file in grid file fort.14' )
goto 900
!
end subroutine SwanReadADCGrid