forked from dnowacki-usgs/swanmod
-
Notifications
You must be signed in to change notification settings - Fork 0
/
SwanGSECorr.f90
282 lines (282 loc) · 12.4 KB
/
SwanGSECorr.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
subroutine SwanGSECorr ( rhs, ac2, cgo, spcdir, idcmin, idcmax, isslow, isstop, trac0 )
!
! --|-----------------------------------------------------------|--
! | 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
!
! 41.00: Marcel Zijlema
!
! Updates
!
! 41.00, February 2009: New subroutine
!
! Purpose
!
! Computes waveage-dependent diffusion terms in x-y space to counteract the garden-sprinkler effect
!
! Modules used
!
use ocpcomm4
use swcomm2
use swcomm3
use swcomm4
use SwanGriddata
use SwanGridobjects
use SwanCompdata
!
implicit none
!
! Argument variables
!
integer, intent(in) :: isslow ! minimum frequency that is propagated within a sweep
integer, intent(in) :: isstop ! maximum frequency that is propagated within a sweep
!
integer, dimension(MSC), intent(in) :: idcmax ! maximum frequency-dependent counter in directional space
integer, dimension(MSC), intent(in) :: idcmin ! minimum frequency-dependent counter in directional space
!
real, dimension(MDC,MSC,nverts), intent(in) :: ac2 ! action density at current time level
real, dimension(MSC,ICMAX), intent(in) :: cgo ! group velocity
real, dimension(MDC,MSC), intent(inout) :: rhs ! right-hand side of system of equations in (sigma,theta) space
real, dimension(MDC,6), intent(in) :: spcdir ! (*,1): spectral direction bins (radians)
! (*,2): cosine of spectral directions
! (*,3): sine of spectral directions
! (*,4): cosine^2 of spectral directions
! (*,5): cosine*sine of spectral directions
! (*,6): sine^2 of spectral directions
real, dimension(MDC,MSC,MTRNP), intent(out) :: trac0 ! explicit part of propagation in present vertex for output purposes
!
! Local variables
!
integer :: icell ! index of present cell
integer :: id ! loop counter over direction bins
integer :: iddum ! counter in directional space for considered sweep
integer, save :: ient = 0 ! number of entries in this subroutine
integer :: is ! loop counter over frequency bins
integer :: ivert ! index of present vertex
integer :: jc ! loop counter
integer :: jcell ! index of next cell
!
integer, dimension(3) :: v ! vertices in present cell
!
real :: area0 ! area of present cell
real :: area1 ! area of next cell
double precision :: carea ! twices the area of centroid dual around present vertex
real :: cslat ! cosine of latitude
real :: dac2dx ! x-derivative of action density
real :: dac2dy ! y-derivative of action density
real :: dcg ! group velocity difference across frequency bin
real :: dgx0 ! x-component of diffusion gradient inside present cell
real :: dgx1 ! x-component of diffusion gradient inside next cell
real :: dgxdx ! x-gradient of x-diffusion gradient component
real :: dgy0 ! y-component of diffusion gradient inside present cell
real :: dgy1 ! y-component of diffusion gradient inside next cell
real :: dgydy ! y-gradient of y-diffusion gradient component
real :: dnn ! waveage-dependent diffusion coefficient normal to propagation direction
real :: dss ! waveage-dependent diffusion coefficient in propagation direction
real :: dxx ! xx-component of diffusion coefficient in Cartesian coordinates
real :: dxy ! xy-component of diffusion coefficient in Cartesian coordinates
real :: dyy ! yy-component of diffusion coefficient in Cartesian coordinates
double precision :: x0 ! x-coordinate of the centroid of present cell
double precision :: x1 ! x-coordinate of the centroid of next cell
double precision :: y0 ! y-coordinate of the centroid of present cell
double precision :: y1 ! y-coordinate of the centroid of next cell
!
type(celltype), dimension(:), pointer :: cell ! datastructure for cells 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,'SwanGSECorr')
!
! point to vertex and cell objects
!
vert => gridobject%vert_grid
cell => gridobject%cell_grid
!
ivert = vs(1)
!
if ( vert(ivert)%atti(VMARKER) == 1 ) return ! no GSE correction in boundary vertex
!
cslat = cos(DEGRAD*(vert(ivert)%attr(VERTY) + YOFFS))
!
do is = isslow, isstop
!
! calculate waveage-dependent diffusion coefficients in polar coordinates
!
if ( is == 1 ) then
dcg = abs(cgo(is+1,1)-cgo(is,1))
elseif ( is == isstop ) then
dcg = abs(cgo(is,1)-cgo(is-1,1))
else
dcg = 0.5 * abs(cgo(is+1,1)-cgo(is-1,1))
endif
!
dss = dcg**2*WAVAGE/12.
dnn = (cgo(is,1)*DDIR)**2 * WAVAGE/12.
!
do iddum = idcmin(is), idcmax(is)
id = mod ( iddum - 1 + MDC , MDC ) + 1
!
! calculate diffusion coefficients in Cartesian coordinates
!
dxx = dss*spcdir(id,4) + dnn*spcdir(id,6)
dyy = dss*spcdir(id,6) + dnn*spcdir(id,4)
dxy = (dss-dnn)*spcdir(id,5)
!
! compute contribution to the diffusion terms in present vertex
!
carea = 0d0
dgxdx = 0.
dgydy = 0.
!
! loop over cells around considered vertex
!
do jc = 1, vert(ivert)%noc
!
! get present cell and its vertices
!
icell = vert(ivert)%cell(jc)%atti(CELLID)
!
v(1) = cell(icell)%atti(CELLV1)
v(2) = cell(icell)%atti(CELLV2)
v(3) = cell(icell)%atti(CELLV3)
!
! determine centroid and area of present cell
!
x0 = cell(icell)%attr(CELLCX)
y0 = cell(icell)%attr(CELLCY)
area0 = cell(icell)%attr(CELLAREA)
!
! determine derivatives of action density inside present cell
!
dac2dx = 0.5*( ac2(id,is,v(1))*(ycugrd(v(2))-ycugrd(v(3))) + &
ac2(id,is,v(2))*(ycugrd(v(3))-ycugrd(v(1))) + &
ac2(id,is,v(3))*(ycugrd(v(1))-ycugrd(v(2))) )/area0
!
dac2dy = 0.5*( ac2(id,is,v(1))*(xcugrd(v(3))-xcugrd(v(2))) + &
ac2(id,is,v(2))*(xcugrd(v(1))-xcugrd(v(3))) + &
ac2(id,is,v(3))*(xcugrd(v(2))-xcugrd(v(1))) )/area0
!
! in case of spherical coordinates, transform back to Cartesian coordinates
!
if ( KSPHER > 0 ) then
!
dac2dx = dac2dx/(cslat * LENDEG)
dac2dy = dac2dy/LENDEG
!
endif
!
! determine diffusion gradients in centroid of present cell
!
dgx0 = dxx*dac2dx + dxy*dac2dy
dgy0 = dxy*dac2dx + dyy*dac2dy
!
! get next cell in counterclockwise direction
!
jcell = vert(ivert)%cell(jc)%atti(NEXTCELL)
!
v(1) = cell(jcell)%atti(CELLV1)
v(2) = cell(jcell)%atti(CELLV2)
v(3) = cell(jcell)%atti(CELLV3)
!
! determine centroid and area of next cell
!
x1 = cell(jcell)%attr(CELLCX)
y1 = cell(jcell)%attr(CELLCY)
area1 = cell(jcell)%attr(CELLAREA)
!
! determine derivatives of action density inside next cell
!
dac2dx = 0.5*( ac2(id,is,v(1))*(ycugrd(v(2))-ycugrd(v(3))) + &
ac2(id,is,v(2))*(ycugrd(v(3))-ycugrd(v(1))) + &
ac2(id,is,v(3))*(ycugrd(v(1))-ycugrd(v(2))) )/area1
!
dac2dy = 0.5*( ac2(id,is,v(1))*(xcugrd(v(3))-xcugrd(v(2))) + &
ac2(id,is,v(2))*(xcugrd(v(1))-xcugrd(v(3))) + &
ac2(id,is,v(3))*(xcugrd(v(2))-xcugrd(v(1))) )/area1
!
! in case of spherical coordinates, transform back to Cartesian coordinates
!
if ( KSPHER > 0 ) then
!
dac2dx = dac2dx/(cslat * LENDEG)
dac2dy = dac2dy/LENDEG
!
endif
!
! determine diffusion gradients in centroid of next cell
!
dgx1 = dxx*dac2dx + dxy*dac2dy
dgy1 = dxy*dac2dx + dyy*dac2dy
!
! compute contribution to area of centroid dual
!
carea = carea + x0*y1 - x1*y0
!
! compute contribution to x-gradient of x-diffusion gradient in centroid dual
!
dgxdx = dgxdx + ( dgx0 + dgx1 ) * real( y1 - y0 )
!
! compute contribution to y-gradient of y-diffusion gradient in centroid dual
!
dgydy = dgydy + ( dgy0 + dgy1 ) * real( x0 - x1 )
!
enddo
!
! add diffusion terms to the right-hand side of the action balance equation
!
if ( carea > 0d0 ) then
!
dgxdx = dgxdx/real(carea)
dgydy = dgydy/real(carea)
!
! in case of spherical coordinates, transform back to Cartesian coordinates
!
if ( KSPHER > 0 ) then
!
dgxdx = dgxdx/(cslat * LENDEG)
dgydy = dgydy/LENDEG
!
endif
!
rhs (id,is ) = rhs (id,is ) + dgxdx + dgydy
trac0(id,is,1) = trac0(id,is,1) - dgxdx - dgydy
!
endif
!
enddo
!
enddo
!
end subroutine SwanGSECorr