forked from sandialabs/StrideSearch
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathTropicalStrideSearch.f90
370 lines (332 loc) · 14 KB
/
TropicalStrideSearch.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
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
module TropicalStrideSearchModule
!> @file TropicalStrideSearch.f90
!! @brief Extends StrideSearch to the application of tropical cyclone detection.
!> @author Peter Bosler SNL
!!
!> @defgroup TropicalStrideSearch TropicalStrideSearch
!! @brief Extends the Stride Search data structure for tropical cyclone detection
!!
!! @{
use StrideSearchModule
use TropicalDataModule
use StormListNodeModule, only : DEG_2_RAD
use TropicalStormListNodeModule
!use OMP_LIB
implicit none
private
public TropicalStrideSearchSector, TropicalSearchSetup, DoTropicalSearch, FinalizeTropicalSector, PrintTropicalSearchInfo
public RemoveMarkedTropicalNodes, MarkTropicalNodesForRemoval, ApplyTropicalLandMask, GetLandMask
type, extends(StrideSearchSector) :: TropicalStrideSearchSector
real, allocatable, dimension(:,:) :: tempWork
real :: vortPslDistThreshold
real :: tempExcessThreshold
real :: tempPslDistThreshold
end type
contains
!> @brief Defines search sectors based on sectorRadius, links sectors to data.
!!
!> @param tSearch StrideSearch object defined by this subroutine
!> @param southernBoundary southernmost latitude (in degrees_north) of search domain
!> @param northernBoundary northernmost latitude (in degrees_north) of search domain
!> @param sectorRadius spatial scale of a typical storm (in km); a maximum of 1 storm will be found in any
!> spherical circle with geodesic radius = sectorRadius
!> @param pslThreshold storms with minimum pressures greater than this value (in Pa) will be ignored
!> @param vortThreshold storms with cyclonic vorticity (1/s) less than this value will be ignored
!> @param windThreshold storms with maximum winds (m/s) less than this value will be ignored
!> @param vortPslDistThreshold maximum allowable distance (in km) separating a storm's pressure minimum from its
!> vorticity maximum
!> @param tempExcessThreshold storms whose maximum vertically averaged temperature does not exceed the sector average of
!> the vertically average temperature by at least this amount will be ignored
!> @param tempPslDistThreshold maximum allowable distance (in km) separating a storms' pressure minimum from the location
!> of its maximum vertically averaged temperature
!> @param tData data read from netcdf file by @ref tropicaldatamodule::readtropicalvariablesattimestep
subroutine TropicalSearchSetup( tSearch, southernBoundary, northernBoundary, sectorRadius, &
pslThreshold, vortThreshold, windThreshold, &
vortPslDistThreshold, tempExcessThreshold, tempPslDistThreshold, &
tData )
type(TropicalStrideSearchSector), intent(out) :: tSearch
real, intent(in) :: southernBoundary, northernBoundary, sectorRadius
real, intent(in) :: pslThreshold, vortThreshold, windThreshold
real, intent(in) :: vortPslDistThreshold, tempExcessThreshold, tempPslDistThreshold
type(TropicalData), intent(in) :: tData
call SearchSetup(tSearch%StrideSearchSector, southernBoundary, northernBoundary, sectorRadius, &
pslThreshold, vortThreshold, windThreshold, tData%StrideSearchData)
tSearch%vortPslDistThreshold = vortPslDistThreshold
tSearch%tempExcessThreshold = tempExcessThreshold
tSearch%tempPslDistThreshold = tempPslDistThreshold
end subroutine
subroutine AllocateTropicalSectorMemory( tSearch, tData, stripIndex )
type(TropicalStrideSearchSector), intent(inout) :: tSearch
type(TropicalData), intent(in) :: tData
integer, intent(in) :: stripIndex
call AllocateSectorMemory(tSearch%StrideSearchSector, tData%StrideSearchData, stripIndex)
allocate(tSearch%tempWork( size(tSearch%myLatIs), size(tSearch%myLonJs) ) )
end subroutine
!> @brief Frees memory used by a TropicalStrideSearch object which is allocated in TropicalSearchSetup.
!> @param tSearch object to be deleted
subroutine FinalizeTropicalSector( tSearch )
type(TropicalStrideSearchSector), intent(inout) :: tSearch
deallocate(tSearch%tempWork)
call FinalizeSector(tSearch%StrideSearchSector)
end subroutine
!> @brief Prints basic info about a TropicalStrideSearch object to the console
!> @param tSearch
subroutine PrintTropicalSearchInfo( tSearch )
type(TropicalStrideSearchSector), intent(in) :: tSearch
call PrintSearchInfo( tSearch%StrideSearchSector )
print *, "vortPslDistThreshold = ", tSearch%vortPslDistThreshold
print *, "tempExcessThreshold = ", tSearch%tempExcessThreshold
print *, "tempPslDistThreshold = ", tSearch%tempPslDistThreshold
end subroutine
!> @brief Performs a stride search of a data set, identifies storms according to the criteria set in TropicalSearchSetup.
!> @param tstormList list of storms found during this Stride Search
!> @param tSearch Stride Search object defined by TropicalSearchSetup
!> @param tData netcdf data output by ReadTropicalVariablesAtTimestep
subroutine DoTropicalSearch( tstormList, tSearch, tData )
type(TropicalStormListNode), pointer, intent(inout) :: tstormList
type(TropicalStrideSearchSector), intent(inout) :: tSearch
type(TropicalData), intent(in) :: tData
!
integer :: i, j, k, ii, jj, kk, stormI, stormJ
real :: stormLon, stormLat, stormPsl, stormVort, stormWind
real :: stormTemp, tempLon, tempLat, stormThick, thickLon, thickLat
real :: vortLon, vortLat
logical :: hasWarmCore, hasThickness
type(TropicalStormListNode), pointer :: tempNode
logical :: foundStorm, crit1, crit2, crit3, crit4
real :: avgT
integer :: nStrips, tid
real :: dataRes
dataRes = 360.0 / tData%nLon
allocate(tempNode)
!
! loop over search sector centers
!
kk = 0
do ii = 1, size(sectorCenterLats)
call AllocateTropicalSectorMemory(tSearch, tData, ii)
do jj = 1, nLonsPerLatLine(ii)
kk = kk + 1
call DefineSectorInData(tSearch%StrideSearchSector, tData%StrideSearchData, ii, kk)
!
! collect data from sector's neighborhood
!
tSearch%pslWork = 1.0e20
tSearch%windWork = 0.0
tSearch%vortWork = 0.0
tSearch%tempWork = 0.0
do j = 1, size(tSearch%myLonJs)
do i = 1, size(tSearch%myLatIs)
if ( tSearch%neighborhood(i,j) ) then
tSearch%pslWork(i,j) = tData%psl(tSearch%myLonJs(j), tSearch%myLatIs(i))
tSearch%windWork(i,j) = tData%wind(tSearch%myLonJs(j), tSearch%myLatIs(i))
tSearch%vortWork(i,j) = sign(1.0, tSearch%myLats(i)) * &
tData%vorticity(tSearch%myLonJs(j), tSearch%myLatIs(i))
tSearch%tempWork(i,j) = tData%vertAvgT(tSearch%myLonJs(j), tSearch%myLatIs(i))
endif
enddo
enddo
!
! apply storm identification criteria
!
crit1 = .FALSE.
crit2 = .FALSE.
crit3 = .FALSE.
crit4 = .FALSE.
if ( maxval(tSearch%vortWork) > tSearch%vortThreshold ) then
crit1 = .TRUE.
stormPsl = minval(tSearch%pslWork)
stormVort = maxval(tSearch%vortWork)
stormWind = maxval(tSearch%windWork)
!avgT = ArithmeticAverageTemp(tSearch)
avgT = SpatialAverageTemp(tSearch, dataRes)
tSearch%tempWork = tSearch%tempWork - avgT
stormTemp = maxval(tSearch%tempWork)
do j = 1, size(tSearch%myLonJs)
do i = 1, size(tSearch%myLatIs)
if (stormPsl == tSearch%pslWork(i,j) ) then
stormLon = tSearch%myLons(j)
stormJ = tSearch%myLonJs(j)
stormLat = tSearch%myLats(i)
stormI = tSearch%myLatIs(i)
endif
if ( stormVort == tSearch%vortWork(i,j) ) then
vortLon = tSearch%myLons(j)
vortLat = tSearch%myLats(i)
endif
if ( stormTemp == tSearch%tempWork(i,j) ) then
tempLon = tSearch%myLons(j)
tempLat = tSearch%myLats(i)
endif
enddo
enddo
if ( SphereDistance( stormLon * DEG_2_RAD, stormLat * DEG_2_RAD, &
vortLon * DEG_2_RAD, vortLat * DEG_2_RAD ) < tSearch%vortPslDistThreshold ) &
crit2 = .TRUE.
if ( stormTemp > tSearch%tempExcessThreshold ) crit3 = .TRUE.
if ( SphereDistance( stormLon * DEG_2_RAD, stormLat * DEG_2_RAD, &
tempLon * DEG_2_RAD, tempLat * DEG_2_RAD ) < tSearch%tempPslDistThreshold ) &
crit4 = .TRUE.
endif
foundStorm = ( ( crit1 .AND. crit2) .AND. ( crit3 .AND. crit4) )
hasWarmCore = ( crit3 .AND. crit4 )
! thickness criteria not used currently, but is required by output format
hasThickness = .TRUE.
stormThick = 10000.0
thickLon = 0.0
thickLat = 0.0
if ( foundStorm ) then
call initializeTropical( tempNode, stormLon, stormLat, stormJ, stormI, stormPsl, stormVort, stormWind, &
vortLon, vortLat, stormTemp, tempLon, tempLat, hasWarmCore, stormThick, &
thickLon, thickLat, hasThickness )
call AddTropicalNodeToList( tstormList, tempNode )
endif
enddo
call FinalizeTropicalSector(tSearch)
enddo
deallocate(tempNode)
end subroutine
!> @brief Marks duplicate detections of the same storm, and storms on our outside the search domain boundaries for removal.
!> @param stormList list of possible storms output by DoTropicalSearch
!> @param tSearch TropicalStrideSearch object output by TropicalSearchSetup
!> @param southernBoundary southernmost latitude (in degrees_north) of search domain
!> @param northernBoundary northernmost latitude (in degrees_north) of search domain
subroutine MarkTropicalNodesForRemoval( stormList, tSearch, southernBoundary, northernBoundary )
type(TropicalStormListNode), pointer, intent(inout) :: stormList
type(TropicalStrideSearchSector), intent(in) :: tSearch
real, intent(in) :: southernBoundary, northernBoundary
type(TropicalStormListNode), pointer :: current, next, query, querynext
current => stormList
do while (associated(current) )
next => current%nextTropical
if ( .NOT. current%removeThisNode ) then
if ( current%lat <= southernBoundary .OR. current%lat >= northernBoundary ) then
current%removeThisNode = .TRUE.
else
query => next
do while ( associated( query ) )
querynext => query%nextTropical
if ( SphereDistance( current%lon*DEG_2_RAD, current%lat*DEG_2_RAD, &
query%lon*DEG_2_RAD, query%lat*DEG_2_RAD ) < tSearch%sectorRadius ) then
if ( query%wind > current%wind ) then
current%removeThisNode = .TRUE.
query%psl = min( current%psl, query%psl)
query%vort = max( current%vort, query%vort)
query%vertAvgT = max(current%vertAvgT, query%vertAvgT)
else
query%removeThisNode = .TRUE.
current%psl = min( current%psl, query%psl)
current%vort = max( current%vort, query%vort)
current%vertAvgT = max(current%vertAvgT, query%vertAvgT)
endif
endif
query => querynext
enddo
endif
endif
current => next
enddo
end subroutine
!> @brief Deletes any TropicalStormListNode marked for removal by another subroutine.
!> @param stormList
subroutine RemoveMarkedTropicalNodes(stormList)
type(TropicalStormListNode), pointer, intent(inout) :: stormList
type(TropicalStormListNode), pointer :: current, next
current => stormList
do while ( associated(current) )
next => current%nextTropical
if ( current%removeThisNode ) then
call RemoveTropicalNodeFromList(stormList, current)
endif
current => next
enddo
end subroutine
!> @brief Marks any storms over land for removal
!> @warning This routine applies to all storms, even those that may have originated over water and moved over land.
!> @param stormlist
subroutine ApplyTropicalLandMask( stormList )
type(TropicalStormListNode), pointer, intent(inout) :: stormList
character(len=256), parameter :: maskfile = "/Users/pabosle/Desktop/StrideSearch/gfdlUtilities/landsea.map"
integer, parameter :: maskNLat = 180, &
maskNLon = 360, &
maskFileRows = 1620, &
maskFileColumns = 40
integer :: mask0(maskFileColumns,maskFileRows), mask(maskNLon,maskNLat)
integer :: i, j
real :: dLon, dLat, lon0, lat0
type(TropicalStormListNode), pointer :: current, next
open(unit=12, file=trim(maskfile), status='OLD', action='READ')
do j= 1, maskFileRows
read(12,'(40I2)') mask0(:,j)
enddo
close(12)
mask = RESHAPE(mask0, (/ maskNLon, maskNLat /) )
mask = CSHIFT( mask, maskNLon/2, 1 )
lon0 = 0.0
lat0 = -90.0 + 0.5
dlon = 360.0 / float(maskNLon)
dLat = -2.0 * lat0 / float(maskNLat - 1)
current => stormList
do while ( associated(current) )
next => current%nextTropical
j = (current%lat - lat0)/dLat + 1.5
i = (current%lon - lon0)/dLon + 1.5
if ( i == 0 ) i = maskNLon
if ( i > maskNLon ) i = i - maskNLon
if ( mask(i,j) /= 0 ) then
current%removeThisNode = .TRUE.
endif
current=>next
enddo
end subroutine
!> @brief Computes an arithemetic average of a sector's vertically averaged temperature
!> @param tSearch contains data to be averaged
!> @return ArithmeticAverageTemp
pure function ArithmeticAverageTemp( tSearch )
real :: ArithmeticAverageTemp
type(TropicalStrideSearchSector), intent(in) :: tSearch
!
real :: tsum
tsum = sum( tSearch%tempWork, MASK=tSearch%neighborhood)
ArithmeticAverageTemp = tsum / count(tSearch%neighborhood)
end function
!> @brief Computes a spatial average of a sector's vertically averaged temperature by approximation
!> with midpoint rule quadrature.
!> @param tSearch sector that contains data to be averaged
!> @param dataRes resolution of the data in degrees
!> @return AvgT @f$ = \frac{sum_{K_{IJ}} \overline{T}(\lambda_j,theta_i)\cos\theta_i \Delta\lambda } { \sum_{K_{IJ}} \cos\theta_i\Delta\lambda@f$
pure function SpatialAverageTemp( tSearch, dataRes )
real :: SpatialAverageTemp
type(TropicalStrideSearchSector), intent(in) :: tSearch
real, intent(in) :: dataRes
!
real :: area
real :: tsum
integer :: i, j
area = 0.0
tsum = 0.0
do j = 1, size(tSearch%myLonJs)
do i = 1, size(tSearch%myLatIs)
if ( tSearch%neighborhood(i,j) ) then
area = area + cos( tSearch%myLats(i) * DEG_2_RAD ) * dataRes * DEG_2_RAD
tsum = tsum + tSearch%tempWork(i,j) * cos( tSearch%myLats(i) * DEG_2_RAD ) * dataRes * DEG_2_RAD
endif
enddo
enddo
end function
!> @brief Returns an integer array used to define land masks
function GetLandMask()
real :: GetLandMask(360,180)
character(len=256), parameter :: maskfile = "/Users/pabosle/Desktop/StrideSearch/gfdlUtilities/landsea.map"
integer :: mask0(40,1620)
integer :: j
open(unit=12, file=maskfile, status='OLD', action='READ')
do j = 1, 1620
read(12, '(40I2)') mask0(:,j)
enddo
close(12)
GetLandMask = RESHAPE(mask0, (/ 360, 180 /) )
GetLandMask = CSHIFT( GetLandMask, 180, 1 )
end function
!> @}
end module