-
Notifications
You must be signed in to change notification settings - Fork 0
/
SwanGradDepthorK.ftn90
390 lines (390 loc) · 13.9 KB
/
SwanGradDepthorK.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
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
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
subroutine SwanGradDepthorK ( dep2, mudl2, spcsig, dhdx, dhdy, dkdx, dkdy, ivert )
!
! --|-----------------------------------------------------------|--
! | 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
!
! 41.60: Marcel Zijlema
!
! Updates
!
! 41.60, July 2015: New subroutine
!
! Purpose
!
! Computes gradients of depth or wave number in vertex
! meant for computing turning rate
!
! Method
!
! Application of the Green-Gauss theorem with the assumption of
! a constant gradient over the controle volume (centroid dual)
!
! Modules used
!
use ocpcomm4
use swcomm2
use swcomm3
use swcomm4
use SwanGriddata
use SwanGridobjects
!
implicit none
!
! Argument variables
!
integer , intent(in) :: ivert ! counter corresponding to current vertex
!
real, dimension(nverts), intent(in) :: dep2 ! water depth at current time level
real , intent(out) :: dhdx ! derivative of depth in x-direction
real , intent(out) :: dhdy ! derivative of depth in y-direction
real, dimension(MSC) , intent(out) :: dkdx ! derivative of wave number in x-direction
real, dimension(MSC) , intent(out) :: dkdy ! derivative of wave number in y-direction
real, dimension(nverts), intent(in) :: mudl2 ! mud layer at current time level
real, dimension(MSC) , intent(in) :: spcsig ! relative frequency bins
!
! Local variables
!
integer :: icell ! index of present cell
integer, save :: ient = 0 ! number of entries in this subroutine
integer :: jc ! loop counter
integer :: jcell ! index of next cell
!
integer, dimension(3) :: v ! vertices in present cell
!
double precision :: area ! twices the area of centroid dual around present vertex
real, dimension(MSC) :: arr ! auxiliary array
real :: cslat ! cosine of latitude
real :: dpmax ! maximum depth found in centroid dual
real :: dpmin ! minimum depth found in centroid dual
real, dimension(3) :: dloc ! local depth at vertices
real, dimension(3) :: dm ! local mud layer at vertices
real :: dmaxc ! maximum depth found per cell of centroid dual
real :: dminc ! minimum depth found per cell of centroid dual
real, parameter :: drat= 5. ! ratio between maximum and minimum depths in centroid dual
real :: h0 ! depth in centroid of present cell
real :: h1 ! depth in centroid of next cell
real, dimension(MSC) :: k0 ! wave number in centroid of present cell
real, dimension(MSC) :: k1 ! wave number in centroid of next cell
real, dimension(MSC,3) :: kloc ! local wave number at vertices
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
!
character(80) :: msgstr ! string to pass message
!
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,'SwanGradDepthorK')
!
dhdx = 0.
dhdy = 0.
!
dkdx = 0.
dkdy = 0.
!
! if no frequency shift and no refraction, return
!
if ( (ITFRE == 0 .or. ICUR == 0) .and. IREFR == 0 ) return
!
! point to vertex and cell objects
!
vert => gridobject%vert_grid
cell => gridobject%cell_grid
!
if ( vert(ivert)%atti(VMARKER) == 1 ) return ! boundary vertex
!
if ( (IREFR /= 0 .and. int(PNUMS(32)) == 0) .or. (ITFRE /= 0 .and. ICUR /= 0) ) then
!
area = 0d0
dpmax = -99999.
dpmin = 99999.
!
! 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)
!
dloc(1) = dep2(v(1))
dloc(2) = dep2(v(2))
dloc(3) = dep2(v(3))
!
if ( dloc(1) <= DEPMIN .or. dloc(2) <= DEPMIN .or. dloc(3) <= DEPMIN ) goto 10
!
dminc = min ( min( dloc(1), dloc(2) ), dloc(3) )
if ( dminc < dpmin ) dpmin = dminc
!
dmaxc = max ( max( dloc(1), dloc(2) ), dloc(3) )
if ( dmaxc > dpmax ) dpmax = dmaxc
!
! determine centroid of present cell
!
x0 = cell(icell)%attr(CELLCX)
y0 = cell(icell)%attr(CELLCY)
!
! determine depth in centroid in present cell
!
h0 = ( dloc(1) + dloc(2) + dloc(3) )/ 3.
!
! 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)
!
dloc(1) = dep2(v(1))
dloc(2) = dep2(v(2))
dloc(3) = dep2(v(3))
!
if ( dloc(1) <= DEPMIN .or. dloc(2) <= DEPMIN .or. dloc(3) <= DEPMIN ) goto 10
!
! determine centroid of next cell
!
x1 = cell(jcell)%attr(CELLCX)
y1 = cell(jcell)%attr(CELLCY)
!
! determine depth in centroid of next cell
!
h1 = ( dloc(1) + dloc(2) + dloc(3) )/ 3.
!
! compute contribution to area of centroid dual
!
area = area + x0*y1 - x1*y0
!
! compute contribution to x-gradient of depth
!
dhdx = dhdx + ( h0 + h1 ) * real( y1 - y0 )
!
! compute contribution to y-gradient of depth
!
dhdy = dhdy + ( h0 + h1 ) * real( x1 - x0 )
!
enddo
!
! if ratio between max and min depth is too large, set gradients to zero and skip to next vertex
!
!if ( dpmax > drat * dpmin ) goto 10
!
! if area is non-positive, give error and go to next vertex
!
if ( .not. area > 0d0 ) then
!NADC write (msgstr, '(a,i5)') ' Area of centroid dual is negative or zero in vertex ', ivert
!NADC call msgerr( 2, trim(msgstr) )
return
endif
!
dhdx = dhdx/real(area)
dhdy = -dhdy/real(area)
!
! in case of spherical coordinates, transform back to Cartesian coordinates
!
if ( KSPHER > 0 ) then
!
cslat = cos(DEGRAD*(vert(ivert)%attr(VERTY) + YOFFS))
!
dhdx = dhdx/(cslat * LENDEG)
dhdy = dhdy/LENDEG
!
endif
!
goto 20
10 dhdx = 0.
dhdy = 0.
!
endif
!
20 if ( IREFR /= 0 .and. int(PNUMS(32)) == 1 ) then
!
area = 0d0
dpmax = -99999.
dpmin = 99999.
!
! 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)
!
dloc(1) = dep2(v(1))
dloc(2) = dep2(v(2))
dloc(3) = dep2(v(3))
!
if ( dloc(1) <= DEPMIN .or. dloc(2) <= DEPMIN .or. dloc(3) <= DEPMIN ) goto 30
!
dminc = min ( min( dloc(1), dloc(2) ), dloc(3) )
if ( dminc < dpmin ) dpmin = dminc
!
dmaxc = max ( max( dloc(1), dloc(2) ), dloc(3) )
if ( dmaxc > dpmax ) dpmax = dmaxc
!
! determine centroid of present cell
!
x0 = cell(icell)%attr(CELLCX)
y0 = cell(icell)%attr(CELLCY)
!
! compute wave numbers for all frequencies
!
call KSCIP1 (MSC, spcsig, dloc(1), kloc(1,1), arr, arr, arr)
call KSCIP1 (MSC, spcsig, dloc(2), kloc(1,2), arr, arr, arr)
call KSCIP1 (MSC, spcsig, dloc(3), kloc(1,3), arr, arr, arr)
!
if ( IMUD == 1 ) then
!
if (VARMUD) then
dm(1) = mudl2(v(1))
dm(2) = mudl2(v(2))
dm(3) = mudl2(v(3))
else
dm(1) = PMUD(1)
dm(2) = PMUD(1)
dm(3) = PMUD(1)
endif
!
call KSCIP2 (MSC, spcsig, dloc(1), kloc(1,1), arr, arr, arr, arr, dm(1))
call KSCIP2 (MSC, spcsig, dloc(2), kloc(1,2), arr, arr, arr, arr, dm(2))
call KSCIP2 (MSC, spcsig, dloc(3), kloc(1,3), arr, arr, arr, arr, dm(3))
!
endif
!
! determine wave number in centroid in present cell
!
k0(:) = ( kloc(:,1) + kloc(:,2) + kloc(:,3) )/ 3.
!
! 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)
!
dloc(1) = dep2(v(1))
dloc(2) = dep2(v(2))
dloc(3) = dep2(v(3))
!
if ( dloc(1) <= DEPMIN .or. dloc(2) <= DEPMIN .or. dloc(3) <= DEPMIN ) goto 30
!
! determine centroid of next cell
!
x1 = cell(jcell)%attr(CELLCX)
y1 = cell(jcell)%attr(CELLCY)
!
! compute wave numbers for all frequencies
!
call KSCIP1 (MSC, spcsig, dloc(1), kloc(1,1), arr, arr, arr)
call KSCIP1 (MSC, spcsig, dloc(2), kloc(1,2), arr, arr, arr)
call KSCIP1 (MSC, spcsig, dloc(3), kloc(1,3), arr, arr, arr)
!
if ( IMUD == 1 ) then
!
if (VARMUD) then
dm(1) = mudl2(v(1))
dm(2) = mudl2(v(2))
dm(3) = mudl2(v(3))
endif
!
call KSCIP2 (MSC, spcsig, dloc(1), kloc(1,1), arr, arr, arr, arr, dm(1))
call KSCIP2 (MSC, spcsig, dloc(2), kloc(1,2), arr, arr, arr, arr, dm(2))
call KSCIP2 (MSC, spcsig, dloc(3), kloc(1,3), arr, arr, arr, arr, dm(3))
!
endif
!
! determine wave number in centroid of next cell
!
k1(:) = ( kloc(:,1) + kloc(:,2) + kloc(:,3) )/ 3.
!
! compute contribution to area of centroid dual
!
area = area + x0*y1 - x1*y0
!
! compute contribution to x-gradient of wave number
!
dkdx(:) = dkdx(:) + ( k0(:) + k1(:) ) * real( y1 - y0 )
!
! compute contribution to y-gradient of wave number
!
dkdy(:) = dkdy(:) + ( k0(:) + k1(:) ) * real( x1 - x0 )
!
enddo
!
! if ratio between max and min depth is too large, set gradients to zero and skip to next vertex
!
!if ( dpmax > drat * dpmin ) goto 30
!
! if area is non-positive, give error and go to next vertex
!
if ( .not. area > 0d0 ) then
!NADC write (msgstr, '(a,i5)') ' Area of centroid dual is negative or zero in vertex ', ivert
!NADC call msgerr( 2, trim(msgstr) )
return
endif
!
dkdx(:) = dkdx(:)/real(area)
dkdy(:) = -dkdy(:)/real(area)
!
! in case of spherical coordinates, transform back to Cartesian coordinates
!
if ( KSPHER > 0 ) then
!
cslat = cos(DEGRAD*(vert(ivert)%attr(VERTY) + YOFFS))
!
dkdx(:) = dkdx(:)/(cslat * LENDEG)
dkdy(:) = dkdy(:)/LENDEG
!
endif
!
return
30 dkdx = 0.
dkdy = 0.
!
endif
!
end subroutine SwanGradDepthorK