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