SCALE-RM
mod_realinput_nicam.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
10 #include "scalelib.h"
12  !-----------------------------------------------------------------------------
13  !
14  !++ used modules
15  !
16  use scale_precision
17  use scale_io
18  use scale_tracer
19  use scale_prc, only: &
20  myrank => prc_myrank, &
21  prc_abort
22 
23  !-----------------------------------------------------------------------------
24  implicit none
25  private
26  !-----------------------------------------------------------------------------
27  !
28  !++ Public procedure
29  !
30  public :: parentatmossetupnicam
31  public :: parentatmosopennicam
32  public :: parentatmosinputnicam
33  public :: parentlandsetupnicam
34  public :: parentlandinputnicam
35  public :: parentoceansetupnicam
36  public :: parentoceanopennicam
37  public :: parentoceaninputnicam
38 
39  !-----------------------------------------------------------------------------
40  !
41  !++ Public parameters & variables
42  !
43  !-----------------------------------------------------------------------------
44  !
45  !++ Private procedure
46  !
47  !-----------------------------------------------------------------------------
48  !
49  !++ Private parameters & variables
50  !
51  real(RP), parameter :: missval_tg = 50.0_rp
52  real(RP), parameter :: missval_strg = 0.0_rp
53  real(RP), parameter :: missval_sst = 273.1506_rp
54 
55 
56  real(SP), allocatable :: read1dx(:)
57  real(SP), allocatable :: read1dy(:)
58  real(SP), allocatable :: read1dz(:)
59  real(SP), allocatable :: read3d (:,:,:)
60 
61  real(SP), allocatable :: read3dl(:,:,:)
62 
63  !-----------------------------------------------------------------------------
64 contains
65  !-----------------------------------------------------------------------------
67  subroutine parentatmossetupnicam( &
68  dims, &
69  timelen, &
70  basename_org )
71  use scale_file, only: &
72  file_get_shape
73  implicit none
74 
75  integer, intent(out) :: dims(6)
76  integer, intent(out) :: timelen
77  character(len=*), intent(in) :: basename_org
78 
79  character(len=H_LONG) :: basename
80  integer :: dims_ncm(4)
81 
82  !---------------------------------------------------------------------------
83 
84  log_info("ParentAtmosSetupNICAM",*) 'Real Case/Atmos Input File Type: NICAM-NETCDF'
85  basename = "ms_pres"//trim(basename_org)
86  call file_get_shape( basename, "ms_pres", & ! (in)
87  dims_ncm(:), & ! (out)
88  single=.true. ) ! (in)
89  timelen = dims_ncm(4)
90 
91  ! full level
92  dims(1) = dims_ncm(3)
93  dims(2) = dims_ncm(1)
94  dims(3) = dims_ncm(2)
95  ! half level
96  ! nicam lat-lon data doesn't have staggered grid system
97  dims(4) = dims(1)
98  dims(5) = dims(2)
99  dims(6) = dims(3)
100 
101  allocate( read1dx( dims(2) ) )
102  allocate( read1dy( dims(3) ) )
103  allocate( read1dz( dims(1) ) )
104  allocate( read3d( dims(1), dims(2), dims(3) ) )
105 
106  return
107  end subroutine parentatmossetupnicam
108 
109  !-----------------------------------------------------------------------------
110  subroutine parentatmosopennicam( &
111  lon_org, &
112  lat_org, &
113  cz_org, &
114  basename_num, &
115  dims )
116  use scale_const, only: &
117  d2r => const_d2r
118  use scale_file, only: &
119  file_open, &
120  file_read
121  implicit none
122 
123  real(RP), intent(out) :: lon_org(:,:)
124  real(RP), intent(out) :: lat_org(:,:)
125  real(RP), intent(out) :: cz_org(:,:,:)
126  character(len=*), intent(in) :: basename_num
127  integer, intent(in) :: dims(6)
128 
129  character(len=H_LONG) :: basename
130 
131  integer :: fid
132  integer :: k, i, j
133 
134  basename = "ms_pres"//trim(basename_num)
135  call file_open( basename, fid, single=.true., postfix="" )
136 
137  call file_read( fid, "lon", read1dx(:) )
138  do j = 1, dims(3)
139  lon_org(:,j) = read1dx(:) * d2r
140  enddo
141 
142  call file_read( fid, "lat", read1dy(:) )
143  do i = 1, dims(2)
144  lat_org(i,:) = read1dy(:) * d2r
145  enddo
146 
147  call file_read( fid, "lev", read1dz(:) )
148  do j = 1, dims(3)
149  do i = 1, dims(2)
150  cz_org(3:,i,j) = read1dz(:)
151  end do
152  end do
153 
154  call file_read( fid, "ms_pres", read3d(:,:,:) )
155  do j = 1, dims(3)
156  do i = 1, dims(2)
157  do k = dims(1)+2, 3, -1
158  ! missing data implies under ground
159  ! For more accurate surface height, we need topograph data
160  ! So surface data is not used at this moment
161  if ( read3d(i,j,k-2) <= 0.0_rp ) then
162 ! cz_org(k,i,j) = max( ( cz_org(k+1,i,j) + cz_org(k,i,j) ) * 0.5_RP, 0.0_RP )
163  cz_org(1:k,i,j) = 0.0_rp
164  exit
165  else
166  cz_org(k,i,j) = read3d(i,j,k)
167  endif
168  enddo
169  enddo
170  enddo
171 
172  return
173  end subroutine parentatmosopennicam
174 
175  !-----------------------------------------------------------------------------
176  subroutine parentatmosinputnicam( &
177  velz_org, &
178  velx_org, &
179  vely_org, &
180  pres_org, &
181  temp_org, &
182  qtrc_org, &
183  basename_num, &
184  dims, &
185  it )
186  use scale_const, only: &
187  rdry => const_rdry, &
188  cpdry => const_cpdry, &
189  p00 => const_pre00, &
190  d2r => const_d2r, &
191  eps => const_eps
192  use scale_file, only: &
193  file_read
194  use scale_atmos_hydrometeor, only: &
195  i_qv
196  implicit none
197 
198  real(RP), intent(out) :: velz_org(:,:,:)
199  real(RP), intent(out) :: velx_org(:,:,:)
200  real(RP), intent(out) :: vely_org(:,:,:)
201  real(RP), intent(out) :: pres_org(:,:,:)
202  real(RP), intent(out) :: temp_org(:,:,:)
203  real(RP), intent(out) :: qtrc_org(:,:,:,:)
204  character(len=*), intent(in) :: basename_num
205  integer, intent(in) :: dims(6)
206  integer, intent(in) :: it
207 
208  real(RP) :: tsfc_org (dims(2),dims(3))
209  real(RP) :: slp_org (dims(2),dims(3))
210  real(RP) :: qvsfc_org(dims(2),dims(3))
211 
212  real(RP) :: pott
213  real(RP) :: RovCP
214  real(RP) :: CPovR
215 
216  integer :: k, i, j
217 
218  character(len=H_LONG) :: basename
219  !---------------------------------------------------------------------------
220 
221  basename = "ms_u"//trim(basename_num)
222  call file_read( basename, "ms_u", read3d(:,:,:), step=it, rankid=myrank, single=.true., postfix="" )
223  do j = 1, dims(3)
224  do i = 1, dims(2)
225  do k = 1, dims(1)
226  if ( abs(abs(read3d(i,j,k))-300.0_sp) < sqrt(eps) ) then ! missing value
227  velx_org(k+2,i,j) = 0.0_rp
228  else
229  velx_org(k+2,i,j) = read3d(i,j,k)
230  end if
231  end do
232  velx_org(1:2,i,j) = 0.0_rp
233  end do
234  end do
235 
236  basename = "ms_v"//trim(basename_num)
237  call file_read( basename, "ms_v", read3d(:,:,:), step=it, rankid=myrank, single=.true., postfix="" )
238  do j = 1, dims(3)
239  do i = 1, dims(2)
240  do k = 1, dims(1)
241  if ( abs(abs(read3d(i,j,k))-300.0_sp) < sqrt(eps) ) then ! missing value
242  vely_org(k+2,i,j) = 0.0_rp
243  else
244  vely_org(k+2,i,j) = read3d(i,j,k)
245  end if
246  end do
247  vely_org(1:2,i,j) = 0.0_rp
248  end do
249  end do
250 
251  velz_org(:,:,:) = 0.0_rp
252 
253 
254  ! The surface height is not available, so t2 is not used, at this moment.
255 ! basename = "ss_t2m"//trim(basename_num)
256 ! call FILE_read( basename, "ss_t2m", tsfc_org(:,:), step=it, rankid=myrank, single=.true., postfix="" )
257 
258  basename = "ms_tem"//trim(basename_num)
259  call file_read( basename, "ms_tem", read3d(:,:,:), step=it, rankid=myrank, single=.true., postfix="" )
260  do j = 1, dims(3)
261  do i = 1, dims(2)
262  do k = dims(1), 1, -1
263  if ( read3d(i,j,k) <= 50.0_rp ) then ! missing value
264  temp_org(1:k+2,i,j) = temp_org(k+3,i,j)
265  exit
266  else
267  temp_org(k+2,i,j) = read3d(i,j,k)
268  end if
269 ! if( k==1 ) temp_org(2,i,j) = tsfc_org(i,j) ! no missing value case
270  end do
271  end do
272  end do
273 
274  rovcp = rdry/cpdry
275  cpovr = cpdry/rdry
276 
277  basename = "ss_slp"//trim(basename_num)
278  call file_read( basename, "ss_slp", slp_org(:,:), step=it, rankid=myrank, single=.true., postfix="" )
279 
280  basename = "ms_pres"//trim(basename_num)
281  call file_read( basename, "ms_pres", read3d(:,:,:), step=it, rankid=myrank, single=.true., postfix="" )
282  do j = 1, dims(3)
283  do i = 1, dims(2)
284  do k = dims(1), 1, -1
285  if ( read3d(i,j,k) < 0.0_rp ) then ! missing value
286  exit
287  else
288  pres_org(k+2,i,j) = read3d(i,j,k)
289  end if
290  end do
291  ! If data has missing value, k is the top layer having missing value,
292  ! otherwise k is zero.
293  pott = temp_org(k+3,i,j) * (p00/pres_org(k+3,i,j))**rovcp ! lowest level
294 ! pres_org(k+2,i,j) = P00 * (temp_org(k+2,i,j)/pott)**CPovR ! surface
295  pres_org(1:k+2,i,j) = slp_org(i,j) ! sea level
296  temp_org(1:k+2,i,j) = pott * (slp_org(i,j)/p00)**rovcp ! sea level
297  end do
298  end do
299 
300  basename = "ss_q2m"//trim(basename_num)
301  call file_read( basename, "ss_q2m", qvsfc_org(:,:), step=it, rankid=myrank, single=.true., postfix="" )
302 
303  qtrc_org(:,:,:,:) = 0.0_rp
304  basename = "ms_qv"//trim(basename_num)
305  call file_read( basename, "ms_qv", read3d(:,:,:), step=it, rankid=myrank, single=.true., postfix="" )
306  do j = 1, dims(3)
307  do i = 1, dims(2)
308  do k = dims(1), 1, -1
309  if ( read3d(i,j,k) < 0.0_rp ) then ! missing value
310  exit
311  else
312  qtrc_org(k+2,i,j,i_qv) = read3d(i,j,k)
313  end if
314  end do
315  qtrc_org(1:k+2,i,j,i_qv) = qvsfc_org(i,j)
316  end do
317  end do
318 
319 
320  return
321  end subroutine parentatmosinputnicam
322 
323  !-----------------------------------------------------------------------------
325  subroutine parentlandsetupnicam( &
326  ldims, &
327  basename_org )
328  use scale_file, only: &
329  file_get_shape
330  implicit none
331 
332  integer, intent(out) :: ldims(3)
333  character(len=*), intent(in) :: basename_org
334 
335  character(len=H_LONG) :: basename
336  integer :: dims_ncm(4)
337 
338  !---------------------------------------------------------------------------
339 
340  log_info("ParentLandSetupNICAM",*) 'Real Case/Land Input File Type: NICAM-NETCDF'
341  basename = "la_tg"//trim(basename_org)
342  call file_get_shape( basename, "la_tg", & ! (in)
343  dims_ncm(:), & ! (out)
344  single=.true. ) ! (in)
345  ! land
346  ldims(1) = dims_ncm(3) ! vertical grid of land model
347  ldims(2) = dims_ncm(1)
348  ldims(3) = dims_ncm(2)
349 
350  if ( .not. allocated(read1dx) ) then
351  allocate( read1dx( ldims(2) ) )
352  allocate( read1dy( ldims(3) ) )
353  end if
354 
355  allocate( read3dl( ldims(2), ldims(3), ldims(1) ) )
356 
357  return
358  end subroutine parentlandsetupnicam
359 
360  !-----------------------------------------------------------------------------
361  subroutine parentlandinputnicam( &
362  tg_org, &
363  strg_org, &
364  lst_org, &
365  llon_org, &
366  llat_org, &
367  lz_org, &
368  topo_org, &
369  lmask_org, &
370  basename_num, &
371  ldims, &
372  use_file_landwater, &
373  it )
374  use scale_const, only: &
375  undef => const_undef, &
376  d2r => const_d2r, &
377  tem00 => const_tem00, &
378  eps => const_eps
379  use scale_file, only: &
380  file_read
381  implicit none
382 
383  real(RP), intent(out) :: tg_org(:,:,:)
384  real(RP), intent(out) :: strg_org(:,:,:)
385  real(RP), intent(out) :: lst_org(:,:)
386  real(RP), intent(out) :: llon_org(:,:)
387  real(RP), intent(out) :: llat_org(:,:)
388  real(RP), intent(out) :: lz_org(:)
389  real(RP), intent(out) :: topo_org(:,:)
390  real(RP), intent(out) :: lmask_org(:,:)
391  character(len=*), intent(in) :: basename_num
392  integer, intent(in) :: ldims(3)
393  logical, intent(in) :: use_file_landwater ! use land water data from files
394  integer, intent(in) :: it
395 
396  ! work
397  integer :: k, i, j
398 
399  character(len=H_LONG) :: basename
400  !---------------------------------------------------------------------------
401 
402  basename = "la_tg"//trim(basename_num)
403  call file_read( basename, "lev", lz_org(:), single=.true., postfix="" )
404 
405  call file_read( basename, "lon", read1dx(:), single=.true., postfix="" )
406  do j = 1, ldims(3)
407  llon_org(:,j) = read1dx(:) * d2r
408  enddo
409 
410  call file_read( basename, "lat", read1dy(:), single=.true., postfix="" )
411  do i = 1, ldims(2)
412  llat_org(i,:) = read1dy(:) * d2r
413  enddo
414 
415  basename = "lsmask"//trim(basename_num)
416  call file_read( basename, "lsmask", lmask_org(:,:), step=it, rankid=myrank, single=.true., allow_missing=.true., missing_value=undef, postfix="" )
417 
418  basename = "la_tg"//trim(basename_num)
419  call file_read( basename, "la_tg", read3dl(:,:,:), step=it, rankid=myrank, single=.true., postfix="" )
420  do j = 1, ldims(3)
421  do i = 1, ldims(2)
422  do k = 1, ldims(1)
423  tg_org(k,i,j) = read3dl(i,j,k)
424  end do
425  end do
426  end do
427 
428  if( use_file_landwater ) then
429  basename = "la_wg"//trim(basename_num)
430  call file_read( basename, "la_wg", read3dl(:,:,:), step=it, rankid=myrank, single=.true., postfix="" )
431  do j = 1, ldims(3)
432  do i = 1, ldims(2)
433  do k = 1, ldims(1)
434  strg_org(k,i,j) = read3dl(i,j,k)
435  end do
436  end do
437  end do
438  end if
439 
440  basename = "ss_tem_sfc"//trim(basename_num)
441  call file_read( basename, "ss_tem_sfc", lst_org(:,:), step=it, rankid=myrank, single=.true., postfix="" )
442 
443  ! replace missing value
444  do j = 1, ldims(3)
445  do i = 1, ldims(2)
446  do k = 1, ldims(1)
447  if ( abs(tg_org(k,i,j) -missval_tg ) < eps ) tg_org(k,i,j) = undef
448  if ( abs(strg_org(k,i,j)-missval_strg) < eps ) strg_org(k,i,j) = undef
449  end do
450  end do
451  end do
452 
453  topo_org = undef
454 
455  return
456  end subroutine parentlandinputnicam
457 
458  !-----------------------------------------------------------------------------
460  subroutine parentoceansetupnicam( &
461  odims, &
462  timelen, &
463  basename_org )
464  use scale_file, only: &
465  file_get_shape
466  implicit none
467 
468  integer, intent(out) :: odims(2)
469  integer, intent(out) :: timelen
470  character(len=*), intent(in) :: basename_org
471 
472  character(len=H_LONG) :: basename
473  integer :: dims_ncm(4)
474 
475  !---------------------------------------------------------------------------
476 
477  log_info("ParentOceanSetupNICAM",*) 'Real Case/Ocean Input File Type: NICAM-NETCDF'
478 
479  basename = "oa_sst"//trim(basename_org)
480  call file_get_shape( basename, "oa_sst", &
481  dims_ncm(:), &
482  single=.true. )
483  odims(1) = dims_ncm(1)
484  odims(2) = dims_ncm(2)
485 
486  timelen = dims_ncm(4)
487 
488  if ( .not. allocated(read1dx) ) then
489  allocate( read1dx( odims(1) ) )
490  allocate( read1dy( odims(2) ) )
491  end if
492 
493  return
494  end subroutine parentoceansetupnicam
495 
496  !-----------------------------------------------------------------------------
497  subroutine parentoceanopennicam( &
498  olon_org, &
499  olat_org, &
500  omask_org, &
501  basename_num, &
502  odims )
503  use scale_const, only: &
504  d2r => const_d2r, &
505  undef => const_undef
506  use scale_file, only: &
507  file_read
508  implicit none
509 
510  real(RP), intent(out) :: olon_org(:,:)
511  real(RP), intent(out) :: olat_org(:,:)
512  real(RP), intent(out) :: omask_org(:,:)
513  character(len=*), intent(in) :: basename_num
514  integer, intent(in) :: odims(2)
515 
516  character(len=H_LONG) :: basename
517 
518  integer :: k, i, j
519 
520  basename = "oa_sst"//trim(basename_num)
521  call file_read( basename, "lon", read1dx(:), single=.true., postfix="" )
522  do j = 1, odims(2)
523  olon_org(:,j) = read1dx(:) * d2r
524  enddo
525 
526  call file_read( basename, "lat", read1dy(:), single=.true., postfix="" )
527  do i = 1, odims(1)
528  olat_org(i,:) = read1dy(:) * d2r
529  enddo
530 
531  basename = "lsmask"//trim(basename_num)
532  call file_read( basename, "lsmask", omask_org(:,:), step=1, rankid=myrank, single=.true., allow_missing=.true., missing_value=undef, postfix="" )
533 
534  return
535  end subroutine parentoceanopennicam
536 
537  !-----------------------------------------------------------------------------
538  subroutine parentoceaninputnicam( &
539  tw_org, &
540  sst_org, &
541  basename_num, &
542  odims, &
543  omask_org, &
544  it )
545  use scale_const, only: &
546  undef => const_undef, &
547  d2r => const_d2r, &
548  tem00 => const_tem00, &
549  eps => const_eps
550  use scale_file, only: &
551  file_read
552  implicit none
553 
554  real(RP), intent(out) :: tw_org(:,:)
555  real(RP), intent(out) :: sst_org(:,:)
556  character(len=*), intent(in) :: basename_num
557  integer, intent(in) :: odims(2)
558  real(RP), intent(in) :: omask_org(:,:)
559  integer, intent(in) :: it
560 
561  ! work
562  integer :: i, j
563 
564  character(len=H_LONG) :: basename
565  !---------------------------------------------------------------------------
566 
567  basename = "oa_sst"//trim(basename_num)
568  call file_read( basename, "oa_sst", sst_org(:,:), step=1, rankid=myrank, single=.true., postfix="" )
569 
570  !basename = "oa_ice"//trim(basename_num)
571  !call FILE_read( basename, "oa_ice", ice_org(:,:), step=1, rankid=myrank, single=.true., postfix="" )
572 
573 
574  ! SST: retrieve SST data around coast
575  do j = 1, odims(2)
576  do i = 1, odims(1)
577  if ( abs(omask_org(i,j)-1.0_rp) < eps ) then ! land data
578  cycle
579  else
580  sst_org(i,j) = ( sst_org(i,j) - missval_sst*omask_org(i,j) ) &
581  / ( 1.0_rp - omask_org(i,j) )
582  end if
583  end do
584  end do
585 
586  tw_org = sst_org
587 
588  return
589  end subroutine parentoceaninputnicam
590 
591 end module mod_realinput_nicam
real(rp), public const_cpdry
specific heat (dry air,constant pressure) [J/kg/K]
Definition: scale_const.F90:56
subroutine, public parentatmosinputnicam(velz_org, velx_org, vely_org, pres_org, temp_org, qtrc_org, basename_num, dims, it)
subroutine, public parentoceanopennicam(olon_org, olat_org, omask_org, basename_num, odims)
real(rp), parameter, public const_tem00
temperature reference (0C) [K]
Definition: scale_const.F90:90
real(rp), public const_d2r
degree to radian
Definition: scale_const.F90:32
subroutine, public file_open(basename, fid, mode, single, aggregate, rankid, postfix)
Definition: scale_file.F90:464
real(rp), public const_rdry
specific gas constant (dry air) [J/kg/K]
Definition: scale_const.F90:55
real(rp), public const_undef
Definition: scale_const.F90:41
module file
Definition: scale_file.F90:15
subroutine, public parentoceansetupnicam(odims, timelen, basename_org)
Ocean Setup.
module TRACER
module REAL input NICAM
module atmosphere / hydrometeor
subroutine, public parentlandinputnicam(tg_org, strg_org, lst_org, llon_org, llat_org, lz_org, topo_org, lmask_org, basename_num, ldims, use_file_landwater, it)
real(rp), public const_pre00
pressure reference [Pa]
Definition: scale_const.F90:88
subroutine, public parentlandsetupnicam(ldims, basename_org)
Land Setup.
subroutine, public parentatmossetupnicam(dims, timelen, basename_org)
Atmos Setup.
module PROCESS
Definition: scale_prc.F90:11
integer, public prc_myrank
process num in local communicator
Definition: scale_prc.F90:89
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:338
module CONSTANT
Definition: scale_const.F90:11
real(rp), public const_eps
small number
Definition: scale_const.F90:33
subroutine, public parentatmosopennicam(lon_org, lat_org, cz_org, basename_num, dims)
module PRECISION
module STDIO
Definition: scale_io.F90:10
subroutine, public parentoceaninputnicam(tw_org, sst_org, basename_num, odims, omask_org, it)