-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathSwanReadGrid.ftn90
131 lines (131 loc) · 3.63 KB
/
SwanReadGrid.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
subroutine SwanReadGrid ( basenm, lenfnm )
!
! --|-----------------------------------------------------------|--
! | 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-2024 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 3 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.
!
! You should have received a copy of the GNU General Public License
! along with this program. If not, see <http://www.gnu.org/licenses/>.
!
!
! Authors
!
! 40.80: Marcel Zijlema
!
! Updates
!
! 40.80, July 2007: New subroutine
!
! Purpose
!
! Reads unstructured grid generated by a grid generator
!
! Method
!
! Reads data from either ADCIRC, Triangle or Easymesh
!
! Modules used
!
use ocpcomm4
use swcomm2
use swcomm3
use SwanGriddata
!
implicit none
!
! Argument variables
!
integer, intent(in) :: lenfnm ! length of file names
character(lenfnm), intent(in) :: basenm ! base name of unstructured grid files
!
! Local variables
!
integer :: i ! loop counter
integer, save :: ient = 0 ! number of entries in this subroutine
!
! Structure
!
! Description of the pseudo code
!
! Source text
!
if (ltrace) call strace (ient,'SwanReadGrid')
!
if ( grid_generator == meth_adcirc ) then
!
! grid is given in ADCIRC grid file (fort.14)
!
call SwanReadADCGrid
!
elseif ( grid_generator == meth_triangle ) then
!
! grid is generated by Triangle
!
call SwanReadTriangleGrid ( basenm, lenfnm )
!
elseif ( grid_generator == meth_easy ) then
!
! grid is generated by Easymesh
!
call SwanReadEasymeshGrid ( basenm, lenfnm )
!
else
!
call msgerr ( 4, 'Unknown grid generator' )
return
!
endif
!
! compute coordinate offsets and reset grid coordinates
!
do i = 1, nverts
if ( .not.LXOFFS ) then
XOFFS = xcugrd(i)
YOFFS = ycugrd(i)
LXOFFS = .true.
xcugrd(i) = 0.
ycugrd(i) = 0.
else
xcugrd(i) = real(xcugrd(i) - dble(XOFFS))
ycugrd(i) = real(ycugrd(i) - dble(YOFFS))
endif
enddo
!
! compute XCGMIN, XCGMAX, YCGMIN, YCGMAX
!
XCGMIN = 1.e9
YCGMIN = 1.e9
XCGMAX = -1.e9
YCGMAX = -1.e9
do i = 1, nverts
if ( xcugrd(i) < XCGMIN ) XCGMIN = xcugrd(i)
if ( ycugrd(i) < YCGMIN ) YCGMIN = ycugrd(i)
if ( xcugrd(i) > XCGMAX ) XCGMAX = xcugrd(i)
if ( ycugrd(i) > YCGMAX ) YCGMAX = ycugrd(i)
enddo
!
! compute lengths of enclosure of computational domain
!
XCLEN = XCGMAX - XCGMIN
YCLEN = YCGMAX - YCGMIN
!
end subroutine SwanReadGrid