forked from NCAR/ccpp-physics
-
Notifications
You must be signed in to change notification settings - Fork 34
/
aerinterp.F90
561 lines (510 loc) · 17.2 KB
/
aerinterp.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
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
!>\file aerinterp.F90
!! This file contains subroutines of reading and interpolating
!! aerosol data for MG microphysics.
!>\ingroup mod_GFS_phys_time_vary
!! This module contain subroutines of reading and interpolating
!! aerosol data for MG microphysics.
module aerinterp
implicit none
private read_netfaer
public :: read_aerdata, setindxaer, aerinterpol,read_aerdataf
contains
logical function netcdf_check(status, errmsg, errflg, why)
use netcdf
implicit none
character(len=*), intent(inout) :: errmsg
integer, intent(out) :: errflg
integer, intent(in) :: status
character(len=*), intent(in) :: why
netcdf_check = (status == NF90_NOERR)
if(netcdf_check) then
errflg = 0
errmsg = ' '
else
errflg = 1
errmsg = trim(why) // ': ' // trim(nf90_strerror(status))
endif
END function netcdf_check
SUBROUTINE read_aerdata (me, master, iflip, idate, errmsg, errflg)
use machine, only: kind_phys, kind_io4, kind_io8
use aerclm_def
use netcdf
!--- in/out
integer, intent(in) :: me, master, iflip, idate(4)
character(len=*), intent(inout) :: errmsg
integer, intent(inout) :: errflg
!--- locals
integer :: ncid, varid, ndims, hmx
integer :: i, j, k, n, ii, imon, klev
character :: fname*50, mn*2, vname*10
logical :: file_exist
integer :: dimids(NF90_MAX_VAR_DIMS)
integer :: dimlen(NF90_MAX_VAR_DIMS)
errflg = 0
errmsg = ' '
!
!! ===================================================================
if (me == master) then
if ( iflip == 0 ) then ! data from toa to sfc
print *, "GFS is top-down"
else
print *, "GFS is bottom-up"
endif
endif
!
!! ===================================================================
!! check if all necessary files exist
!! ===================================================================
do imon = 1, 12
write(mn,'(i2.2)') imon
fname=trim("aeroclim.m"//mn//".nc")
inquire (file = fname, exist = file_exist)
if (.not. file_exist) then
errmsg = 'Error in read_aerdata: file ' // trim(fname) // ' not found'
errflg = 1
return
endif
enddo
!
!! ===================================================================
!! fetch dim spec and lat/lon from m01 data set
!! ===================================================================
fname=trim("aeroclim.m"//'01'//".nc")
ncid = -1
if(.not.netcdf_check(nf90_open(fname , nf90_NOWRITE, ncid), &
errmsg, errflg, 'open '//trim(fname))) then
return
endif
vname = trim(specname(1))
varid = -1
if(.not.netcdf_check(nf90_inq_varid(ncid, vname, varid), &
errmsg, errflg, 'find id of '//trim(vname)//' var')) then
return
endif
ndims = 0
if(.not.netcdf_check(nf90_inquire_variable(ncid, varid, ndims=ndims, dimids=dimids), &
errmsg, errflg, 'inquire details about '//trim(vname)//' var')) then
return
endif
do i=1,ndims
if(.not.netcdf_check(nf90_inquire_dimension(ncid, dimids(i), len=dimlen(i)), &
errmsg, errflg, 'inquire details about dimension')) then
return
endif
enddo
! specify latsaer, lonsaer, hmx
lonsaer = dimlen(1)
latsaer = dimlen(2)
levsw = dimlen(3)
if(me==master) then
print *, 'MERRA2 dim: ',dimlen(1:ndims)
endif
! allocate arrays
if (.not. allocated(aer_lat)) then
allocate(aer_lat(latsaer))
allocate(aer_lon(lonsaer))
endif
! construct lat/lon array
varid = -1
if(.not.netcdf_check(nf90_inq_varid(ncid, 'lat', varid), &
errmsg, errflg, 'find id of lat var')) then
return
endif
aer_lat = 0
if(.not.netcdf_check(nf90_get_var(ncid, varid, aer_lat, (/ 1, 1, 1 /), (/latsaer, 1, 1/)), &
errmsg, errflg, 'read lat var')) then
return
endif
varid = -1
if(.not.netcdf_check(nf90_inq_varid(ncid, 'lon', varid), &
errmsg, errflg, 'find id of lon var')) then
return
endif
aer_lon = 0
if(.not.netcdf_check(nf90_get_var(ncid, varid, aer_lon, (/ 1, 1, 1 /), (/lonsaer, 1, 1/)), &
errmsg, errflg, 'read lon var')) then
return
endif
if(.not.netcdf_check(nf90_close(ncid), errmsg, errflg, 'close '//trim(fname))) then
return
endif
END SUBROUTINE read_aerdata
!
!**********************************************************************
SUBROUTINE read_aerdataf ( me, master, iflip, idate, FHOUR, errmsg, errflg)
use machine, only: kind_phys, kind_io4, kind_io8
use aerclm_def
!--- in/out
integer, intent(in) :: me, master, iflip, idate(4)
character(len=*), intent(inout) :: errmsg
integer, intent(inout) :: errflg
real(kind=kind_phys), intent(in) :: fhour
!--- locals
integer :: i, j, k, n, ii, imon, klev, n1, n2
logical :: file_exist
integer IDAT(8),JDAT(8)
real(kind=kind_phys) rjday
real(8) RINC(5)
integer jdow, jdoy, jday
real(4) rinc4(5)
integer w3kindreal,w3kindint
integer, allocatable :: invardims(:)
!
if (.not. allocated(aerin)) then
allocate(aerin(iamin:iamax,jamin:jamax,levsaer,ntrcaerm,timeaer))
allocate(aer_pres(iamin:iamax,jamin:jamax,levsaer,timeaer))
endif
! allocate local working arrays
!! found interpolation months
IDAT = 0
IDAT(1) = IDATE(4)
IDAT(2) = IDATE(2)
IDAT(3) = IDATE(3)
IDAT(5) = IDATE(1)
RINC = 0.
RINC(2) = FHOUR
call w3kind(w3kindreal,w3kindint)
if(w3kindreal == 4) then
rinc4 = rinc
CALL W3MOVDAT(RINC4,IDAT,JDAT)
else
CALL W3MOVDAT(RINC,IDAT,JDAT)
endif
!
jdow = 0
jdoy = 0
jday = 0
call w3doxdat(jdat,jdow,jdoy,jday)
rjday = jdoy + jdat(5) / 24.
IF (RJDAY < aer_time(1)) RJDAY = RJDAY+365.
!
n2 = 13
do j=2, 12
if (rjday < aer_time(j)) then
n2 = j
exit
endif
enddo
n1 = n2 - 1
if (n2 > 12) n2 = n2 -12
!! ===================================================================
call read_netfaer(n1, iflip, 1, errmsg, errflg)
if(errflg/=0) return
call read_netfaer(n2, iflip, 2, errmsg, errflg)
if(errflg/=0) return
!! ===================================================================
n1sv=n1
n2sv=n2
!---
END SUBROUTINE read_aerdataf
!
SUBROUTINE setindxaer(npts,dlat,jindx1,jindx2,ddy,dlon, &
iindx1,iindx2,ddx,me,master)
!
USE MACHINE, ONLY: kind_phys
use aerclm_def, only: aer_lat, jaero=>latsaer, &
aer_lon, iaero=>lonsaer
!
implicit none
!
integer me, master
integer npts, JINDX1(npts),JINDX2(npts),IINDX1(npts),IINDX2(npts)
real(kind=kind_phys) dlat(npts),DDY(npts),dlon(npts),DDX(npts)
!
integer i,j
DO J=1,npts
jindx2(j) = jaero + 1
do i=1,jaero
if (dlat(j) < aer_lat(i)) then
jindx2(j) = i
exit
endif
enddo
jindx1(j) = max(jindx2(j)-1,1)
jindx2(j) = min(jindx2(j),jaero)
if (jindx2(j) .ne. jindx1(j)) then
DDY(j) = (dlat(j) - aer_lat(jindx1(j))) &
/ (aer_lat(jindx2(j)) - aer_lat(jindx1(j)))
else
ddy(j) = 1.0
endif
ENDDO
DO J=1,npts
iindx2(j) = iaero + 1
do i=1,iaero
if (dlon(j) < aer_lon(i)) then
iindx2(j) = i
exit
endif
enddo
iindx1(j) = max(iindx2(j)-1,1)
iindx2(j) = min(iindx2(j),iaero)
if (iindx2(j) .ne. iindx1(j)) then
ddx(j) = (dlon(j) - aer_lon(iindx1(j))) &
/ (aer_lon(iindx2(j)) - aer_lon(iindx1(j)))
else
ddx(j) = 1.0
endif
ENDDO
RETURN
END SUBROUTINE setindxaer
!
!**********************************************************************
!**********************************************************************
!
SUBROUTINE aerinterpol( me,master,nthrds,npts,IDATE,FHOUR,iflip, jindx1,jindx2, &
ddy,iindx1,iindx2,ddx,lev,prsl,aerout, errmsg,errflg)
!
use machine, only: kind_phys, kind_io4, kind_io8
use aerclm_def
implicit none
integer, intent(inout) :: errflg
character(*), intent(inout) :: errmsg
integer, intent(in) :: iflip
integer i1,i2, iday,j,j1,j2,l,npts,nc,n1,n2,lev,k,i,ii, klev
real(kind=kind_phys) fhour,temj, tx1, tx2,temi, tem
real(kind=kind_phys), dimension(npts) :: temij,temiy,temjx,ddxy
!
integer JINDX1(npts), JINDX2(npts), iINDX1(npts), iINDX2(npts)
integer me,idate(4), master, nthrds
integer IDAT(8),JDAT(8)
!
real(kind=kind_phys) DDY(npts), ddx(npts),ttt
real(kind=kind_phys) aerout(npts,lev,ntrcaer)
real(kind=kind_phys) aerpm(npts,levsaer,ntrcaer)
real(kind=kind_phys) prsl(npts,lev), aerpres(npts,levsaer)
real(kind=kind_phys) rjday
integer jdow, jdoy, jday
real(8) RINC(5)
real(4) rinc4(5)
integer w3kindreal,w3kindint
!
errflg = 0
errmsg = ' '
IDAT = 0
IDAT(1) = IDATE(4)
IDAT(2) = IDATE(2)
IDAT(3) = IDATE(3)
IDAT(5) = IDATE(1)
RINC = 0.
RINC(2) = FHOUR
call w3kind(w3kindreal,w3kindint)
if(w3kindreal == 4) then
rinc4 = rinc
CALL W3MOVDAT(RINC4,IDAT,JDAT)
else
CALL W3MOVDAT(RINC,IDAT,JDAT)
endif
!
jdow = 0
jdoy = 0
jday = 0
call w3doxdat(jdat,jdow,jdoy,jday)
rjday = jdoy + jdat(5) / 24.
IF (RJDAY < aer_time(1)) RJDAY = RJDAY+365.
!
n2 = 13
do j=2, 12
if (rjday < aer_time(j)) then
n2 = j
exit
endif
enddo
n1 = n2 - 1
if (n2 > 12) n2 = n2 -12
! need to read a new month
if (n1.ne.n1sv) then
#ifdef DEBUG
if (me == master) write(*,*)"read in a new month MERRA2", n2
#endif
DO ii = 1, ntrcaerm
do j = jamin, jamax
do k = 1, levsaer
do i = iamin, iamax
aerin(i,j,k,ii,1) = aerin(i,j,k,ii,2)
enddo !i-loop (lon)
enddo !k-loop (lev)
enddo !j-loop (lat)
ENDDO ! ii-loop (ntracaerm)
!! ===================================================================
call read_netfaer(n2, iflip, 2, errmsg, errflg)
n1sv=n1
n2sv=n2
end if
!
tx1 = (aer_time(n2) - rjday) / (aer_time(n2) - aer_time(n1))
tx2 = 1.0 - tx1
if (n2 > 12) n2 = n2 -12
do j=1,npts
TEMJ = 1.0 - DDY(J)
TEMI = 1.0 - DDX(J)
temij(j) = TEMI*TEMJ
temiy(j) = TEMI*DDY(j)
temjx(j) = TEMJ*DDX(j)
ddxy(j) = DDX(j)*DDY(J)
enddo
#ifndef __GFORTRAN__
!$OMP parallel num_threads(nthrds) default(none) &
!$OMP shared(npts,ntrcaer,aerin,aer_pres,prsl) &
!$OMP shared(ddx,ddy,jindx1,jindx2,iindx1,iindx2) &
!$OMP shared(aerpm,aerpres,aerout,lev,nthrds) &
!$OMP shared(temij,temiy,temjx,ddxy) &
!$OMP private(l,j,k,ii,i1,i2,j1,j2,tem) &
!$OMP copyin(tx1,tx2) firstprivate(tx1,tx2)
!$OMP do
#endif
DO L=1,levsaer
DO J=1,npts
J1 = JINDX1(J)
J2 = JINDX2(J)
I1 = IINDX1(J)
I2 = IINDX2(J)
DO ii=1,ntrcaer
aerpm(j,L,ii) = &
tx1*(TEMIJ(j)*aerin(I1,J1,L,ii,1)+DDXY(j)*aerin(I2,J2,L,ii,1) &
+TEMIY(j)*aerin(I1,J2,L,ii,1)+temjx(j)*aerin(I2,J1,L,ii,1))&
+tx2*(TEMIJ(j)*aerin(I1,J1,L,ii,2)+DDXY(j)*aerin(I2,J2,L,ii,2) &
+TEMIY(j)*aerin(I1,J2,L,ii,2)+temjx(j)*aerin(I2,J1,L,ii,2))
ENDDO
aerpres(j,L) = &
tx1*(TEMIJ(j)*aer_pres(I1,J1,L,1)+DDXY(j)*aer_pres(I2,J2,L,1) &
+TEMIY(j)*aer_pres(I1,J2,L,1)+temjx(j)*aer_pres(I2,J1,L,1))&
+tx2*(TEMIJ(j)*aer_pres(I1,J1,L,2)+DDXY(j)*aer_pres(I2,J2,L,2) &
+TEMIY(j)*aer_pres(I1,J2,L,2)+temjx(j)*aer_pres(I2,J1,L,2))
ENDDO
ENDDO
#ifndef __GFORTRAN__
!$OMP end do
! don't flip, input is the same direction as GFS (bottom-up)
!$OMP do
#endif
DO J=1,npts
DO L=1,lev
if(prsl(j,L) >= aerpres(j,1)) then
DO ii=1, ntrcaer
aerout(j,L,ii) = aerpm(j,1,ii) !! sfc level
ENDDO
else if(prsl(j,L) <= aerpres(j,levsaer)) then
DO ii=1, ntrcaer
aerout(j,L,ii) = aerpm(j,levsaer,ii) !! toa top
ENDDO
else
DO k=1, levsaer-1 !! from sfc to toa
IF(prsl(j,L) <= aerpres(j,k) .and. prsl(j,L)>aerpres(j,k+1)) then
i1 = k
i2 = min(k+1,levsaer)
exit
ENDIF
ENDDO
tem = 1.0 / (aerpres(j,i1) - aerpres(j,i2))
tx1 = (prsl(j,L) - aerpres(j,i2)) * tem
tx2 = (aerpres(j,i1) - prsl(j,L)) * tem
DO ii = 1, ntrcaer
aerout(j,L,ii) = aerpm(j,i1,ii)*tx1 + aerpm(j,i2,ii)*tx2
ENDDO
endif
ENDDO !L-loop
ENDDO !J-loop
#ifndef __GFORTRAN__
!$OMP end do
!$OMP end parallel
#endif
RETURN
END SUBROUTINE aerinterpol
subroutine read_netfaer(nf, iflip,nt, errmsg,errflg)
use machine, only: kind_phys, kind_io4, kind_io8
use aerclm_def
use netcdf
integer, intent(in) :: iflip, nf, nt
integer, intent(inout) :: errflg
character(*), intent(inout) :: errmsg
integer :: ncid, varid, i,j,k,ii,klev
character :: fname*50, mn*2, vname*10
real(kind=kind_io4),allocatable,dimension(:,:,:) :: buff
real(kind=kind_io4),allocatable,dimension(:,:,:,:):: buffx
real(kind=kind_io4),allocatable,dimension(:,:) :: pres_tmp
!! ===================================================================
allocate (buff(lonsaer, latsaer, levsw))
allocate (pres_tmp(lonsaer, levsw))
allocate (buffx(lonsaer, latsaer, levsw, 1))
errflg = 0
errmsg = ' '
buff = 0
pres_tmp = 0
buffx = 0
write(mn,'(i2.2)') nf
fname=trim("aeroclim.m"//mn//".nc")
ncid = -1
if(.not.netcdf_check(nf90_open(fname , nf90_NOWRITE, ncid), &
errmsg, errflg, 'open '//trim(fname))) then
return
endif
! ====> construct 3-d pressure array (Pa)
varid = -1
if(.not.netcdf_check(nf90_inq_varid(ncid, "DELP", varid), &
errmsg, errflg, 'find id of DELP var')) then
return
endif
if(.not.netcdf_check(nf90_get_var(ncid, varid, buff), &
errmsg, errflg, 'read DELP var')) then
return
endif
do j = jamin, jamax
do i = iamin, iamax
! constract pres_tmp (top-down), note input is top-down
pres_tmp(i,1) = 0.
do k=2, levsw
pres_tmp(i,k) = pres_tmp(i,k-1)+buff(i,j,k)
enddo !k-loop
enddo !i-loop (lon)
! extract pres_tmp to fill aer_pres (in Pa)
do k = 1, levsaer
if ( iflip == 0 ) then ! data from toa to sfc
klev = k
else ! data from sfc to top
klev = ( levsw - k ) + 1
endif
do i = iamin, iamax
aer_pres(i,j,k,nt) = 1.d0*pres_tmp(i,klev)
enddo !i-loop (lon)
enddo !k-loop (lev)
enddo !j-loop (lat)
! ====> construct 4-d aerosol array (kg/kg)
! merra2 data is top down
! for GFS, iflip 0: toa to sfc; 1: sfc to toa
DO ii = 1, ntrcaerm
vname=trim(specname(ii))
varid = -1
if(.not.netcdf_check(nf90_inq_varid(ncid, vname, varid), &
errmsg, errflg, 'get id of '//trim(vname)//' var')) then
return
endif
if(.not.netcdf_check(nf90_get_var(ncid, varid, buffx), &
errmsg, errflg, 'read '//trim(vname)//' var')) then
return
endif
do j = jamin, jamax
do k = 1, levsaer
! input is from toa to sfc
if ( iflip == 0 ) then ! data from toa to sfc
klev = k
else ! data from sfc to top
klev = ( levsw - k ) + 1
endif
do i = iamin, iamax
aerin(i,j,k,ii,nt) = 1.d0*buffx(i,j,klev,1)
if(aerin(i,j,k,ii,nt) < 0 .or. aerin(i,j,k,ii,nt) > 1.) then
aerin(i,j,k,ii,nt) = 1.e-15
endif
enddo !i-loop (lon)
enddo !k-loop (lev)
enddo !j-loop (lat)
ENDDO ! ii-loop (ntracaerm)
! close the file
if(.not.netcdf_check(nf90_close(ncid), errmsg, errflg, 'close '//trim(fname))) then
return
endif
deallocate (buff, pres_tmp)
deallocate (buffx)
END SUBROUTINE read_netfaer
end module aerinterp