forked from NOAA-PSL/stochastic_physics
-
Notifications
You must be signed in to change notification settings - Fork 0
/
four_to_grid_stochy.F
280 lines (256 loc) · 10.9 KB
/
four_to_grid_stochy.F
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
!>@brief The module 'four_to_grid_mod' contains the subroute four_to_grid
module four_to_grid_mod
use spectral_layout_mod, only: num_parthds_stochy => ompthreads
implicit none
contains
!>@brief The subroutine 'epslon_stochy' calculate coeffients for use in spectral space
!>@details This code is taken from the legacy spectral GFS
subroutine four_to_grid(syn_gr_a_1,syn_gr_a_2,
& lon_dim_coef,lon_dim_grid,lons_lat,lot)
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
use kinddef
implicit none
!!
real(kind=kind_dbl_prec) syn_gr_a_1(lon_dim_coef,lot)
real(kind=kind_dbl_prec) syn_gr_a_2(lon_dim_grid,lot)
integer lon_dim_coef
integer lon_dim_grid
integer lons_lat
integer lot
!________________________________________________________
#ifdef MKL
integer*8 plan
#else
real(kind=kind_dbl_prec) aux1crs(42002)
real(kind=kind_dbl_prec) scale_ibm
integer ibmsign
integer init
#endif
integer lot_thread
integer num_threads
integer nvar_thread_max
integer nvar_1
integer nvar_2
integer thread
#ifdef MKL
include "fftw/fftw3.f"
integer NULL
#else
external dcrft
external scrft
#endif
!________________________________________________________
num_threads = min(num_parthds_stochy,lot)
nvar_thread_max = (lot+num_threads-1)/num_threads
if ( kind_dbl_prec == 8 ) then !------------------------------------
#ifdef MKL
!$omp parallel do num_threads(num_threads)
!$omp+shared(syn_gr_a_1,syn_gr_a_2,lons_lat)
!$omp+shared(lon_dim_coef,lon_dim_grid)
!$omp+shared(lot,num_threads,nvar_thread_max)
!$omp+private(thread,nvar_1,nvar_2,lot_thread,plan)
#else
!$omp parallel do num_threads(num_threads)
!$omp+shared(syn_gr_a_1,syn_gr_a_2,lons_lat)
!$omp+shared(lon_dim_coef,lon_dim_grid)
!$omp+shared(lot,num_threads,nvar_thread_max)
!$omp+shared(ibmsign,scale_ibm)
!$omp+private(thread,nvar_1,nvar_2,lot_thread,init,aux1crs)
#endif
do thread=1,num_threads ! start of thread loop ..............
nvar_1=(thread-1)*nvar_thread_max + 1
nvar_2=min(nvar_1+nvar_thread_max-1,lot)
lot_thread=nvar_2 - nvar_1 + 1
if (nvar_2 >= nvar_1) then
#ifdef MKL
!call dfftw_plan_many_dft_c2r(
! plan, 1, N, m, &
! X, NULL, 1, dimx, &
! Y, NULL, 1, dimy, &
! fftw_flag)
call dfftw_plan_many_dft_c2r( &
& plan, 1, lons_lat, lot_thread, &
& syn_gr_a_1, NULL, 1, size(syn_gr_a_1,dim=1), &
& syn_gr_a_2, NULL, 1, size(syn_gr_a_2,dim=1), &
& FFTW_ESTIMATE)
call dfftw_execute(plan)
call dfftw_destroy_plan(plan)
#else
init = 1
ibmsign = -1
scale_ibm = 1.0d0
call dcrft(init,
& syn_gr_a_1(1,nvar_1) ,lon_dim_coef/2,
& syn_gr_a_2(1,nvar_1) ,lon_dim_grid,
& lons_lat,lot_thread,ibmsign,scale_ibm,
& aux1crs,22000,
& aux1crs(22001),20000)
init = 0
call dcrft(init,
& syn_gr_a_1(1,nvar_1) ,lon_dim_coef/2,
& syn_gr_a_2(1,nvar_1) ,lon_dim_grid,
& lons_lat,lot_thread,ibmsign,scale_ibm,
& aux1crs,22000,
& aux1crs(22001),20000)
#endif
endif
enddo ! fin thread loop ......................................
else !------------------------------------------------------------
#ifdef MKL
!$omp parallel do num_threads(num_threads)
!$omp+shared(syn_gr_a_1,syn_gr_a_2,lons_lat)
!$omp+shared(lon_dim_coef,lon_dim_grid)
!$omp+shared(lot,num_threads,nvar_thread_max)
!$omp+private(thread,nvar_1,nvar_2,lot_thread,plan)
#else
!$omp parallel do num_threads(num_threads)
!$omp+shared(syn_gr_a_1,syn_gr_a_2,lons_lat)
!$omp+shared(lon_dim_coef,lon_dim_grid)
!$omp+shared(lot,num_threads,nvar_thread_max)
!$omp+shared(ibmsign,scale_ibm)
!$omp+private(thread,nvar_1,nvar_2,lot_thread,init,aux1crs)
#endif
do thread=1,num_threads ! start of thread loop ..............
nvar_1 = (thread-1)*nvar_thread_max + 1
nvar_2 = min(nvar_1+nvar_thread_max-1,lot)
lot_thread = nvar_2 - nvar_1 + 1
if (nvar_2 >= nvar_1) then
#ifdef MKL
!call sfftw_plan_many_dft_c2r(
! plan, 1, N, m, &
! X, NULL, 1, dimx, &
! Y, NULL, 1, dimy, &
! fftw_flag)
call sfftw_plan_many_dft_c2r( &
& plan, 1, lons_lat, lot_thread, &
& syn_gr_a_1, NULL, 1, size(syn_gr_a_1,dim=1), &
& syn_gr_a_2, NULL, 1, size(syn_gr_a_2,dim=1), &
& FFTW_ESTIMATE)
call sfftw_execute(plan)
call sfftw_destroy_plan(plan)
#else
init = 1
ibmsign = -1
scale_ibm = 1.0d0
call scrft(init,
& syn_gr_a_1(1,nvar_1) ,lon_dim_coef/2,
& syn_gr_a_2(1,nvar_1) ,lon_dim_grid,
& lons_lat,lot_thread,ibmsign,scale_ibm,
& aux1crs,22000,
& aux1crs(22001),20000,
& aux1crs(22001),0)
init = 0
call scrft(init,
& syn_gr_a_1(1,nvar_1) ,lon_dim_coef/2,
& syn_gr_a_2(1,nvar_1) ,lon_dim_grid,
& lons_lat,lot_thread,ibmsign,scale_ibm,
& aux1crs,22000,
& aux1crs(22001),20000,
& aux1crs(22001),0)
#endif
endif
enddo ! fin thread loop ......................................
endif !-----------------------------------------------------------
!!
return
end
subroutine grid_to_four(anl_gr_a_2,anl_gr_a_1,
& lon_dim_grid,lon_dim_coef,lons_lat,lot)
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
use kinddef
implicit none
!!
real(kind=kind_dbl_prec) anl_gr_a_2(lon_dim_grid,lot)
real(kind=kind_dbl_prec) anl_gr_a_1(lon_dim_coef,lot)
integer lon_dim_grid
integer lon_dim_coef
integer lons_lat
integer lot
!________________________________________________________
real(kind=kind_dbl_prec) aux1crs(42002)
real(kind=kind_dbl_prec) scale_ibm,rone
integer ibmsign
integer init
integer lot_thread
integer num_threads
integer nvar_thread_max
integer nvar_1,nvar_2
integer thread
!________________________________________________________
#ifdef MKL
write(0,*) "ERROR in grid_to_four: srcft and drcft ",
& " must be replaced with MKL's FFTW calls. ABORT."
call sleep(5)
stop
#endif
num_threads=min(num_parthds_stochy,lot)
nvar_thread_max=(lot+num_threads-1)/num_threads
if ( kind_dbl_prec == 8 ) then !------------------------------------
!$omp parallel do num_threads(num_threads)
!$omp+shared(anl_gr_a_1,anl_gr_a_2,lons_lat)
!$omp+shared(lon_dim_coef,lon_dim_grid)
!$omp+shared(lot,num_threads,nvar_thread_max)
!$omp+shared(ibmsign,scale_ibm,rone)
!$omp+private(thread,nvar_1,nvar_2,lot_thread,init,aux1crs)
do thread=1,num_threads ! start of thread loop ..............
nvar_1 = (thread-1)*nvar_thread_max + 1
nvar_2 = min(nvar_1+nvar_thread_max-1,lot)
if (nvar_2 >= nvar_1) then
lot_thread = nvar_2 - nvar_1 + 1
init = 1
ibmsign = 1
rone = 1.0d0
scale_ibm = rone/lons_lat
call drcft(init,
& anl_gr_a_2(1,nvar_1), lon_dim_grid,
& anl_gr_a_1(1,nvar_1), lon_dim_coef/2,
& lons_lat,lot_thread,ibmsign,scale_ibm,
& aux1crs,22000,
& aux1crs(22001),20000)
init = 0
call drcft(init,
& anl_gr_a_2(1,nvar_1), lon_dim_grid,
& anl_gr_a_1(1,nvar_1), lon_dim_coef/2,
& lons_lat,lot_thread,ibmsign,scale_ibm,
& aux1crs,22000,
& aux1crs(22001),20000)
endif
enddo ! fin thread loop ......................................
else !------------------------------------------------------------
!$omp parallel do num_threads(num_threads)
!$omp+shared(anl_gr_a_1,anl_gr_a_2,lons_lat)
!$omp+shared(lon_dim_coef,lon_dim_grid)
!$omp+shared(lot,num_threads,nvar_thread_max)
!$omp+shared(ibmsign,scale_ibm,rone)
!$omp+private(thread,nvar_1,nvar_2,lot_thread,init,aux1crs)
do thread=1,num_threads ! start of thread loop ..............
nvar_1 = (thread-1)*nvar_thread_max + 1
nvar_2 = min(nvar_1+nvar_thread_max-1,lot)
if (nvar_2 >= nvar_1) then
lot_thread=nvar_2 - nvar_1 + 1
init = 1
ibmsign = 1
rone = 1.0d0
scale_ibm = rone/lons_lat
call srcft(init,
& anl_gr_a_2(1,nvar_1), lon_dim_grid,
& anl_gr_a_1(1,nvar_1), lon_dim_coef/2,
& lons_lat,lot_thread,ibmsign,scale_ibm,
& aux1crs,22000,
& aux1crs(22001),20000,
& aux1crs(22001),0)
init = 0
call srcft(init,
& anl_gr_a_2(1,nvar_1), lon_dim_grid,
& anl_gr_a_1(1,nvar_1), lon_dim_coef/2,
& lons_lat,lot_thread,ibmsign,scale_ibm,
& aux1crs,22000,
& aux1crs(22001),20000,
& aux1crs(22001),0)
endif
enddo ! fin thread loop ......................................
endif !-----------------------------------------------------------
!!
return
end
end module four_to_grid_mod