-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathutils_mod.F90
542 lines (444 loc) · 18.9 KB
/
utils_mod.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
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
module utils_mod
use netcdf
use ocnvars, only : vardefs
implicit none
private
interface getfield
module procedure getfield2d
module procedure getfield3d
end interface getfield
interface packarrays
module procedure packarrays2d
module procedure packarrays3d
end interface packarrays
interface remap
module procedure remap1d
module procedure remap2d
module procedure remap3d
end interface remap
interface getvecpair
module procedure getvecpair2d
module procedure getvecpair3d
end interface getvecpair
interface dumpnc
module procedure dumpnc2d
module procedure dumpnc3d
end interface dumpnc
public getfield
public packarrays
public remap
public dumpnc
logical :: debug = .true.
contains
!----------------------------------------------------------
! pack 2D fields into arrays by mapping type
!----------------------------------------------------------
subroutine packarrays2d(filesrc, wgtsdir, cosrot, sinrot, vars, dims, nflds, fields)
character(len=*), intent(in) :: filesrc,wgtsdir
real, intent(in) :: cosrot(:),sinrot(:)
type(vardefs), intent(in) :: vars(:)
integer, intent(in) :: dims(:)
integer, intent(in) :: nflds
real, intent(out) :: fields(:,:)
! local variables
integer :: n, nn
real, allocatable, dimension(:,:) :: vecpair
character(len=20) :: subname = 'packarrays2d'
fields=0.0
if (debug)print '(a)','enter '//trim(subname)
! obtain vector pairs
do n = 1,nflds
if (trim(vars(n)%var_grid) == 'Cu') then
allocate(vecpair(dims(1)*dims(2),2)); vecpair = 0.0
call getvecpair(trim(filesrc), trim(wgtsdir), cosrot, sinrot, &
trim(vars(n)%input_var_name), trim(vars(n)%var_grid), &
trim(vars(n)%var_pair), trim(vars(n)%var_pair_grid), &
dims=(/dims(1),dims(2)/), vecpair=vecpair)
end if
end do
! create packed array
nn = 0
do n = 1,nflds
if (len_trim(vars(n)%var_pair) == 0) then
nn = nn + 1
call getfield(trim(filesrc), trim(vars(n)%input_var_name), dims=(/dims(1),dims(2)/), &
field=fields(:,nn))
else ! fill with vector pairs
nn = nn+1
if (trim(vars(n)%var_grid) == 'Cu')fields(:,nn) = vecpair(:,1)
if (trim(vars(n)%var_grid) == 'Cv')fields(:,nn) = vecpair(:,2)
end if
end do
if (debug)print '(a)','exit '//trim(subname)
end subroutine packarrays2d
!----------------------------------------------------------
! pack 3D fields into arrays by mapping type
!----------------------------------------------------------
subroutine packarrays3d(filesrc, wgtsdir, cosrot, sinrot, vars, dims, nflds, fields)
character(len=*), intent(in) :: filesrc,wgtsdir
real, intent(in) :: cosrot(:),sinrot(:)
type(vardefs), intent(in) :: vars(:)
integer, intent(in) :: dims(:)
integer, intent(in) :: nflds
real, intent(out) :: fields(:,:,:)
! local variables
integer :: n, nn
real, allocatable, dimension(:,:,:) :: vecpair
character(len=20) :: subname = 'packarrays3d'
fields=0.0
if (debug)print '(a)','enter '//trim(subname)
! obtain vector pairs
do n = 1,dims(3)
if (trim(vars(n)%var_grid) == 'Cu') then
allocate(vecpair(dims(1)*dims(2),dims(3),2)); vecpair = 0.0
call getvecpair(trim(filesrc), trim(wgtsdir), cosrot, sinrot, &
trim(vars(n)%input_var_name), trim(vars(n)%var_grid), &
trim(vars(n)%var_pair), trim(vars(n)%var_pair_grid), &
dims=(/dims(1),dims(2),dims(3)/), vecpair=vecpair)
end if
end do
! create packed array
nn = 0
do n = 1,nflds
if (len_trim(vars(n)%var_pair) == 0) then
nn = nn + 1
call getfield(trim(filesrc), trim(vars(n)%input_var_name), dims=(/dims(1),dims(2),dims(3)/), &
field=fields(:,:,nn))
else ! fill with vector pairs
nn = nn+1
if (trim(vars(n)%var_grid) == 'Cu')fields(:,:,nn) = vecpair(:,:,1)
if (trim(vars(n)%var_grid) == 'Cv')fields(:,:,nn) = vecpair(:,:,2)
end if
end do
if (debug)print '(a)','exit '//trim(subname)
end subroutine packarrays3d
!----------------------------------------------------------
! obtain 2D vector pairs mapped to Ct and rotated to EW
!----------------------------------------------------------
subroutine getvecpair2d(fname, wdir, cosrot, sinrot, vname1, vgrid1, &
vname2, vgrid2, dims, vecpair)
character(len=*), intent(in) :: fname
character(len=*), intent(in) :: wdir
real, intent(in) :: cosrot(:), sinrot(:)
character(len=*), intent(in) :: vname1, vgrid1, vname2, vgrid2
integer, intent(in) :: dims(:)
real, intent(out) :: vecpair(:,:)
! local variables
integer :: ii
real, dimension(dims(1)*dims(2)) :: urot, vrot
character(len=240) :: wgtsfile
character(len=20) :: subname = 'getvecpair2d'
if (debug)print '(a)','enter '//trim(subname)
wgtsfile = trim(wdir)//'tripole.mx025.'//vgrid1//'.to.Ct.bilinear.nc'
call getfield(fname, vname1, dims=dims, field=vecpair(:,1), wgts=trim(wgtsfile))
wgtsfile = trim(wdir)//'tripole.mx025.'//vgrid2//'.to.Ct.bilinear.nc'
call getfield(fname, vname2, dims=dims, field=vecpair(:,2), wgts=trim(wgtsfile))
urot = 0.0; vrot = 0.0
do ii = 1,dims(1)*dims(2)
urot(ii) = vecpair(ii,1)*cosrot(ii) + vecpair(ii,2)*sinrot(ii)
vrot(ii) = vecpair(ii,2)*cosrot(ii) - vecpair(ii,1)*sinrot(ii)
end do
vecpair(:,1) = urot(:)
vecpair(:,2) = vrot(:)
if (debug) print '(a)','exit '//trim(subname)
end subroutine getvecpair2d
!----------------------------------------------------------
! obtain 3D vector pairs, mapped to Ct and rotated to EW
!----------------------------------------------------------
subroutine getvecpair3d(fname, wdir, cosrot, sinrot, vname1, vgrid1, &
vname2, vgrid2, dims, vecpair)
character(len=*), intent(in) :: fname
character(len=*), intent(in) :: wdir
real, intent(in) :: cosrot(:), sinrot(:)
character(len=*), intent(in) :: vname1, vgrid1, vname2, vgrid2
integer, intent(in) :: dims(:)
real, intent(out) :: vecpair(:,:,:)
! local variables
integer :: ii,k
real, dimension(dims(1)*dims(2)) :: urot, vrot
character(len=240) :: wgtsfile
character(len=20) :: subname = 'getfield3d'
if (debug)print '(a)','enter '//trim(subname)
wgtsfile = trim(wdir)//'tripole.mx025.'//vgrid1//'.to.Ct.bilinear.nc'
call getfield(fname, vname1, dims=dims, field=vecpair(:,:,1), wgts=trim(wgtsfile))
wgtsfile = trim(wdir)//'tripole.mx025.'//vgrid2//'.to.Ct.bilinear.nc'
call getfield(fname, vname2, dims=dims, field=vecpair(:,:,2), wgts=trim(wgtsfile))
do k = 1,dims(3)
urot = 0.0; vrot = 0.0
do ii = 1,dims(1)*dims(2)
urot(ii)= vecpair(ii,k,1)*cosrot(ii) + vecpair(ii,k,2)*sinrot(ii)
vrot(ii)= vecpair(ii,k,2)*cosrot(ii) - vecpair(ii,k,1)*sinrot(ii)
end do
vecpair(:,k,1) = urot(:)
vecpair(:,k,2) = vrot(:)
end do
if (debug) print '(a)','exit '//trim(subname)
end subroutine getvecpair3d
!----------------------------------------------------------
! obtain a 2D field and return a 1-D vector array
!----------------------------------------------------------
subroutine getfield2d(fname, vname, dims, field, wgts)
character(len=*), intent(in) :: fname, vname
integer, intent(in) :: dims(:)
real, intent(out) :: field(:)
character(len=*), optional, intent(in) :: wgts
! local variable
integer :: ncid, varid, rc
real :: fval
real, allocatable :: a2d(:,:)
real, allocatable :: atmp(:)
character(len=20) :: subname = 'getfield2d'
if (debug)print '(a)','enter '//trim(subname)//' variable '//vname
allocate(a2d(dims(1),dims(2))); a2d = 0.0
allocate(atmp(dims(1)*dims(2))); atmp = 0.0
rc = nf90_open(fname, nf90_nowrite, ncid)
call handle_err(rc,' nf90_open '//fname)
rc = nf90_inq_varid(ncid, vname, varid)
call handle_err(rc,' get variable ID '// vname)
rc = nf90_get_var(ncid, varid, a2d)
call handle_err(rc,' get variable'// vname)
rc = nf90_get_att(ncid, varid, '_FillValue', fval)
rc = nf90_close(ncid)
atmp(:) = reshape(a2d, (/dims(1)*dims(2)/))
if(present(wgts)) then
where(atmp .eq. fval)atmp = 0.0
call remap(trim(wgts), src_field=atmp, dst_field=field)
else
field = atmp
end if
if (debug) print '(a)','exit '//trim(subname)//' variable '//vname
end subroutine getfield2d
!----------------------------------------------------------
! obtain a 3D field and return a 2-D vector array
!----------------------------------------------------------
subroutine getfield3d(fname, vname, dims, field, wgts)
character(len=*), intent(in) :: fname, vname
integer, intent(in) :: dims(:)
real, intent(out) :: field(:,:)
character(len=*), optional, intent(in) :: wgts
! local variable
integer :: ncid, varid, rc
real :: fval
real, allocatable :: a3d(:,:,:)
real, allocatable :: atmp(:,:)
character(len=20) :: subname = 'getfield3d'
if (debug)print '(a)','enter '//trim(subname)//' variable '//vname
allocate(a3d(dims(1),dims(2),dims(3))); a3d = 0.0
allocate(atmp(dims(1)*dims(2),dims(3))); atmp = 0.0
rc = nf90_open(fname, nf90_nowrite, ncid)
call handle_err(rc,' nf90_open '//fname)
rc = nf90_inq_varid(ncid, vname, varid)
call handle_err(rc,' get variable ID '// vname)
rc = nf90_get_var(ncid, varid, a3d)
call handle_err(rc,' get variable'// vname)
rc = nf90_get_att(ncid, varid, '_FillValue', fval)
rc = nf90_close(ncid)
atmp(:,:) = reshape(a3d, (/dims(1)*dims(2),dims(3)/))
if(present(wgts)) then
where(atmp .eq. fval)atmp = 0.0
call remap(trim(wgts), dim3=dims(3), src_field=atmp, dst_field=field)
else
field = atmp
end if
if (debug) print '(a)','exit '//trim(subname)//' variable '//vname
end subroutine getfield3d
!----------------------------------------------------------
! remap a 1-D vector array
!----------------------------------------------------------
subroutine remap1d(fname, src_field, dst_field)
character(len=*), intent(in) :: fname
real, intent(in) :: src_field(:)
real, intent(out) :: dst_field(:)
! local variables
integer :: ncid, rc, id
integer :: i,ii,jj
integer :: n_a, n_b, n_s
! not real, but get_var(ncid, id, col) seg faults if col is int
! see https://github.com/Unidata/netcdf-fortran/issues/413
real(kind=4), allocatable, dimension(:) :: col, row
real(kind=8), allocatable, dimension(:) :: S
character(len=20) :: subname = 'remap2d'
if (debug)print '(a)','enter '//trim(subname)
! retrieve the weights
rc = nf90_open(trim(fname), nf90_nowrite, ncid)
call handle_err(rc,' nf90_open '//fname)
rc = nf90_inq_dimid(ncid, 'n_s', id)
rc = nf90_inquire_dimension(ncid, id, len=n_s)
rc = nf90_inq_dimid(ncid, 'n_a', id)
rc = nf90_inquire_dimension(ncid, id, len=n_a)
rc = nf90_inq_dimid(ncid, 'n_b', id)
rc = nf90_inquire_dimension(ncid, id, len=n_b)
allocate(col(1:n_s)); col = 0.0
allocate(row(1:n_s)); row = 0.0
allocate( S(1:n_s)); S = 0.0
rc = nf90_inq_varid(ncid, 'col', id)
rc = nf90_get_var(ncid, id, col)
rc = nf90_inq_varid(ncid, 'row', id)
rc = nf90_get_var(ncid, id, row)
rc = nf90_inq_varid(ncid, 'S', id)
rc = nf90_get_var(ncid, id, S)
rc = nf90_close(ncid)
dst_field = 0.0
do i = 1,n_s
ii = row(i); jj = col(i)
dst_field(ii) = dst_field(ii) + S(i)*real(src_field(jj),8)
enddo
if (debug) print '(a)','exit '//trim(subname)
end subroutine remap1d
!----------------------------------------------------------
! remap a packed field of either nflds or nlevs
!----------------------------------------------------------
subroutine remap2d(fname, dim2, src_field, dst_field)
character(len=*), intent(in) :: fname
integer, intent(in) :: dim2
real, intent(in) :: src_field(:,:)
real, intent(out) :: dst_field(:,:)
! local variables
integer :: ncid, rc, id
integer :: i,ii,jj
integer :: n_a, n_b, n_s
! not real, but get_var(ncid, id, col) seg faults if col is int
! see https://github.com/Unidata/netcdf-fortran/issues/413
real(kind=4), allocatable, dimension(:) :: col, row
real(kind=8), allocatable, dimension(:) :: S
character(len=20) :: subname = 'remap3d'
if (debug)print '(a)','enter '//trim(subname)//' weights = '//trim(fname)
! retrieve the weights
rc = nf90_open(trim(fname), nf90_nowrite, ncid)
call handle_err(rc,' nf90_open '//fname)
rc = nf90_inq_dimid(ncid, 'n_s', id)
rc = nf90_inquire_dimension(ncid, id, len=n_s)
rc = nf90_inq_dimid(ncid, 'n_a', id)
rc = nf90_inquire_dimension(ncid, id, len=n_a)
rc = nf90_inq_dimid(ncid, 'n_b', id)
rc = nf90_inquire_dimension(ncid, id, len=n_b)
allocate(col(1:n_s)); col = 0.0
allocate(row(1:n_s)); row = 0.0
allocate( S(1:n_s)); S = 0.0
rc = nf90_inq_varid(ncid, 'col', id)
rc = nf90_get_var(ncid, id, col)
rc = nf90_inq_varid(ncid, 'row', id)
rc = nf90_get_var(ncid, id, row)
rc = nf90_inq_varid(ncid, 'S', id)
rc = nf90_get_var(ncid, id, S)
rc = nf90_close(ncid)
dst_field = 0.0
do i = 1,n_s
ii = row(i); jj = col(i)
dst_field(ii,:) = dst_field(ii,:) + S(i)*real(src_field(jj,:),8)
enddo
if (debug) print '(a)','exit '//trim(subname)
end subroutine remap2d
!----------------------------------------------------------
! remap a field packed array of nk levels and nflds fields
!----------------------------------------------------------
subroutine remap3d(fname, nk, nflds, src_field, dst_field)
character(len=*), intent(in) :: fname
integer, intent(in) :: nk, nflds
real, intent(in) :: src_field(:,:,:)
real, intent(out) :: dst_field(:,:,:)
! local variables
integer :: ncid, rc, id
integer :: i,ii,jj
integer :: n_a, n_b, n_s
! not real, but get_var(ncid, id, col) seg faults if col is int
! see https://github.com/Unidata/netcdf-fortran/issues/413
real(kind=4), allocatable, dimension(:) :: col, row
real(kind=8), allocatable, dimension(:) :: S
character(len=20) :: subname = 'remap4d'
if (debug)print '(a)','enter '//trim(subname)//' weights = '//trim(fname)
! retrieve the weights
rc = nf90_open(trim(fname), nf90_nowrite, ncid)
call handle_err(rc,' nf90_open '//fname)
rc = nf90_inq_dimid(ncid, 'n_s', id)
rc = nf90_inquire_dimension(ncid, id, len=n_s)
rc = nf90_inq_dimid(ncid, 'n_a', id)
rc = nf90_inquire_dimension(ncid, id, len=n_a)
rc = nf90_inq_dimid(ncid, 'n_b', id)
rc = nf90_inquire_dimension(ncid, id, len=n_b)
allocate(col(1:n_s)); col = 0.0
allocate(row(1:n_s)); row = 0.0
allocate( S(1:n_s)); S = 0.0
rc = nf90_inq_varid(ncid, 'col', id)
rc = nf90_get_var(ncid, id, col)
rc = nf90_inq_varid(ncid, 'row', id)
rc = nf90_get_var(ncid, id, row)
rc = nf90_inq_varid(ncid, 'S', id)
rc = nf90_get_var(ncid, id, S)
rc = nf90_close(ncid)
dst_field = 0.0
do i = 1,n_s
ii = row(i); jj = col(i)
dst_field(ii,:,:) = dst_field(ii,:,:) + S(i)*real(src_field(jj,:,:),8)
enddo
if (debug) print '(a)','exit '//trim(subname)
end subroutine remap3d
!----------------------------------------------------------
! write a bare netcdf file of a 2D packed field
!----------------------------------------------------------
subroutine dumpnc2d(fname, vname, dims, nflds, field)
character(len=*), intent(in) :: fname, vname
integer, intent(in) :: dims(:)
integer, intent(in) :: nflds
real, intent(in) :: field(:,:)
! local variable
integer :: ncid, varid, rc, idimid, jdimid, fdimid
real, allocatable :: a3d(:,:,:)
character(len=20) :: subname = 'dumpnc2d'
if (debug)print '(a)','enter '//trim(subname)//' variable '//vname
print *,dims
print *,size(field,1),size(field,2)
allocate(a3d(dims(1),dims(2),nflds)); a3d = 0.0
rc = nf90_create(trim(fname), nf90_clobber, ncid)
rc = nf90_def_dim(ncid, 'nx', dims(1), idimid)
rc = nf90_def_dim(ncid, 'ny', dims(2), jdimid)
rc = nf90_def_dim(ncid, 'nf', dims(3), fdimid)
rc = nf90_def_var(ncid, vname, nf90_float, (/idimid,jdimid,fdimid/), varid)
rc = nf90_enddef(ncid)
a3d(:,:,:) = reshape(field(1:dims(1)*dims(2),1:nflds), (/dims(1),dims(2),nflds/))
rc = nf90_put_var(ncid, varid, a3d)
rc = nf90_close(ncid)
if (debug)print '(a)','exit '//trim(subname)//' variable '//vname
end subroutine dumpnc2d
!----------------------------------------------------------
! write a bare netcdf file of a packed 3D field
!----------------------------------------------------------
subroutine dumpnc3d(fname, vname, dims, nflds, field)
character(len=*), intent(in) :: fname, vname
integer, intent(in) :: dims(:)
integer, intent(in) :: nflds
real, intent(in) :: field(:,:,:)
! local variable
integer :: n, ncid, varid, rc, idimid, jdimid, kdimid, fdimid
real, allocatable :: a4d(:,:,:,:)
character(len=20) :: subname = 'dumpnc3d'
if (debug)print '(a)','enter '//trim(subname)//' variable '//vname
print *,dims
print *,size(field,1),size(field,2),size(field,3)
allocate(a4d(dims(1),dims(2),dims(3),nflds)); a4d = 0.0
rc = nf90_create(trim(fname), nf90_clobber, ncid)
rc = nf90_def_dim(ncid, 'nx', dims(1), idimid)
rc = nf90_def_dim(ncid, 'ny', dims(2), jdimid)
rc = nf90_def_dim(ncid, 'nk', dims(3), kdimid)
rc = nf90_def_dim(ncid, 'nf', nflds, fdimid)
rc = nf90_def_var(ncid, vname, nf90_float, (/idimid,jdimid,kdimid,fdimid/), varid)
rc = nf90_enddef(ncid)
do n = 1,nflds
a4d(:,:,:,n) = reshape(field(1:dims(1)*dims(2),1:dims(3),n), (/dims(1),dims(2),dims(3)/))
end do
rc = nf90_put_var(ncid, varid, a4d)
rc = nf90_close(ncid)
if (debug)print '(a)','exit '//trim(subname)//' variable '//vname
end subroutine dumpnc3d
!----------------------------------------------------------
! handle netcdf errors
!----------------------------------------------------------
subroutine handle_err(ierr,string)
integer , intent(in) :: ierr
character(len=*), intent(in) :: string
if (ierr /= nf90_noerr) then
write(ndse,*) "*** ERROR ***: ",trim(string),':',trim(nf90_strerror(ierr))
stop
end if
end subroutine handle_err
end module utils_mod