SCALE-RM
scale_file_history_cartesC.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
10 !-------------------------------------------------------------------------------
11 #include "scalelib.h"
13  !-----------------------------------------------------------------------------
14  !
15  !++ Used modules
16  !
17  use scale_precision
18  use scale_io
19  use scale_prof
24  use scale_prc, only: &
25  prc_abort
26  !-----------------------------------------------------------------------------
27  implicit none
28  private
29  !-----------------------------------------------------------------------------
30  !
31  !++ Public procedures
32  !
35 
36  !-----------------------------------------------------------------------------
37  !
38  !++ Public parameters & variables
39  !
40  !-----------------------------------------------------------------------------
41  !
42  !++ Private procedures
43  !
44  private :: file_history_cartesc_set_dims
45  private :: file_history_cartesc_set_axes
46  private :: file_history_cartesc_set_axes_attributes
47 
48  !-----------------------------------------------------------------------------
49  !
50  !++ Private parameters & variables
51  !
52  integer, parameter :: nzs = 3
53  character(len=8), parameter :: zs(nzs) = (/ "model ", &
54  "z ", &
55  "pressure" /)
56 
57  integer :: FILE_HISTORY_CARTESC_MODEL_nlayer = -1
58  integer :: FILE_HISTORY_CARTESC_PRES_nlayer = 0
59  real(RP), allocatable :: FILE_HISTORY_CARTESC_PRES_val(:)
60 
61  integer :: im, jm, km
62  integer :: ims, ime
63  integer :: jms, jme
64  integer :: imh, jmh
65  integer :: imsh, jmsh
66 
67  integer :: FILE_HISTORY_CARTESC_STARTDATE(6)
68  real(DP) :: FILE_HISTORY_CARTESC_STARTSUBSEC
69 
70  logical :: FILE_HISTORY_CARTESC_BOUNDARY = .false.
71 
72  logical :: pres_set = .false.
73 
74  !-----------------------------------------------------------------------------
75 contains
76  !-----------------------------------------------------------------------------
79  use scale_file_h, only: &
81  use scale_file_history, only: &
87  use scale_prc, only: &
90  use scale_prc_cartesc, only: &
93  prc_twod, &
94  prc_has_w, &
95  prc_has_s
96  use scale_time, only: &
97  time_nowdate, &
99  time_nowstep, &
100  time_dtsec, &
102  use scale_calendar, only: &
104  use scale_interp_vert, only: &
106  implicit none
107 
108  integer, parameter :: nlayer_max = 300
109  real(rp) :: file_history_cartesc_pres(nlayer_max)
110 
111  namelist / param_file_history_cartesc / &
112  file_history_cartesc_model_nlayer, &
113  file_history_cartesc_pres_nlayer, &
114  file_history_cartesc_pres, &
115  file_history_cartesc_boundary
116 
117  character(len=H_MID) :: file_history_cartesc_h_title = 'SCALE-RM FILE_HISTORY_CARTESC OUTPUT'
118  character(len=H_MID) :: file_history_cartesc_t_since
119 
120  character(len=FILE_HSHORT) :: calendar
121  real(dp) :: start_daysec
122  integer :: ierr
123  integer :: k
124  !---------------------------------------------------------------------------
125 
126  log_newline
127  log_info("FILE_HISTORY_CARTESC_setup",*) 'Setup'
128 
129  file_history_cartesc_pres(:) = 0.0_rp
130 
131  !--- read namelist
132  rewind(io_fid_conf)
133  read(io_fid_conf,nml=param_file_history_cartesc,iostat=ierr)
134  if( ierr < 0 ) then !--- missing
135  log_info("FILE_HISTORY_CARTESC_setup",*) 'Not found namelist. Default used.'
136  elseif( ierr > 0 ) then !--- fatal error
137  log_error("FILE_HISTORY_CARTESC_setup",*) 'Not appropriate names in namelist PARAM_FILE_HISTORY_CARTESC. Check!'
138  call prc_abort
139  endif
140  log_nml(param_file_history_cartesc)
141 
142  if ( file_history_cartesc_model_nlayer < 0 ) then
143  file_history_cartesc_model_nlayer = kmax
144  end if
145 
146 
147  ! check pressure coordinate
148  if ( file_history_cartesc_pres_nlayer > 0 ) then
149  if ( file_history_cartesc_pres_nlayer > nlayer_max ) then
150  log_error("FILE_HISTORY_CARTESC_setup",'(a,i3)') 'FILE_HISTORY_CARTESC_PRES_nlayer must be <= ', nlayer_max
151  call prc_abort
152  end if
153  allocate( file_history_cartesc_pres_val(file_history_cartesc_pres_nlayer) )
154 
155  do k = 1, file_history_cartesc_pres_nlayer
156  if ( file_history_cartesc_pres(k) <= 0.0_rp ) then
157  log_error("FILE_HISTORY_CARTESC_setup",'(a,i3,f7.1)') 'Invalid value found in pressure coordinate! (k,value)=', k, file_history_cartesc_pres(k)
158  call prc_abort
159  elseif ( file_history_cartesc_pres(k+1) >= file_history_cartesc_pres(k) ) then
160  log_error("FILE_HISTORY_CARTESC_setup",'(a,i3,2f7.1)') 'The value of pressure coordinate must be descending order! ', &
161  '(k,value[k],value[k+1])=', k, file_history_cartesc_pres(k), file_history_cartesc_pres(k+1)
162  call prc_abort
163  endif
164  file_history_cartesc_pres_val(k) = file_history_cartesc_pres(k) * 100.0_rp ! [hPa->Pa]
165  enddo
166 
167  call interp_vert_alloc_pres( file_history_cartesc_pres_nlayer, ka, ia, ja ) ! [IN]
168  else
169  log_info("FILE_HISTORY_CARTESC_setup",*) 'FILE_HISTORY_CARTESC_PRES_nlayer is not set.'
170  log_info("FILE_HISTORY_CARTESC_setup",*) 'Output with pressure coordinate is disabled'
171  endif
172 
173 
174 
175  file_history_cartesc_startdate(:) = time_nowdate
176  file_history_cartesc_startsubsec = time_nowsubsec
177 
178  start_daysec = time_startdaysec
179  if ( time_nowdate(1) > 0 ) then
180  write(file_history_cartesc_t_since,'(I4.4,5(A1,I2.2))') time_nowdate(1), &
181  '-', time_nowdate(2), &
182  '-', time_nowdate(3), &
183  ' ', time_nowdate(4), &
184  ':', time_nowdate(5), &
185  ':', time_nowdate(6)
186  start_daysec = time_nowsubsec
187  else
188  file_history_cartesc_t_since = ''
189  endif
190 
191  if ( file_history_cartesc_boundary ) then
192  ims = isb
193  ime = ieb
194  jms = jsb
195  jme = jeb
196 
197  imsh = ims
198  jmsh = jms
199 
200  im = imaxb
201  jm = jmaxb
202  imh = im
203  jmh = jm
204  else
205  ims = is
206  ime = ie
207  jms = js
208  jme = je
209 
210  if ( prc_twod .OR. prc_has_w .OR. prc_periodic_x ) then
211  imsh = ims
212  else
213  imsh = ims - 1 ! including i = IS-1
214  endif
215  if ( prc_has_s .OR. prc_periodic_y ) then
216  jmsh = jms
217  else
218  jmsh = jms - 1 ! include j = JS-1
219  endif
220 
221  im = ime - ims + 1
222  jm = jme - jms + 1
223  imh = ime - imsh + 1
224  jmh = jme - jmsh + 1
225  endif
226 
227  ! get calendar name
228  call calendar_get_name( calendar )
229 
230  call file_history_setup( file_history_cartesc_h_title, & ! [IN]
231  h_source, h_institute, & ! [IN]
232  start_daysec, time_dtsec, & ! [IN]
233  time_since = file_history_cartesc_t_since, & ! [IN]
234  calendar = calendar, & ! [IN]
235  default_zcoord = 'model', & ! [IN]
236  myrank = prc_myrank ) ! [IN]
237 
239 
240  call file_history_cartesc_set_dims
241 
242  call file_history_cartesc_set_axes
243 
244  call file_history_cartesc_set_axes_attributes
245 
247  file_history_truncate_2d => file_history_cartesc_truncate_2d
249 
250  return
251  end subroutine file_history_cartesc_setup
252 
253  !-----------------------------------------------------------------------------
255  !-----------------------------------------------------------------------------
256  subroutine file_history_cartesc_set_pres( &
257  PRES, &
258  PRESH, &
259  SFC_PRES )
260  use scale_interp_vert, only: &
262  implicit none
263 
264  real(rp), intent(in) :: pres (:,:,:) ! pressure at the full level [Pa]
265  real(rp), intent(in) :: presh (:,:,:) ! pressure at the half level [Pa]
266  real(rp), intent(in) :: sfc_pres( :,:) ! surface pressure [Pa]
267  !---------------------------------------------------------------------------
268 
269  if ( file_history_cartesc_pres_nlayer > 0 ) then
270  call interp_vert_setcoef_pres( file_history_cartesc_pres_nlayer, & ! [IN]
271  ka, ks, ke, & ! [IN]
272  ia, isb, ieb, & ! [IN]
273  ja, jsb, jeb, & ! [IN]
274  pres(:,:,:), & ! [IN]
275  presh(:,:,:), & ! [IN]
276  sfc_pres(:,:) , & ! [IN]
277  file_history_cartesc_pres_val(:) ) ! [IN]
278  pres_set = .true.
279  endif
280 
281  return
282  end subroutine file_history_cartesc_set_pres
283 
284  ! private routines
285 
286  !-----------------------------------------------------------------------------
288  !-----------------------------------------------------------------------------
289  subroutine file_history_cartesc_set_dims
290  use scale_prc, only: &
291  prc_myrank
292  use scale_prc_cartesc, only: &
293  prc_2drank, &
294  prc_has_w, &
295  prc_has_s, &
296  prc_periodic_x, &
298  use scale_file_history, only: &
300  use scale_mapprojection, only: &
302  implicit none
303 
304  character(len=H_SHORT) :: mapping
305 
306  character(len=H_SHORT) :: dims(3,3)
307 
308  integer :: start(3,3), count(3,3)
309  integer :: xs, xc, ys, yc
310  integer :: xoff, yoff
311  !---------------------------------------------------------------------------
312 
313  ! get start and count for x and y
314  if ( file_history_cartesc_boundary ) then
315  xs = isgb
316  ys = jsgb
317  xc = imaxb
318  yc = jmaxb
319  xoff = 0
320  yoff = 0
321  else
322  ! for the case the shared-file contains no halos
323  xs = prc_2drank(prc_myrank,1) * imax + 1 ! no IHALO
324  xc = imax
325  ys = prc_2drank(prc_myrank,2) * jmax + 1 ! no JHALO
326  yc = jmax
327  if ( prc_periodic_x ) then
328  xoff = 0
329  else
330  xoff = 1
331  end if
332  if ( prc_periodic_y ) then
333  yoff = 0
334  else
335  yoff = 1
336  end if
337  end if
338 
339  ! get mapping name
340  mapping = mapprojection_mappinginfo%mapping_name
341 
342 
343  ! 0D
344  dims(1,1) = ""
345  start(1,1) = 1
346  count(1,1) = 1
347  call file_history_set_dim( "", 0, 1, dims(:,:), zs(:), start(:,:), count(:,:) ) ! [IN]
348 
349  ! Vertical 1D
350  start(1,1) = 1
351  dims(1,1) = "z"
352  count(1,1) = file_history_cartesc_model_nlayer
353  call file_history_set_dim( "Z", 1, 1, dims(:,:), zs(:), start(:,:), count(:,:) ) ! [IN]
354  dims(1,1) = "zh"
355  count(1,1) = file_history_cartesc_model_nlayer + 1
356  call file_history_set_dim( "ZH", 1, 1, dims(:,:), zs(:), start(:,:), count(:,:) ) ! [IN]
357 
358  dims(1,1) = "oz"
359  count(1,1) = okmax
360  call file_history_set_dim( "OZ", 1, 1, dims(:,:), zs(:), start(:,:), count(:,:) ) ! [IN]
361  dims(1,1) = "ozh"
362  count(1,1) = okmax + 1
363  call file_history_set_dim( "OZH", 1, 1, dims(:,:), zs(:), start(:,:), count(:,:) ) ! [IN]
364 
365  dims(1,1) = "lz"
366  count(1,1) = lkmax
367  call file_history_set_dim( "LZ", 1, 1, dims(:,:), zs(:), start(:,:), count(:,:) ) ! [IN]
368  dims(1,1) = "lzh"
369  count(1,1) = lkmax + 1
370  call file_history_set_dim( "LZH", 1, 1, dims(:,:), zs(:), start(:,:), count(:,:) ) ! [IN]
371 
372  dims(1,1) = "uz"
373  count(1,1) = ukmax
374  call file_history_set_dim( "UZ", 1, 1, dims(:,:), zs(:), start(:,:), count(:,:) ) ! [IN]
375  dims(1,1) = "uzh"
376  count(1,1) = ukmax + 1
377  call file_history_set_dim( "UZH", 1, 1, dims(:,:), zs(:), start(:,:), count(:,:) ) ! [IN]
378 
379  ! X, Y
380  start(1,:) = xs
381  start(2,:) = ys
382  dims(1,:) = 'lon'
383  dims(2,:) = 'lat'
384  count(1,:) = xc
385  count(2,:) = yc
386  call file_history_set_dim( "XY", 2, 1, dims(:,:), zs(:), start(:,:), count(:,:), mapping=mapping, area="cell_area", location="face" ) ! [IN]
387 
388  start(3,:) = 1
389  dims(3,:) = (/ "height ", "z ", "pressure " /)
390  count(3,:) = (/ file_history_cartesc_model_nlayer, file_history_cartesc_model_nlayer, file_history_cartesc_pres_nlayer /)
391  call file_history_set_dim( "ZXY", 3, nzs, dims(:,:), zs(:), start(:,:), count(:,:), mapping=mapping, & ! [IN]
392  area="cell_area", area_x="cell_area_xyz_x", area_y="cell_area_xyz_y", volume="cell_volume", location="face" ) ! [IN]
393  dims(3,:) = (/ "height_xyw", "zh ", "pressure " /)
394  count(3,:) = (/ file_history_cartesc_model_nlayer+1, file_history_cartesc_model_nlayer+1, file_history_cartesc_pres_nlayer /)
395  call file_history_set_dim( "ZHXY", 3, nzs, dims(:,:), zs(:), start(:,:), count(:,:), mapping=mapping, area="cell_area", volume="cell_volume_xyw", location="face" )
396 
397  if ( okmax > 0 ) then
398  dims(3,1) = "oz"
399  count(3,1) = okmax
400  call file_history_set_dim( "OXY", 3, 1, dims(:,:), zs(:), start(:,:), count(:,:), mapping=mapping, area="cell_area", volume="cell_volume_xyo", location="face", grid="ocean" ) ! [IN]
401  dims(3,1) = "ozh"
402  count(3,1) = okmax + 1
403  call file_history_set_dim( "OHXY", 3, 1, dims(:,:), zs(:), start(:,:), count(:,:), mapping=mapping, area="cell_area", location="face", grid="ocean" ) ! [IN]
404  endif
405 
406  if ( lkmax > 0 ) then
407  dims(3,1) = "lz"
408  count(3,1) = lkmax
409  call file_history_set_dim( "LXY", 3, 1, dims(:,:), zs(:), start(:,:), count(:,:), mapping=mapping, area="cell_area", volume="cell_volume_xyl", location="face", grid="land" ) ! [IN]
410  dims(3,1) = "lzh"
411  count(3,1) = lkmax + 1
412  call file_history_set_dim( "LHXY", 3, 1, dims(:,:), zs(:), start(:,:), count(:,:), mapping=mapping, area="cell_area", location="face", grid="land" ) ! [IN]
413  endif
414 
415  if ( ukmax > 0 ) then
416  dims(3,1) = "uz"
417  count(3,1) = ukmax
418  call file_history_set_dim( "UXY", 3, 1, dims(:,:), zs(:), start(:,:), count(:,:), mapping=mapping, area="cell_area", volume="cell_volume_xyu", grid="urban" ) ! [IN]
419  dims(3,1) = "uzh"
420  count(3,1) = ukmax + 1
421  call file_history_set_dim( "UHXY", 3, 1, dims(:,:), zs(:), start(:,:), count(:,:), mapping=mapping, area="cell_area", location="face", grid="urban" ) ! [IN]
422  endif
423 
424  ! XH, Y
425  dims(1,:) = 'lon_uy'
426  dims(2,:) = 'lat_uy'
427  if ( prc_has_w ) then
428  start(1,:) = xs+xoff
429  count(1,:) = xc
430  else
431  start(1,:) = xs
432  count(1,:) = xc+xoff
433  endif
434  call file_history_set_dim( "XHY", 2, 1, dims(:,:), zs(:), start(:,:), count(:,:), mapping=mapping, area="cell_area_uy", location="edge1" ) ! [IN]
435 
436  dims(3,:) = (/ "height_uyz", "z ", "pressure " /)
437  count(3,:) = (/ file_history_cartesc_model_nlayer, file_history_cartesc_model_nlayer, file_history_cartesc_pres_nlayer /)
438  call file_history_set_dim( "ZXHY", 3, nzs, dims(:,:), zs(:), start(:,:), count(:,:), mapping=mapping, area="cell_area_uy", volume="cell_volume_uyz", location="edge1" ) ! [IN]
439  dims(3,:) = (/ "height_uyw", "zh ", "pressure " /)
440  count(3,:) = (/ file_history_cartesc_model_nlayer+1, file_history_cartesc_model_nlayer+1, file_history_cartesc_pres_nlayer /)
441  call file_history_set_dim( "ZHXHY", 3, nzs, dims(:,:), zs(:), start(:,:), count(:,:), mapping=mapping, area="cell_area_uy", location="edge1" ) ! [IN]
442 
443  ! X, YH
444  dims(1,:) = 'lon_xv'
445  dims(2,:) = 'lat_xv'
446  start(1,:) = xs
447  count(1,:) = xc
448  if ( prc_has_s ) then
449  start(2,:) = ys+yoff
450  count(2,:) = yc
451  else
452  start(2,:) = ys
453  count(2,:) = yc+yoff
454  endif
455  call file_history_set_dim( "XYH", 2, 1, dims(:,:), zs(:), start(:,:), count(:,:), mapping=mapping, area="cell_area_xv", location="edge2" ) ! [IN]
456 
457  dims(3,:) = (/ "height_xvz", "z ", "pressure " /)
458  count(3,:) = (/ file_history_cartesc_model_nlayer, file_history_cartesc_model_nlayer, file_history_cartesc_pres_nlayer /)
459  call file_history_set_dim( "ZXYH", 3, nzs, dims(:,:), zs(:), start(:,:), count(:,:), mapping=mapping, area="cell_area_xv", volume="cell_volume_xvz", location="edge2" ) ! [IN]
460  dims(3,:) = (/ "height_xvw", "zh ", "pressure " /)
461  count(3,:) = (/ file_history_cartesc_model_nlayer+1, file_history_cartesc_model_nlayer+1, file_history_cartesc_pres_nlayer /)
462  call file_history_set_dim( "ZHXYH", 3, nzs, dims(:,:), zs(:), start(:,:), count(:,:), mapping=mapping, area="cell_area_xv", location="edge2" ) ! [IN]
463 
464  ! XH, YH
465  dims(1,:) = 'lon_uv'
466  dims(2,:) = 'lat_uv'
467  if ( prc_has_w ) then
468  start(1,:) = xs+xoff
469  count(1,:) = xc
470  else
471  start(1,:) = xs
472  count(1,:) = xc+xoff
473  endif
474  call file_history_set_dim( "XHYH", 2, 1, dims(:,:), zs(:), start(:,:), count(:,:), mapping=mapping, location="face" ) ! [IN]
475 
476  dims(3,:) = (/ "height_uvz", "z ", "pressure " /)
477  count(3,:) = (/ file_history_cartesc_model_nlayer, file_history_cartesc_model_nlayer, file_history_cartesc_pres_nlayer /)
478  call file_history_set_dim( "ZXHYH", 3, nzs, dims(:,:), zs(:), start(:,:), count(:,:), mapping=mapping, & ![IN]
479  area="cell_area_uv", area_x="cell_area_uvz_x", area_y="cell_area_uvz_y", location="node" ) ! [IN]
480  dims(3,:) = (/ "height_uvw", "zh ", "pressure " /)
481  count(3,:) = (/ file_history_cartesc_model_nlayer+1, file_history_cartesc_model_nlayer+1, file_history_cartesc_pres_nlayer /)
482  call file_history_set_dim( "ZHXHYH", 3, nzs, dims(:,:), zs(:), start(:,:), count(:,:), mapping=mapping, location="node" ) ! [IN]
483 
484  return
485  end subroutine file_history_cartesc_set_dims
486 
487  !-----------------------------------------------------------------------------
490  src, &
491  dim_type, zcoord, fill_halo, &
492  dst )
493  implicit none
494 
495  real(RP), intent(in) :: src(:)
496  character(len=*), intent(in) :: dim_type
497  character(len=*), intent(in) :: zcoord
498  logical, intent(in) :: fill_halo ! ignored
499 
500  real(DP), intent(out) :: dst(:)
501 
502  integer :: ksize
503  integer :: kstart
504  integer :: k
505  !---------------------------------------------------------------------------
506 
507  ! select dimension
508  select case ( dim_type )
509  case ('Z')
510  ksize = file_history_cartesc_model_nlayer
511  kstart = ks
512  case ('ZH')
513  ksize = file_history_cartesc_model_nlayer+1
514  kstart = ks-1
515  case ('OZ')
516  ksize = okmax
517  kstart = oks
518  case ('OZH')
519  ksize = okmax+1
520  kstart = oks-1
521  case ('LZ')
522  ksize = lkmax
523  kstart = lks
524  case ('LZH')
525  ksize = lkmax+1
526  kstart = lks-1
527  case ('UZ')
528  ksize = ukmax
529  kstart = uks
530  case ('UZH')
531  ksize = ukmax+1
532  kstart = uks-1
533  case default
534  log_error("FILE_HISTORY_CARTESC_truncate_1D",*) 'dim_type is invalid: ', trim(dim_type)
535  call prc_abort
536  end select
537 
538  do k = 1, ksize
539  dst(k) = src(kstart+k-1)
540  enddo
541 
542  return
543  end subroutine file_history_cartesc_truncate_1d
544 
545  !-----------------------------------------------------------------------------
547  subroutine file_history_cartesc_truncate_2d( &
548  src, &
549  dim_type, zcoord, fill_halo, &
550  dst )
551  use scale_file_h, only: &
552  rmiss => file_rmiss
553  implicit none
554 
555  real(RP), intent(in) :: src(:,:)
556  character(len=*), intent(in) :: dim_type
557  character(len=*), intent(in) :: zcoord ! ignored
558  logical, intent(in) :: fill_halo
559 
560  real(DP), intent(out) :: dst(:)
561 
562  integer :: isize, jsize
563  integer :: istart, jstart
564  integer :: i, j
565  !---------------------------------------------------------------------------
566 
567  ! select dimension
568  select case( dim_type )
569  case ( 'XY', 'XYH' )
570  isize = im
571  istart = ims
572  case ( 'XHY', 'XHYH' )
573  isize = imh
574  istart = imsh
575  case default
576  log_error("FILE_HISTORY_CARTESC_truncate_2D",*) 'dim_type is invalid: ', trim(dim_type)
577  call prc_abort
578  end select
579 
580  select case ( dim_type )
581  case ( 'XY', 'XHY' )
582  jsize = jm
583  jstart = jms
584  case ( 'XYH', 'XHYH' )
585  jsize = jmh
586  jstart = jmsh
587  case default
588  log_error("FILE_HISTORY_CARTESC_truncate_2D",*) 'dim_type is invalid: ', trim(dim_type)
589  call prc_abort
590  end select
591 
592  !$omp parallel do
593  do j = 1, jsize
594  do i = 1, isize
595  dst((j-1)*isize+i) = src(istart+i-1,jstart+j-1)
596  enddo
597  enddo
598 
599  if ( fill_halo ) then
600  ! W halo
601  do j = 1, jsize
602  do i = 1, is-istart
603  dst((j-1)*isize+i) = rmiss
604  enddo
605  enddo
606  ! E halo
607  do j = 1, jsize
608  do i = ie-istart+2, ime-istart+1
609  dst((j-1)*isize+i) = rmiss
610  enddo
611  enddo
612  ! S halo
613  do j = 1, js-jstart
614  do i = 1, isize
615  dst((j-1)*isize+i) = rmiss
616  enddo
617  enddo
618  ! N halo
619  do j = je-jstart+2, jme-jstart+1
620  do i = 1, isize
621  dst((j-1)*isize+i) = rmiss
622  enddo
623  enddo
624 
625  end if
626 
627  return
628  end subroutine file_history_cartesc_truncate_2d
629 
630  !-----------------------------------------------------------------------------
633  src, &
634  dim_type, zcoord, fill_halo, &
635  dst )
636  use scale_file_h, only: &
637  rmiss => file_rmiss
638  use scale_const, only: &
639  undef => const_undef
640  use scale_atmos_grid_cartesc, only: &
643  use scale_atmos_grid_cartesc_real, only: &
646  use scale_interp_vert, only: &
652  implicit none
653 
654  real(RP), intent(in) :: src(:,:,:)
655  character(len=*), intent(in) :: dim_type
656  character(len=*), intent(in) :: zcoord
657  logical, intent(in) :: fill_halo
658  real(DP), intent(out) :: dst(:)
659 
660  real(RP) :: src_Z(KA,IA,JA)
661  real(RP) :: src_P(FILE_HISTORY_CARTESC_PRES_nlayer,IA,JA)
662 
663  logical :: atmos
664 
665  integer :: isize, jsize, ksize
666  integer :: istart, jstart, kstart
667  integer :: i, j, k
668  !---------------------------------------------------------------------------
669 
670  ! select dimension
671  if ( index( dim_type, 'XH' ) == 0 ) then
672  isize = im
673  istart = ims
674  else
675  isize = imh
676  istart = imsh
677  end if
678 
679  if ( index( dim_type, 'YH' ) == 0 ) then
680  jsize = jm
681  jstart = jms
682  else
683  jsize = jmh
684  jstart = jmsh
685  end if
686 
687  select case( dim_type(1:1) )
688  case ( 'Z' )
689  ksize = file_history_cartesc_model_nlayer
690  kstart = ks
691  atmos = .true.
692  case('O')
693  ksize = okmax
694  kstart = oks
695  atmos = .false.
696  case('L')
697  ksize = lkmax
698  kstart = lks
699  atmos = .false.
700  case('U')
701  ksize = ukmax
702  kstart = uks
703  atmos = .false.
704  case default
705  log_error("FILE_HISTORY_CARTESC_truncate_3D",*) 'dim_type is invalid: ', trim(dim_type)
706  call prc_abort
707  end select
708  if ( dim_type(2:2) == 'H' ) then
709  ksize = ksize + 1
710  kstart = kstart - 1
711  end if
712 
713 
714  if ( ksize == file_history_cartesc_model_nlayer .and. zcoord == "z" .and. interp_available .and. atmos ) then ! z*->z interpolation (full level)
715 
716  call prof_rapstart('FILE_O_interp', 2)
717  call interp_vert_xi2z( ka, ks, ke, ia, isb, ieb, ja, jsb, jeb, &
718  atmos_grid_cartesc_cz(:), & ! [IN]
719  atmos_grid_cartesc_real_cz(:,:,:), & ! [IN]
720  src(:,:,:), & ! [IN]
721  src_z(:,:,:) ) ! [OUT]
722  call prof_rapend ('FILE_O_interp', 2)
723 
724  !$omp parallel do
725  do k = 1, ksize
726  do j = 1, jsize
727  do i = 1, isize
728  dst((k-1)*jsize*isize+(j-1)*isize+i) = src_z(kstart+k-1,istart+i-1,jstart+j-1)
729  enddo
730  enddo
731  enddo
732 
733  else if( ksize == file_history_cartesc_model_nlayer+1 .and. zcoord == "z" .and. interp_available .and. atmos ) then ! z*->z interpolation (half level)
734 
735 
736  call prof_rapstart('FILE_O_interp', 2)
737  call interp_vert_xih2zh( ka, ks, ke, ia, isb, ieb, ja, jsb, jeb, &
738  atmos_grid_cartesc_fz(:), & ! [IN]
739  atmos_grid_cartesc_real_fz(:,:,:), & ! [IN]
740  src(:,:,:), & ! [IN]
741  src_z(:,:,:) ) ! [OUT]
742  call prof_rapend ('FILE_O_interp', 2)
743 
744  !$omp parallel do
745  do k = 1, ksize
746  do j = 1, jsize
747  do i = 1, isize
748  dst((k-1)*jsize*isize+(j-1)*isize+i) = src_z(kstart+k-1,istart+i-1,jstart+j-1)
749  enddo
750  enddo
751  enddo
752 
753  elseif( zcoord == "pressure" .and. ksize == file_history_cartesc_model_nlayer .and. atmos ) then ! z*->p interpolation (full level)
754  ksize = file_history_cartesc_pres_nlayer
755  if ( ksize == 0 ) then
756  log_error("FILE_HISTORY_CARTESC_truncate_3D",*) 'FILE_HISTORY_CARTESC_PRES_nlayer must be set to output variable with the pressure coordinate'
757  call prc_abort
758  end if
759 
760  if ( pres_set ) then
761 
762  call prof_rapstart('FILE_O_interp', 2)
763  call interp_vert_xi2p( file_history_cartesc_pres_nlayer, & ! [IN]
764  ka, ks, ke, & ! [IN]
765  ia, isb, ieb, & ! [IN]
766  ja, jsb, jeb, & ! [IN]
767  src(:,:,:), & ! [IN]
768  src_p(:,:,:) ) ! [OUT]
769  call prof_rapend ('FILE_O_interp', 2)
770 
771  !$omp parallel do
772  do k = 1, ksize
773  do j = 1, jsize
774  do i = 1, isize
775  dst((k-1)*jsize*isize+(j-1)*isize+i) = src_p(k,istart+i-1,jstart+j-1)
776  enddo
777  enddo
778  enddo
779 
780  else
781 
782  !$omp parallel do
783  do k = 1, ksize
784  do j = 1, jsize
785  do i = 1, isize
786  dst((k-1)*jsize*isize+(j-1)*isize+i) = undef
787  enddo
788  enddo
789  enddo
790 
791  end if
792 
793  elseif( zcoord == "pressure" .and. ksize == file_history_cartesc_model_nlayer+1 .and. atmos ) then ! z*->p interpolation (half level)
794  ksize = file_history_cartesc_pres_nlayer
795  if ( ksize == 0 ) then
796  log_error("FILE_HISTORY_CARTESC_truncate_3D",*) 'FILE_HISTORY_CARTESC_PRES_nlayer must be set to output variable with the pressure coordinate'
797  call prc_abort
798  end if
799 
800  if ( pres_set ) then
801 
802  call prof_rapstart('FILE_O_interp', 2)
803  call interp_vert_xih2p( file_history_cartesc_pres_nlayer, & ! [IN]
804  ka, ks, ke, & ! [IN]
805  ia, isb, ieb, & ! [IN]
806  ja, jsb, jeb, & ! [IN]
807  src(:,:,:), & ! [IN]
808  src_p(:,:,:) ) ! [OUT]
809  call prof_rapend ('FILE_O_interp', 2)
810 
811  !$omp parallel do
812  do k = 1, ksize
813  do j = 1, jsize
814  do i = 1, isize
815  dst((k-1)*jsize*isize+(j-1)*isize+i) = src_p(k,istart+i-1,jstart+j-1)
816  enddo
817  enddo
818  enddo
819 
820  else
821 
822  !$omp parallel do
823  do k = 1, ksize
824  do j = 1, jsize
825  do i = 1, isize
826  dst((k-1)*jsize*isize+(j-1)*isize+i) = undef
827  enddo
828  enddo
829  enddo
830 
831  end if
832 
833  else ! no interpolation
834 
835  do k = 1, ksize
836  do j = 1, jsize
837  do i = 1, isize
838  dst((k-1)*jsize*isize+(j-1)*isize+i) = src(kstart+k-1,istart+i-1,jstart+j-1)
839  enddo
840  enddo
841  enddo
842 
843  endif
844 
845  if ( fill_halo ) then
846  ! W halo
847  do k = 1, ksize
848  do j = 1, jsize
849  do i = 1, is-istart
850  dst((k-1)*jsize*isize+(j-1)*isize+i) = rmiss
851  enddo
852  enddo
853  enddo
854  ! E halo
855  do k = 1, ksize
856  do j = 1, jsize
857  do i = ie-istart+2, ime-istart+1
858  dst((k-1)*jsize*isize+(j-1)*isize+i) = rmiss
859  enddo
860  enddo
861  enddo
862  ! S halo
863  do k = 1, ksize
864  do j = 1, js-jstart
865  do i = 1, isize
866  dst((k-1)*jsize*isize+(j-1)*isize+i) = rmiss
867  enddo
868  enddo
869  enddo
870  ! N halo
871  do k = 1, ksize
872  do j = je-jstart+2, jme-jstart+1
873  do i = 1, isize
874  dst((k-1)*jsize*isize+(j-1)*isize+i) = rmiss
875  enddo
876  enddo
877  enddo
878  endif
879 
880  return
881  end subroutine file_history_cartesc_truncate_3d
882 
883  !-----------------------------------------------------------------------------
885  ! only register the axis and coordinate variables into internal buffers
886  ! The actual write happens later when calling FILE_HISTORY_CARTESC_write
887  subroutine file_history_cartesc_set_axes
888  use scale_file_history, only: &
891  file_history_set_associatedcoordinate
892  use scale_const, only: &
893  undef => const_undef, &
894  d2r => const_d2r
895  use scale_prc, only: &
896  prc_myrank
897  use scale_prc_cartesc, only: &
898  prc_2drank, &
899  prc_num_x, &
900  prc_num_y, &
901  prc_periodic_x, &
902  prc_periodic_y, &
903  prc_twod, &
904  prc_has_w, &
905  prc_has_s
906  use scale_atmos_grid_cartesc, only: &
937  use scale_ocean_grid_cartesc, only: &
941  use scale_land_grid_cartesc, only: &
945  use scale_urban_grid_cartesc, only: &
949  use scale_atmos_grid_cartesc_real, only: &
961  areauy => atmos_grid_cartesc_real_areauy, &
962  areaxv => atmos_grid_cartesc_real_areaxv, &
963  areazuy_x => atmos_grid_cartesc_real_areazuy_x, &
964  areazxv_y => atmos_grid_cartesc_real_areazxv_y, &
965  areawuy_x => atmos_grid_cartesc_real_areawuy_x, &
966  areawxv_y => atmos_grid_cartesc_real_areawxv_y, &
967  areazxy_x => atmos_grid_cartesc_real_areazxy_x, &
968  areazuv_y => atmos_grid_cartesc_real_areazuv_y, &
969  areazuv_x => atmos_grid_cartesc_real_areazuv_x, &
970  areazxy_y => atmos_grid_cartesc_real_areazxy_y, &
972  volwxy => atmos_grid_cartesc_real_volwxy, &
973  volzuy => atmos_grid_cartesc_real_volzuy, &
975  use scale_ocean_grid_cartesc_real, only: &
977  use scale_land_grid_cartesc_real, only: &
979  use scale_urban_grid_cartesc_real, only: &
981  use scale_topography, only: &
983  use scale_landuse, only: &
985  implicit none
986 
987  real(RP) :: AXIS (imh,jmh,0:FILE_HISTORY_CARTESC_MODEL_nlayer)
988  real(RP) :: AXISO(im, jm, OKMAX)
989  real(RP) :: AXISL(im, jm, LKMAX)
990  real(RP) :: AXISU(im, jm, UKMAX)
991  character(len=2) :: AXIS_name(3)
992 
993  integer :: rankidx(2)
994  integer :: start(3,4)
995  integer :: startX, startY, startZ
996  integer :: startXH, startYH
997  integer :: XAG, YAG
998  integer :: XAGH, YAGH
999 
1000  real(RP) :: z_bnds(2,KA), zh_bnds(2,0:KA)
1001  real(RP) :: oz_bnds(2,OKA), ozh_bnds(2,0:OKA)
1002  real(RP) :: lz_bnds(2,LKA), lzh_bnds(2,0:LKA)
1003  real(RP) :: uz_bnds(2,UKA), uzh_bnds(2,0:UKA)
1004  real(RP) :: x_bnds(2,IA), xh_bnds(2,0:IA)
1005  real(RP) :: y_bnds(2,JA), yh_bnds(2,0:JA)
1006 
1007  real(RP) :: FDXG(0:IAG), FDYG(0:JAG)
1008  real(RP) :: FDX(0:IA), FDY(0:JA)
1009 
1010  integer :: k, i, j
1011  !---------------------------------------------------------------------------
1012 
1013  rankidx(1) = prc_2drank(prc_myrank,1)
1014  rankidx(2) = prc_2drank(prc_myrank,2)
1015 
1016  ! For parallel I/O, some variables are written by a subset of processes.
1017  ! 1. Only PRC_myrank 0 writes all z axes
1018  ! 2. Only south-most processes (rankidx(2) == 0) write x axes
1019  ! rankidx(1) == 0 writes west HALO
1020  ! rankidx(1) == PRC_NUM_X-1 writes east HALO
1021  ! others writes without HALO
1022  ! 3. Only west-most processes (rankidx(1) == 0) write y axes
1023  ! rankidx(1) == 0 writes south HALO
1024  ! rankidx(1) == PRC_NUM_Y-1 writes north HALO
1025  ! others writes without HALO
1026 
1027  if ( file_history_aggregate ) then
1028 
1029  startz = 1
1030 
1031  if ( file_history_cartesc_boundary ) then
1032  startx = isgb ! global subarray starting index
1033  starty = jsgb ! global subarray starting index
1034  startxh = startx
1035  startyh = starty
1036  xag = iagb
1037  yag = jagb
1038  xagh = xag
1039  yagh = yag
1040  else
1041  startx = prc_2drank(prc_myrank,1) * imax + 1
1042  starty = prc_2drank(prc_myrank,2) * jmax + 1
1043  startxh = startx
1044  startyh = starty
1045  xag = imaxg
1046  yag = jmaxg
1047  xagh = xag
1048  yagh = yag
1049 
1050  if ( .NOT. prc_periodic_x ) then
1051  xagh = xagh + 1
1052  if( prc_has_w ) startxh = startxh + 1
1053  endif
1054 
1055  if ( .NOT. prc_periodic_y ) then
1056  yagh = yagh + 1
1057  if( prc_has_s ) startyh = startyh + 1
1058  endif
1059  endif
1060 
1061  ! for shared-file parallel I/O, only a part of rank writes variables
1062  if ( prc_myrank > 0 ) then ! only rank 0 writes Z axes
1063  startz = -1
1064  endif
1065  if ( rankidx(2) > 0 ) then ! only south-most processes write
1066  startx = -1
1067  startxh = -1
1068  endif
1069  if ( rankidx(1) > 0 ) then ! only west-most processes write
1070  starty = -1
1071  startyh = -1
1072  endif
1073  else
1074  startz = 1
1075  startx = 1
1076  starty = 1
1077  startxh = startx
1078  startyh = starty
1079  xag = im
1080  yag = jm
1081  xagh = imh
1082  yagh = jmh
1083  endif
1084 
1085  ! bounds
1086  do k = ks, ke
1087  z_bnds(1,k) = atmos_grid_cartesc_fz(k-1)
1088  z_bnds(2,k) = atmos_grid_cartesc_fz(k )
1089  end do
1090  do k = ks-1, ke
1091  zh_bnds(1,k) = atmos_grid_cartesc_cz(k )
1092  zh_bnds(2,k) = atmos_grid_cartesc_cz(k+1)
1093  end do
1094 
1095  if ( oka > 0 ) then
1096  do k = oks, oke
1097  oz_bnds(1,k) = ocean_grid_cartesc_fz(k-1)
1098  oz_bnds(2,k) = ocean_grid_cartesc_fz(k )
1099  end do
1100  ozh_bnds(1,oks-1) = ocean_grid_cartesc_fz(oks-1)
1101  do k = oks-1, oke-1
1102  ozh_bnds(2,k ) = ocean_grid_cartesc_cz(k+1)
1103  ozh_bnds(1,k+1) = ocean_grid_cartesc_cz(k+1)
1104  end do
1105  ozh_bnds(2,oke) = ocean_grid_cartesc_fz(oke)
1106  end if
1107 
1108  if ( lka > 0 ) then
1109  do k = lks, lke
1110  lz_bnds(1,k) = land_grid_cartesc_fz(k-1)
1111  lz_bnds(2,k) = land_grid_cartesc_fz(k )
1112  end do
1113  lzh_bnds(1,lks-1) = land_grid_cartesc_fz(lks-1)
1114  do k = lks-1, lke-1
1115  lzh_bnds(2,k ) = land_grid_cartesc_cz(k+1)
1116  lzh_bnds(1,k+1) = land_grid_cartesc_cz(k+1)
1117  end do
1118  lzh_bnds(2,lke) = land_grid_cartesc_fz(lke)
1119  end if
1120 
1121  if ( uka > 0 ) then
1122  do k = uks, uke
1123  uz_bnds(1,k) = urban_grid_cartesc_fz(k-1)
1124  uz_bnds(2,k) = urban_grid_cartesc_fz(k )
1125  end do
1126  uzh_bnds(1,uks-1) = urban_grid_cartesc_fz(uks-1)
1127  do k = uks-1, uke-1
1128  uzh_bnds(2,k ) = urban_grid_cartesc_cz(k+1)
1129  uzh_bnds(1,k+1) = urban_grid_cartesc_cz(k+1)
1130  end do
1131  uzh_bnds(2,uke) = urban_grid_cartesc_fz(uke)
1132  end if
1133 
1134  do i = ims, ime
1135  x_bnds(1,i) = atmos_grid_cartesc_fx(i-1)
1136  x_bnds(2,i) = atmos_grid_cartesc_fx(i )
1137  end do
1138  if ( imsh == 0 ) then
1139  xh_bnds(1,0) = atmos_grid_cartesc_fx(0)
1140  else
1141  xh_bnds(1,imsh) = atmos_grid_cartesc_cx(imsh)
1142  end if
1143  do i = imsh, ime-1
1144  xh_bnds(2,i ) = atmos_grid_cartesc_cx(i+1)
1145  xh_bnds(1,i+1) = atmos_grid_cartesc_cx(i+1)
1146  end do
1147  if ( ime == ia ) then
1148  xh_bnds(2,ime) = atmos_grid_cartesc_fx(ia)
1149  else
1150  xh_bnds(2,ime) = atmos_grid_cartesc_cx(ime+1)
1151  end if
1152 
1153  do j = jms, jme
1154  y_bnds(1,j) = atmos_grid_cartesc_fy(j-1)
1155  y_bnds(2,j) = atmos_grid_cartesc_fy(j )
1156  end do
1157  if ( jmsh == 0 ) then
1158  yh_bnds(1,jmsh) = atmos_grid_cartesc_fy(jmsh)
1159  else
1160  yh_bnds(1,jmsh) = atmos_grid_cartesc_cy(jmsh)
1161  end if
1162  do j = jmsh, jme-1
1163  yh_bnds(2,j ) = atmos_grid_cartesc_cy(j+1)
1164  yh_bnds(1,j+1) = atmos_grid_cartesc_cy(j+1)
1165  end do
1166  if ( jme == ja ) then
1167  yh_bnds(2,jme) = atmos_grid_cartesc_fy(jme)
1168  else
1169  yh_bnds(2,jme) = atmos_grid_cartesc_cy(jme+1)
1170  end if
1171 
1172 
1173  do i = 1, iag-1
1174  fdxg(i) = atmos_grid_cartesc_fdxg(i)
1175  enddo
1176  fdxg(0) = undef
1177  fdxg(iag) = undef
1178  do j = 1, jag-1
1179  fdyg(j) = atmos_grid_cartesc_fdyg(j)
1180  enddo
1181  fdyg(0) = undef
1182  fdyg(jag) = undef
1183 
1184  do i = 1, ia-1
1185  fdx(i) = atmos_grid_cartesc_fdx(i)
1186  enddo
1187  fdx(0) = fdxg(is_ing-ihalo-1)
1188  fdx(ia) = fdxg(ie_ing+ihalo )
1189  do j = 1, ja-1
1190  fdy(j) = atmos_grid_cartesc_fdy(j)
1191  enddo
1192  fdy(0) = fdyg(js_ing-jhalo-1)
1193  fdy(ja) = fdyg(je_ing+jhalo )
1194 
1195 
1196  ! for the shared-file I/O method, the axes are global (gsize)
1197  ! for one-file-per-process I/O method, the axes size is equal to the local buffer size
1198 
1199  call file_history_set_axis( 'z', 'Z', 'm', 'z', atmos_grid_cartesc_cz(ks :ks+file_history_cartesc_model_nlayer-1), &
1200  bounds=z_bnds(:,ks :ks+file_history_cartesc_model_nlayer-1), gsize=file_history_cartesc_model_nlayer, start=startz )
1201  call file_history_set_axis( 'zh', 'Z (half level)', 'm', 'zh', atmos_grid_cartesc_fz(ks-1:ks+file_history_cartesc_model_nlayer-1), &
1202  bounds=zh_bnds(:,ks-1:ks+file_history_cartesc_model_nlayer-1), gsize=file_history_cartesc_model_nlayer+1, start=startz )
1203 
1204  if ( file_history_cartesc_pres_nlayer > 0 ) then
1205  call file_history_set_axis( 'pressure', 'Pressure', 'hPa', 'pressure', file_history_cartesc_pres_val(:)/100.0_rp, &
1206  gsize=file_history_cartesc_pres_nlayer, start=startz, down=.true. )
1207  endif
1208 
1209  if ( oka > 0 ) then
1210  call file_history_set_axis( 'oz', 'OZ', 'm', 'oz', ocean_grid_cartesc_cz(oks :oke), &
1211  bounds=oz_bnds(:,oks :oke), gsize=okmax , start=startz, down=.true. )
1212  call file_history_set_axis( 'ozh', 'OZ (half level)', 'm', 'ozh', ocean_grid_cartesc_fz(oks-1:oke), &
1213  bounds=ozh_bnds(:,oks-1:oke), gsize=okmax+1, start=startz, down=.true. )
1214  endif
1215 
1216  if ( lka > 0 ) then
1217  call file_history_set_axis( 'lz', 'LZ', 'm', 'lz', land_grid_cartesc_cz(lks :lke), &
1218  bounds=lz_bnds(:,lks :lke), gsize=lkmax , start=startz, down=.true. )
1219  call file_history_set_axis( 'lzh', 'LZ (half level)', 'm', 'lzh', land_grid_cartesc_fz(lks-1:lke), &
1220  bounds=lzh_bnds(:,lks-1:lke), gsize=lkmax+1, start=startz, down=.true. )
1221  endif
1222 
1223  if ( uka > 0 ) then
1224  call file_history_set_axis( 'uz', 'UZ', 'm', 'uz', urban_grid_cartesc_cz(uks :uke), &
1225  bounds=uz_bnds(:,uks :uke), gsize=ukmax , start=startz, down=.true. )
1226  call file_history_set_axis( 'uzh', 'UZ (half level)', 'm', 'uzh', urban_grid_cartesc_fz(uks-1:uke), &
1227  bounds=uzh_bnds(:,uks-1:uke), gsize=ukmax+1, start=startz, down=.true. )
1228  endif
1229 
1230  call file_history_set_axis( 'x', 'X', 'm', 'x', atmos_grid_cartesc_cx(ims :ime), &
1231  bounds=x_bnds(:,ims :ime), gsize=xag , start=startx )
1232  call file_history_set_axis( 'xh', 'X (half level)', 'm', 'xh', atmos_grid_cartesc_fx(imsh:ime), &
1233  bounds=xh_bnds(:,imsh:ime), gsize=xagh, start=startxh )
1234 
1235  call file_history_set_axis( 'y', 'Y', 'm', 'y', atmos_grid_cartesc_cy(jms :jme), &
1236  bounds=y_bnds(:,jms :jme), gsize=yag , start=starty )
1237  call file_history_set_axis( 'yh', 'Y (half level)', 'm', 'yh', atmos_grid_cartesc_fy(jmsh:jme), &
1238  bounds=yh_bnds(:,jmsh:jme), gsize=yagh, start=startyh )
1239 
1240  ! axes below always include halos when written to file regardless of PRC_PERIODIC_X/PRC_PERIODIC_Y
1241  call file_history_set_axis( 'CZ', 'Atmos Grid Center Position Z', 'm', 'CZ', atmos_grid_cartesc_cz, gsize=ka, start=startz )
1242  call file_history_set_axis( 'FZ', 'Atmos Grid Face Position Z', 'm', 'FZ', atmos_grid_cartesc_fz, gsize=ka+1, start=startz )
1243  call file_history_set_axis( 'CDZ', 'Grid Cell length Z', 'm', 'CZ', atmos_grid_cartesc_cdz, gsize=ka, start=startz )
1244  call file_history_set_axis( 'FDZ', 'Grid distance Z', 'm', 'FDZ', atmos_grid_cartesc_fdz, gsize=ka-1, start=startz )
1245  call file_history_set_axis( 'CBFZ', 'Boundary factor Center Z', '1', 'CZ', atmos_grid_cartesc_cbfz, gsize=ka, start=startz )
1246  call file_history_set_axis( 'FBFZ', 'Boundary factor Face Z', '1', 'FZ', atmos_grid_cartesc_fbfz, gsize=ka+1, start=startz )
1247 
1248  if ( okmax > 0 ) then
1249  call file_history_set_axis( 'OCZ', 'Ocean Grid Center Position Z', 'm', 'OCZ', ocean_grid_cartesc_cz, gsize=okmax, start=startz, down=.true. )
1250  call file_history_set_axis( 'OFZ', 'Ocean Grid Face Position Z', 'm', 'OFZ', ocean_grid_cartesc_fz, gsize=okmax+1, start=startz, down=.true. )
1251  call file_history_set_axis( 'OCDZ', 'Ocean Grid Cell length Z', 'm', 'OCZ', ocean_grid_cartesc_cdz, gsize=okmax, start=startz )
1252  endif
1253 
1254  if ( lkmax > 0 ) then
1255  call file_history_set_axis( 'LCZ', 'Land Grid Center Position Z', 'm', 'LCZ', land_grid_cartesc_cz, gsize=lkmax, start=startz, down=.true. )
1256  call file_history_set_axis( 'LFZ', 'Land Grid Face Position Z', 'm', 'LFZ', land_grid_cartesc_fz, gsize=lkmax+1, start=startz, down=.true. )
1257  call file_history_set_axis( 'LCDZ', 'Land Grid Cell length Z', 'm', 'LCZ', land_grid_cartesc_cdz, gsize=lkmax, start=startz )
1258  endif
1259 
1260  if ( ukmax > 0 ) then
1261  call file_history_set_axis( 'UCZ', 'Urban Grid Center Position Z', 'm', 'UCZ', urban_grid_cartesc_cz, gsize=ukmax, start=startz, down=.true. )
1262  call file_history_set_axis( 'UFZ', 'Urban Grid Face Position Z', 'm', 'UFZ', urban_grid_cartesc_fz, gsize=ukmax+1, start=startz, down=.true. )
1263  call file_history_set_axis( 'UCDZ', 'Urban Grid Cell length Z', 'm', 'UCZ', urban_grid_cartesc_cdz, gsize=ukmax, start=startz )
1264  endif
1265 
1266  if ( file_history_aggregate ) then
1267  call file_history_set_axis( 'CX', 'Atmos Grid Center Position X', 'm', 'CX', atmos_grid_cartesc_cxg, gsize=iag, start=startz )
1268  call file_history_set_axis( 'CY', 'Atmos Grid Center Position Y', 'm', 'CY', atmos_grid_cartesc_cyg, gsize=jag, start=startz )
1269  call file_history_set_axis( 'FX', 'Atmos Grid Face Position X', 'm', 'FX', atmos_grid_cartesc_fxg, gsize=iag+1, start=startz )
1270  call file_history_set_axis( 'FY', 'Atmos Grid Face Position Y', 'm', 'FY', atmos_grid_cartesc_fyg, gsize=jag+1, start=startz )
1271  call file_history_set_axis( 'CDX', 'Grid Cell length X', 'm', 'CX', atmos_grid_cartesc_cdxg, gsize=iag, start=startz )
1272  call file_history_set_axis( 'CDY', 'Grid Cell length Y', 'm', 'CY', atmos_grid_cartesc_cdyg, gsize=jag, start=startz )
1273  call file_history_set_axis( 'FDX', 'Grid distance X', 'm', 'FX', fdxg, gsize=iag+1, start=startz )
1274  call file_history_set_axis( 'FDY', 'Grid distance Y', 'm', 'FY', fdyg, gsize=jag+1, start=startz )
1275  call file_history_set_axis( 'CBFX', 'Boundary factor Center X', '1', 'CX', atmos_grid_cartesc_cbfxg, gsize=iag, start=startz )
1276  call file_history_set_axis( 'CBFY', 'Boundary factor Center Y', '1', 'CY', atmos_grid_cartesc_cbfyg, gsize=jag, start=startz )
1277  call file_history_set_axis( 'FBFX', 'Boundary factor Face X', '1', 'FX', atmos_grid_cartesc_fbfxg, gsize=iag+1, start=startz )
1278  call file_history_set_axis( 'FBFY', 'Boundary factor Face Y', '1', 'FY', atmos_grid_cartesc_fbfyg, gsize=jag+1, start=startz )
1279  else
1280  call file_history_set_axis( 'CX', 'Atmos Grid Center Position X', 'm', 'CX', atmos_grid_cartesc_cx )
1281  call file_history_set_axis( 'CY', 'Atmos Grid Center Position Y', 'm', 'CY', atmos_grid_cartesc_cy )
1282  call file_history_set_axis( 'FX', 'Atmos Grid Face Position X', 'm', 'FX', atmos_grid_cartesc_fx )
1283  call file_history_set_axis( 'FY', 'Atmos Grid Face Position Y', 'm', 'FY', atmos_grid_cartesc_fy )
1284  call file_history_set_axis( 'CDX', 'Grid Cell length X', 'm', 'CX', atmos_grid_cartesc_cdx )
1285  call file_history_set_axis( 'CDY', 'Grid Cell length Y', 'm', 'CY', atmos_grid_cartesc_cdy )
1286  call file_history_set_axis( 'FDX', 'Grid distance X', 'm', 'FX', fdx )
1287  call file_history_set_axis( 'FDY', 'Grid distance Y', 'm', 'FY', fdy )
1288  call file_history_set_axis( 'CBFX', 'Boundary factor Center X', '1', 'CX', atmos_grid_cartesc_cbfx )
1289  call file_history_set_axis( 'CBFY', 'Boundary factor Center Y', '1', 'CY', atmos_grid_cartesc_cbfy )
1290  call file_history_set_axis( 'FBFX', 'Boundary factor Face X', '1', 'FX', atmos_grid_cartesc_fbfx )
1291  call file_history_set_axis( 'FBFY', 'Boundary factor Face Y', '1', 'FY', atmos_grid_cartesc_fbfy )
1292  endif
1293 
1294  call file_history_set_axis('CXG', 'Grid Center Position X (global)', 'm', 'CXG', atmos_grid_cartesc_cxg, gsize=iag, start=startz )
1295  call file_history_set_axis('CYG', 'Grid Center Position Y (global)', 'm', 'CYG', atmos_grid_cartesc_cyg, gsize=jag, start=startz )
1296  call file_history_set_axis('FXG', 'Grid Face Position X (global)', 'm', 'FXG', atmos_grid_cartesc_fxg, gsize=iag+1, start=startz )
1297  call file_history_set_axis('FYG', 'Grid Face Position Y (global)', 'm', 'FYG', atmos_grid_cartesc_fyg, gsize=jag+1, start=startz )
1298  call file_history_set_axis('CDXG', 'Grid Cell length X (global)', 'm', 'CXG', atmos_grid_cartesc_cdxg, gsize=iag, start=startz )
1299  call file_history_set_axis('CDYG', 'Grid Cell length Y (global)', 'm', 'CYG', atmos_grid_cartesc_cdyg, gsize=jag, start=startz )
1300  call file_history_set_axis('FDXG', 'Grid distance X (global)', 'm', 'FDXG', fdxg, gsize=iag+1, start=startz )
1301  call file_history_set_axis('FDYG', 'Grid distance Y (global)', 'm', 'FDYG', fdyg, gsize=jag+1, start=startz )
1302  call file_history_set_axis('CBFXG', 'Boundary factor Center X (global)', '1', 'CXG', atmos_grid_cartesc_cbfxg, gsize=iag, start=startz )
1303  call file_history_set_axis('CBFYG', 'Boundary factor Center Y (global)', '1', 'CYG', atmos_grid_cartesc_cbfyg, gsize=jag, start=startz )
1304  call file_history_set_axis('FBFXG', 'Boundary factor Face X (global)', '1', 'FXG', atmos_grid_cartesc_fbfxg, gsize=iag+1, start=startz )
1305  call file_history_set_axis('FBFYG', 'Boundary factor Face Y (global)', '1', 'FYG', atmos_grid_cartesc_fbfyg, gsize=jag+1, start=startz )
1306 
1307 
1308 
1309  ! associate coordinates
1310  if ( file_history_aggregate ) then
1311  if ( file_history_cartesc_boundary ) then
1312  start(1,:) = isgb ! global subarray starting index
1313  start(2,:) = jsgb ! global subarray starting index
1314  else
1315  start(1,:) = prc_2drank(prc_myrank,1) * imax + 1 ! no IHALO
1316  start(2,:) = prc_2drank(prc_myrank,2) * jmax + 1 ! no JHALO
1317  if ( (.NOT. prc_periodic_x) .AND. prc_has_w ) then
1318  start(1,2) = start(1,2) + 1
1319  start(1,4) = start(1,4) + 1
1320  endif
1321  if ( (.NOT. prc_periodic_y) .AND. prc_has_s ) then
1322  start(2,3) = start(2,3) + 1
1323  start(2,4) = start(2,4) + 1
1324  endif
1325  endif
1326  start(3,:) = 1
1327  else
1328  start(:,:) = 1
1329  endif
1330 
1331  do k = 1, file_history_cartesc_model_nlayer
1332  axis(1:im,1:jm,k) = atmos_grid_cartesc_real_cz(k+ks-1,ims:ime,jms:jme)
1333  enddo
1334  axis_name(1:3) = (/'x ', 'y ', 'z '/)
1335  call file_history_set_associatedcoordinate( 'height', 'height above ground level', &
1336  'm', axis_name(1:3), axis(1:im,1:jm,1:), start=start(:,1) )
1337 
1338  do k = 0, file_history_cartesc_model_nlayer
1339  axis(1:im,1:jm,k) = atmos_grid_cartesc_real_fz(k+ks-1,ims:ime,jms:jme)
1340  enddo
1341  axis_name(1:3) = (/'x ', 'y ', 'zh'/)
1342  call file_history_set_associatedcoordinate( 'height_xyw', 'height above ground level (half level xyw)', &
1343  'm' , axis_name(1:3), axis(1:im,1:jm,0:), start=start(:,1) )
1344 
1345  if ( prc_twod ) then
1346  !$omp parallel do
1347  do k = 1, file_history_cartesc_model_nlayer
1348  do j = 1, jm
1349  axis(1,j,k) = atmos_grid_cartesc_real_cz(k+ks-1,is,jms+j-1)
1350  enddo
1351  enddo
1352  else
1353  !$omp parallel do
1354  do k = 1, file_history_cartesc_model_nlayer
1355  do j = 1, jm
1356  do i = 1, min(imh,ia-imsh)
1357  axis(i,j,k) = ( atmos_grid_cartesc_real_cz(k+ks-1,imsh+i-1,jms+j-1) + atmos_grid_cartesc_real_cz(k+ks-1,imsh+i,jms+j-1) ) * 0.5_rp
1358  enddo
1359  enddo
1360  enddo
1361  if ( imh == ia-imsh+1 ) then
1362  !$omp parallel do
1363  do k = 1, file_history_cartesc_model_nlayer
1364  do j = 1, jm
1365  axis(imh,j,k) = atmos_grid_cartesc_real_cz(k+ks-1,imsh+imh-1,jms+j-1)
1366  enddo
1367  enddo
1368  end if
1369  end if
1370  axis_name(1:3) = (/'xh', 'y ', 'z '/)
1371  call file_history_set_associatedcoordinate( 'height_uyz', 'height above ground level (half level uyz)', &
1372  'm', axis_name(1:3), axis(1:imh,1:jm,1:), start=start(:,2) )
1373 
1374  do k = 1, file_history_cartesc_model_nlayer
1375  do j = 1, min(jmh,ja-jmsh)
1376  do i = 1, im
1377  axis(i,j,k) = ( atmos_grid_cartesc_real_cz(k+ks-1,ims+i-1,jmsh+j-1) + atmos_grid_cartesc_real_cz(k+ks-1,ims+i-1,jmsh+j) ) * 0.5_rp
1378  enddo
1379  enddo
1380  enddo
1381  if ( jmh == ja-jmsh+1 ) then
1382  do k = 1, file_history_cartesc_model_nlayer
1383  do i = 1, im
1384  axis(i,jmh,k) = atmos_grid_cartesc_real_cz(k+ks-1,ims+i-1,jmsh+jmh-1)
1385  enddo
1386  enddo
1387  endif
1388  axis_name(1:3) = (/'x ', 'yh', 'z '/)
1389  call file_history_set_associatedcoordinate( 'height_xvz', 'height above ground level (half level xvz)', &
1390  'm', axis_name(1:3), axis(1:im,1:jmh,1:), start=start(:,3) )
1391 
1392  do k = 1, file_history_cartesc_model_nlayer
1393  do j = 1, min(jmh,ja-jmsh)
1394  do i = 1, min(imh,ia-imsh)
1395  axis(i,j,k) = ( atmos_grid_cartesc_real_cz(k+ks-1,imsh+i-1,jmsh+j-1) + atmos_grid_cartesc_real_cz(k+ks-1,imsh+i ,jmsh+j-1) &
1396  + atmos_grid_cartesc_real_cz(k+ks-1,imsh+i-1,jmsh+j ) + atmos_grid_cartesc_real_cz(k+ks-1,imsh+i ,jmsh+j ) ) * 0.25_rp
1397  enddo
1398  enddo
1399  enddo
1400  if ( jmh == ja-jmsh+1 ) then
1401  do k = 1, file_history_cartesc_model_nlayer
1402  do i = 1, min(imh,ia-imsh)
1403  axis(i,jmh,k) = ( atmos_grid_cartesc_real_cz(k+ks-1,imsh+i-1,jmsh+jmh-1) + atmos_grid_cartesc_real_cz(k+ks-1,imsh+i,jmsh+jmh-1) ) * 0.5_rp
1404  enddo
1405  enddo
1406  endif
1407  if ( imh == ia-imsh+1 ) then
1408  do k = 1, file_history_cartesc_model_nlayer
1409  do j = 1, min(jmh,ja-jmsh)
1410  axis(imh,j,k) = ( atmos_grid_cartesc_real_cz(k+ks-1,imsh+imh-1,jmsh+j-1) + atmos_grid_cartesc_real_cz(k+ks-1,imsh+imh-1,jmsh+j) ) * 0.5_rp
1411  enddo
1412  enddo
1413  endif
1414  if ( imh == ia-imsh+1 .AND. jmh == ja-jmsh+1 ) then
1415  do k = 1, file_history_cartesc_model_nlayer
1416  axis(imh,jmh,k) = atmos_grid_cartesc_real_cz(k+ks-1,imsh+imh-1,jmsh+jmh-1)
1417  enddo
1418  endif
1419  axis_name(1:3) = (/'xh', 'yh', 'z '/)
1420  call file_history_set_associatedcoordinate( 'height_uvz', 'height above ground level (half level uvz)', &
1421  'm', axis_name(1:3), axis(1:imh,1:jmh,1:), start=start(:,4) )
1422 
1423  do k = 0, file_history_cartesc_model_nlayer
1424  do j = 1, jm
1425  do i = 1, min(imh,ia-imsh)
1426  axis(i,j,k) = ( atmos_grid_cartesc_real_fz(k+ks-1,imsh+i-1,jms+j-1) + atmos_grid_cartesc_real_fz(k+ks-1,imsh+i,jms+j-1) ) * 0.5_rp
1427  enddo
1428  enddo
1429  enddo
1430  if ( imh == ia-imsh+1 ) then
1431  do k = 0, file_history_cartesc_model_nlayer
1432  do j = 1, jm
1433  axis(imh,j,k) = atmos_grid_cartesc_real_fz(k+ks-1,imsh+imh-1,jms+j-1)
1434  enddo
1435  enddo
1436  endif
1437  axis_name(1:3) = (/'xh', 'y ', 'zh'/)
1438  call file_history_set_associatedcoordinate( 'height_uyw', 'height above ground level (half level uyw)', &
1439  'm', axis_name(1:3), axis(1:imh,1:jm,0:), start=start(:,2) )
1440 
1441  do k = 0, file_history_cartesc_model_nlayer
1442  do j = 1, min(jmh,ja-jmsh)
1443  do i = 1, im
1444  axis(i,j,k) = ( atmos_grid_cartesc_real_fz(k+ks-1,ims+i-1,jmsh+j-1) + atmos_grid_cartesc_real_fz(k+ks-1,ims+i-1,jmsh+j) ) * 0.5_rp
1445  enddo
1446  enddo
1447  enddo
1448  if ( jmh == ja-jmsh+1 ) then
1449  do k = 0, file_history_cartesc_model_nlayer
1450  do i = 1, im
1451  axis(i,jmh,k) = atmos_grid_cartesc_real_fz(k+ks-1,ims+i-1,jmsh+jmh-1)
1452  enddo
1453  enddo
1454  endif
1455  axis_name(1:3) = (/'x ', 'yh', 'zh'/)
1456  call file_history_set_associatedcoordinate( 'height_xvw', 'height above ground level (half level xvw)', &
1457  'm', axis_name(1:3), axis(1:im,1:jmh,0:), start=start(:,3) )
1458 
1459  do k = 0, file_history_cartesc_model_nlayer
1460  do j = 1, min(jmh,ja-jmsh)
1461  do i = 1, min(imh,ia-imsh)
1462  axis(i,j,k) = ( atmos_grid_cartesc_real_fz(k+ks-1,imsh+i-1,jmsh+j-1) + atmos_grid_cartesc_real_fz(k+ks-1,imsh+i ,jmsh+j-1) &
1463  + atmos_grid_cartesc_real_fz(k+ks-1,imsh+i-1,jmsh+j ) + atmos_grid_cartesc_real_fz(k+ks-1,imsh+i ,jmsh+j ) ) * 0.25_rp
1464  enddo
1465  enddo
1466  enddo
1467  if ( jmh == ja-jmsh+1 ) then
1468  do k = 0, file_history_cartesc_model_nlayer
1469  do i = 1, min(imh,ia-imsh)
1470  axis(i,jmh,k) = ( atmos_grid_cartesc_real_fz(k+ks-1,imsh+i-1,jmsh+jmh-1) + atmos_grid_cartesc_real_fz(k+ks-1,imsh+i,jmsh+jmh-1) ) * 0.5_rp
1471  enddo
1472  enddo
1473  endif
1474  if ( imh == ia-imsh+1 ) then
1475  do k = 0, file_history_cartesc_model_nlayer
1476  do j = 1, min(jmh,ja-jmsh)
1477  axis(imh,j,k) = ( atmos_grid_cartesc_real_fz(k+ks-1,imsh+imh-1,jmsh+j-1) + atmos_grid_cartesc_real_fz(k+ks-1,imsh+imh-1,jmsh+j) ) * 0.5_rp
1478  enddo
1479  enddo
1480  endif
1481  if ( imh == ia-imsh+1 .AND. jm == ja-jmsh+1 ) then
1482  do k = 0, file_history_cartesc_model_nlayer
1483  axis(imh,jmh,k) = atmos_grid_cartesc_real_fz(k+ks-1,imsh+imh-1,jmsh+jmh-1)
1484  enddo
1485  endif
1486  axis_name(1:3) = (/'xh', 'yh', 'zh'/)
1487  call file_history_set_associatedcoordinate( 'height_uvw', 'height above ground level (half level uvw)', 'm', &
1488  axis_name(1:3), axis(1:imh,1:jmh,0:), start=start(:,4) )
1489 
1490  axis(1:im,1:jm,1) = atmos_grid_cartesc_real_lon(ims:ime,jms:jme) / d2r
1491  axis_name(1:2) = (/'x ', 'y '/)
1492  call file_history_set_associatedcoordinate( 'lon', 'longitude', 'degrees_east', &
1493  axis_name(1:2), axis(1:im,1:jm,1), start=start(:,1) )
1494 
1495  axis(1:imh,1:jm,1) = atmos_grid_cartesc_real_lonuy(imsh:ime,jms:jme) / d2r
1496  axis_name(1:2) = (/'xh', 'y '/)
1497  call file_history_set_associatedcoordinate( 'lon_uy', 'longitude (half level uy)', 'degrees_east', &
1498  axis_name(1:2), axis(1:imh,1:jm,1), start=start(:,2) )
1499 
1500  axis(1:im,1:jmh,1) = atmos_grid_cartesc_real_lonxv(ims:ime,jmsh:jme) / d2r
1501  axis_name(1:2) = (/'x ', 'yh'/)
1502  call file_history_set_associatedcoordinate( 'lon_xv', 'longitude (half level xv)', 'degrees_east', &
1503  axis_name(1:2), axis(1:im,1:jmh,1), start=start(:,3) )
1504 
1505  axis(1:imh,1:jmh,1) = atmos_grid_cartesc_real_lonuv(imsh:ime,jmsh:jme) / d2r
1506  axis_name(1:2) = (/'xh', 'yh'/)
1507  call file_history_set_associatedcoordinate( 'lon_uv', 'longitude (half level uv)', 'degrees_east', &
1508  axis_name(1:2), axis(1:imh,1:jmh,1), start=start(:,4) )
1509 
1510  axis(1:im,1:jm,1) = atmos_grid_cartesc_real_lat(ims:ime,jms:jme) / d2r
1511  axis_name(1:2) = (/'x ', 'y '/)
1512  call file_history_set_associatedcoordinate( 'lat', 'latitude', 'degrees_north', &
1513  axis_name(1:2), axis(1:im,1:jm,1), start=start(:,1) )
1514 
1515  axis(1:imh,1:jm,1) = atmos_grid_cartesc_real_latuy(imsh:ime,jms:jme) / d2r
1516  axis_name(1:2) = (/'xh', 'y '/)
1517  call file_history_set_associatedcoordinate( 'lat_uy', 'latitude (half level uy)', 'degrees_north', &
1518  axis_name(1:2), axis(1:imh,1:jm,1), start=start(:,2) )
1519 
1520  axis(1:im,1:jmh,1) = atmos_grid_cartesc_real_latxv(ims:ime,jmsh:jme) / d2r
1521  axis_name(1:2) = (/'x ', 'yh'/)
1522  call file_history_set_associatedcoordinate( 'lat_xv', 'latitude (half level xv)', 'degrees_north', &
1523  axis_name(1:2), axis(1:im,1:jmh,1), start=start(:,3) )
1524 
1525  axis(1:imh,1:jmh,1) = atmos_grid_cartesc_real_latuv(imsh:ime,jmsh:jme) / d2r
1526  axis_name(1:2) = (/'xh', 'yh'/)
1527  call file_history_set_associatedcoordinate( 'lat_uv', 'latitude (half level uv)', 'degrees_north', &
1528  axis_name(1:2), axis(1:imh,1:jmh,1), start=start(:,4) )
1529 
1530  axis_name(1:2) = (/'x ', 'y '/)
1531  call file_history_set_associatedcoordinate( 'topo', 'topography', 'm', axis_name(1:2), &
1532  topography_zsfc(ims:ime,jms:jme), start=start(:,1) )
1533 
1534  axis_name(1:2) = (/'x ', 'y '/)
1535  call file_history_set_associatedcoordinate( 'lsmask', 'fraction for land-sea mask', '1', axis_name(1:2), &
1536  landuse_frac_land(ims:ime,jms:jme), start=start(:,1) )
1537 
1538  axis_name(1:2) = (/'x ','y '/)
1539  call file_history_set_associatedcoordinate( 'cell_area', 'area of grid cell', 'm2', axis_name(1:2), &
1540  area(ims:ime,jms:jme), start=start(:,1) )
1541  axis_name(1:2) = (/'xh','y '/)
1542  call file_history_set_associatedcoordinate( 'cell_area_uy', 'area of grid cell (half level uy)', 'm2', axis_name(1:2), &
1543  areauy(imsh:ime,jms:jme), start=start(:,2) )
1544  axis_name(1:2) = (/'x ','yh'/)
1545  call file_history_set_associatedcoordinate( 'cell_area_xv', 'area of grid cell (half level xv)', 'm2', axis_name(1:2), &
1546  areaxv(ims:ime,jmsh:jme), start=start(:,3) )
1547 
1548  do k = 1, file_history_cartesc_model_nlayer
1549  do j = 1, jm
1550  do i = 1, imh
1551  axis(i,j,k) = areazuy_x(ks+k-1,imsh+i-1,jms+j-1)
1552  end do
1553  end do
1554  end do
1555  axis_name = (/'xh', 'y ', 'z '/)
1556  call file_history_set_associatedcoordinate( 'cell_area_uyz_x', 'area of grid cell face (half level uyz, normal x)', 'm2', &
1557  axis_name(1:3), axis(1:imh,1:jm,1:), start=start(:,2) )
1558  do k = 1, file_history_cartesc_model_nlayer
1559  do j = 1, jmh
1560  do i = 1, im
1561  axis(i,j,k) = areazxv_y(ks+k-1,ims+i-1,jmsh+j-1)
1562  end do
1563  end do
1564  end do
1565  axis_name = (/'x ', 'yh', 'z '/)
1566  call file_history_set_associatedcoordinate( 'cell_area_xvz_y', 'area of grid cell face (half level xvz, normal y)', 'm2', &
1567  axis_name(1:3), axis(1:im,1:jmh,1:), start=start(:,3) )
1568  do k = 0, file_history_cartesc_model_nlayer
1569  do j = 1, jmh
1570  do i = 1, imh
1571  axis(i,j,k) = areawuy_x(ks+k-1,imsh+i-1,jmsh+j-1)
1572  end do
1573  end do
1574  end do
1575  axis_name = (/'xh', 'y ', 'zh'/)
1576  call file_history_set_associatedcoordinate( 'cell_area_uyw_x', 'area of grid cell face (half level uyw, normal x)', 'm2', &
1577  axis_name(1:3), axis(1:imh,1:jmh,0:), start=start(:,2) )
1578  do k = 0, file_history_cartesc_model_nlayer
1579  do j = 1, jmh
1580  do i = 1, im
1581  axis(i,j,k) = areawxv_y(ks+k-1,ims+i-1,jmsh+j-1)
1582  end do
1583  end do
1584  end do
1585  axis_name = (/'x ', 'yh', 'zh'/)
1586  call file_history_set_associatedcoordinate( 'cell_area_xvw_y', 'area of grid cell face (half level xvw, normal y)', 'm2', &
1587  axis_name(1:3), axis(1:im,1:jmh,0:), start=start(:,3) )
1588  do k = 1, file_history_cartesc_model_nlayer
1589  do j = 1, jm
1590  do i = 1, im
1591  axis(i,j,k) = areazxy_x(ks+k-1,ims+i-1,jms+j-1)
1592  end do
1593  end do
1594  end do
1595  axis_name = (/'x ', 'y ', 'z '/)
1596  call file_history_set_associatedcoordinate( 'cell_area_xyz_x', 'area of grid cell face (half level xyz, normal x)', 'm2', &
1597  axis_name(1:3), axis(1:im,1:jm,1:), start=start(:,1) )
1598  do k = 1, file_history_cartesc_model_nlayer
1599  do j = 1, jmh
1600  do i = 1, imh
1601  axis(i,j,k) = areazuv_y(ks+k-1,imsh+i-1,jmsh+j-1)
1602  end do
1603  end do
1604  end do
1605  axis_name = (/'xh', 'yh', 'z '/)
1606  call file_history_set_associatedcoordinate( 'cell_area_uvz_y', 'area of grid cell face (half level uvz, normal y)', 'm2', &
1607  axis_name(1:3), axis(1:imh,1:jmh,1:), start=start(:,4) )
1608  do k = 1, file_history_cartesc_model_nlayer
1609  do j = 1, jmh
1610  do i = 1, imh
1611  axis(i,j,k) = areazuv_x(ks+k-1,imsh+i-1,jmsh+j-1)
1612  end do
1613  end do
1614  end do
1615  axis_name = (/'xh', 'yh', 'z '/)
1616  call file_history_set_associatedcoordinate( 'cell_area_uvz_x', 'area of grid cell face (half level uvz, normal x)', 'm2', &
1617  axis_name(1:3), axis(1:imh,1:jmh,1:), start=start(:,4) )
1618  do k = 1, file_history_cartesc_model_nlayer
1619  do j = 1, jm
1620  do i = 1, im
1621  axis(i,j,k) = areazxy_y(ks+k-1,ims+i-1,jms+j-1)
1622  end do
1623  end do
1624  end do
1625  axis_name = (/'x ', 'y ', 'z '/)
1626  call file_history_set_associatedcoordinate( 'cell_area_xyz_y', 'area of grid cell face (half level xyz, normal y)', 'm2', &
1627  axis_name(1:3), axis(1:im,1:jm,1:), start=start(:,1) )
1628 
1629  do k = 1, file_history_cartesc_model_nlayer
1630  do j = 1, jm
1631  do i = 1, im
1632  axis(i,j,k) = vol(ks+k-1,ims+i-1,jms+j-1)
1633  end do
1634  end do
1635  end do
1636  axis_name = (/ 'x ', 'y ', 'z '/)
1637  call file_history_set_associatedcoordinate( 'cell_volume', 'volume of grid cell', 'm3', &
1638  axis_name(1:3), axis(1:im,1:jm,1:), start=start(:,1) )
1639  do k = 0, file_history_cartesc_model_nlayer
1640  do j = 1, jm
1641  do i = 1, im
1642  axis(i,j,k) = volwxy(ks+k-1,ims+i-1,jms+j-1)
1643  end do
1644  end do
1645  end do
1646  axis_name = (/'x ', 'y ', 'zh'/)
1647  call file_history_set_associatedcoordinate( 'cell_volume_xyw', 'volume of grid cell (half level xyw)', 'm3', &
1648  axis_name(1:3), axis(1:im,1:jm,0:), start=start(:,1) )
1649  do k = 1, file_history_cartesc_model_nlayer
1650  do j = 1, jm
1651  do i = 1, imh
1652  axis(i,j,k) = volzuy(ks+k-1,imsh+i-1,jms+j-1)
1653  end do
1654  end do
1655  end do
1656  axis_name = (/'xh', 'y ', 'z '/)
1657  call file_history_set_associatedcoordinate( 'cell_volume_uyz', 'volume of grid cell (half level uyz)', 'm3', &
1658  axis_name(1:3), axis(1:imh,1:jm,1:), start=start(:,2) )
1659  do k = 1, file_history_cartesc_model_nlayer
1660  do j = 1, jmh
1661  do i = 1, im
1662  axis(i,j,k) = volzxv(ks+k-1,ims+i-1,jmsh+j-1)
1663  end do
1664  end do
1665  end do
1666  axis_name = (/'x ', 'yh', 'z '/)
1667  call file_history_set_associatedcoordinate( 'cell_volume_xvz', 'volume of grid cell (half level xvz)', 'm3', &
1668  axis_name(1:3), axis(1:im,1:jmh,1:), start=start(:,3) )
1669 
1670  if ( okmax > 0 ) then
1671  do k = 1, okmax
1672  do j = 1, jm
1673  do i = 1, im
1674  axiso(i,j,k) = volo(oks+k-1,ims+i-1,jms+j-1)
1675  end do
1676  end do
1677  end do
1678  axis_name = (/'x ', 'y ', 'oz'/)
1679  call file_history_set_associatedcoordinate( 'cell_volume_xyo', 'volume of grid cell', 'm3', &
1680  axis_name(1:3), axiso(:,:,:), start=start(:,1) )
1681  endif
1682 
1683  if ( lkmax > 0 ) then
1684  do k = 1, lkmax
1685  do j = 1, jm
1686  do i = 1, im
1687  axisl(i,j,k) = voll(lks+k-1,ims+i-1,jms+j-1)
1688  end do
1689  end do
1690  end do
1691  axis_name = (/'x ', 'y ', 'lz'/)
1692  call file_history_set_associatedcoordinate( 'cell_volume_xyl', 'volume of grid cell', 'm3', &
1693  axis_name(1:3), axisl(:,:,:), start=start(:,1) )
1694  endif
1695 
1696  if ( ukmax > 0 ) then
1697  do k = 1, ukmax
1698  do j = 1, jm
1699  do i = 1, im
1700  axisu(i,j,k) = volu(uks+k-1,ims+i-1,jms+j-1)
1701  end do
1702  end do
1703  end do
1704  axis_name = (/'x ', 'y ', 'uz'/)
1705  call file_history_set_associatedcoordinate( 'cell_volume_xyu', 'volume of grid cell', 'm3', &
1706  axis_name(1:3), axisu(:,:,:), start=start(:,1) )
1707  endif
1708 
1709  return
1710  end subroutine file_history_cartesc_set_axes
1711 
1712  !-----------------------------------------------------------------------------
1713  subroutine file_history_cartesc_set_axes_attributes
1714  use scale_atmos_grid_cartesc, only: &
1716  use scale_calendar, only: &
1718  use scale_file_history, only: &
1720  file_history_set_attribute
1721  use scale_prc, only: &
1722  prc_myrank
1723  use scale_const, only: &
1724  undef => const_undef
1725  use scale_prc_cartesc, only: &
1726  prc_2drank, &
1727  prc_num_x, &
1728  prc_num_y, &
1729  prc_periodic_x, &
1730  prc_periodic_y, &
1731  prc_has_w, &
1732  prc_has_e, &
1733  prc_has_s, &
1734  prc_has_n
1735  use scale_file, only: &
1737  use scale_file_cartesc, only: &
1738  axisattinfo
1739  use scale_mapprojection, only: &
1740  minfo => mapprojection_mappinginfo
1741  implicit none
1742 
1743  character(len=34) :: tunits
1744  character(len=H_SHORT) :: calendar
1745 
1746  type(axisattinfo) :: ainfo(4) ! x, xh, y, yh
1747  !---------------------------------------------------------------------------
1748 
1749  call file_history_set_attribute( "global", "Conventions", "CF-1.6" ) ! [IN]
1750 
1751  call file_history_set_attribute( "global", "grid_name", atmos_grid_cartesc_name ) ! [IN]
1752 
1753  if ( file_history_aggregate ) then
1754  call file_history_set_attribute( "global", "scale_cartesC_prc_rank_x", (/0/) ) ! [IN]
1755  call file_history_set_attribute( "global", "scale_cartesC_prc_rank_y", (/0/) ) ! [IN]
1756 
1757  call file_history_set_attribute( "global", "scale_cartesC_prc_num_x", (/1/) ) ! [IN]
1758  call file_history_set_attribute( "global", "scale_cartesC_prc_num_y", (/1/) ) ! [IN]
1759  else
1760  call file_history_set_attribute( "global", "scale_cartesC_prc_rank_x", (/prc_2drank(prc_myrank,1)/) ) ! [IN]
1761  call file_history_set_attribute( "global", "scale_cartesC_prc_rank_y", (/prc_2drank(prc_myrank,2)/) ) ! [IN]
1762 
1763  call file_history_set_attribute( "global", "scale_cartesC_prc_num_x", (/prc_num_x/) ) ! [IN]
1764  call file_history_set_attribute( "global", "scale_cartesC_prc_num_y", (/prc_num_y/) ) ! [IN]
1765  endif
1766 
1767  call file_history_set_attribute( "global", "scale_cartesC_prc_periodic_z", .false. ) ! [IN]
1768  call file_history_set_attribute( "global", "scale_cartesC_prc_periodic_x", prc_periodic_x ) ! [IN]
1769  call file_history_set_attribute( "global", "scale_cartesC_prc_periodic_y", prc_periodic_y ) ! [IN]
1770 
1771  call file_history_set_attribute( "global", "scale_atmos_grid_cartesC_index_imaxg", (/imaxg/) ) ! [IN]
1772  call file_history_set_attribute( "global", "scale_atmos_grid_cartesC_index_jmaxg", (/jmaxg/) ) ! [IN]
1773 
1774  call file_history_set_attribute( "global", "scale_atmos_grid_cartesC_index_kmax", (/kmax/) ) ! [IN]
1775  if ( okmax > 0 ) call file_history_set_attribute( "global", "scale_ocean_grid_cartesC_index_kmax", (/okmax/) ) ! [IN]
1776  if ( lkmax > 0 ) call file_history_set_attribute( "global", "scale_land_grid_cartesC_index_kmax", (/lkmax/) ) ! [IN]
1777  if ( ukmax > 0 ) call file_history_set_attribute( "global", "scale_urban_grid_cartesC_index_kmax", (/ukmax/) ) ! [IN]
1778 
1779  call file_history_set_attribute( "global", "scale_atmos_grid_cartesC_index_khalo", (/khalo/) ) ! [IN]
1780  call file_history_set_attribute( "global", "scale_atmos_grid_cartesC_index_ihalo", (/ihalo/) ) ! [IN]
1781  call file_history_set_attribute( "global", "scale_atmos_grid_cartesC_index_jhalo", (/jhalo/) ) ! [IN]
1782 
1783  call calendar_get_name( calendar )
1784  if ( calendar /= "" ) call file_history_set_attribute( "global", "calendar", calendar )
1785  call file_get_cftunits( file_history_cartesc_startdate(:), tunits )
1786  call file_history_set_attribute( "global", "time_units", tunits )
1787  call file_history_set_attribute( "global", "time_start", (/file_history_cartesc_startsubsec/) )
1788 
1789  if ( prc_periodic_x ) then
1790  ainfo(1)%periodic = .true.
1791  ainfo(2)%periodic = .true.
1792  else
1793  ainfo(1)%periodic = .false.
1794  ainfo(2)%periodic = .false.
1795  endif
1796 
1797  if ( prc_periodic_y ) then
1798  ainfo(3)%periodic = .true.
1799  ainfo(4)%periodic = .true.
1800  else
1801  ainfo(3)%periodic = .false.
1802  ainfo(4)%periodic = .false.
1803  endif
1804 
1805  ! for x
1806  if ( prc_periodic_x .OR. .NOT. file_history_cartesc_boundary ) then
1807  ainfo(1)%size_global (1) = imax * prc_num_x
1808  ainfo(1)%start_global(1) = is_ing - ihalo
1809  ainfo(1)%halo_global (1) = 0 ! west side
1810  ainfo(1)%halo_global (2) = 0 ! east side
1811  ainfo(1)%halo_local (1) = 0 ! west side
1812  ainfo(1)%halo_local (2) = 0 ! east side
1813  else
1814  ainfo(1)%size_global (1) = iag
1815  ainfo(1)%start_global(1) = isga
1816  ainfo(1)%halo_global (1) = ihalo ! west side
1817  ainfo(1)%halo_global (2) = ihalo ! east side
1818  ainfo(1)%halo_local (1) = ihalo ! west side
1819  ainfo(1)%halo_local (2) = ihalo ! east side
1820  if ( .not. file_history_aggregate ) then
1821  if( prc_has_w ) ainfo(1)%halo_local(1) = 0
1822  if( prc_has_e ) ainfo(1)%halo_local(2) = 0
1823  end if
1824  endif
1825 
1826  ! for xh
1827  ainfo(2) = ainfo(1)
1828  if ( .NOT. prc_periodic_x .AND. .NOT. file_history_cartesc_boundary ) then
1829  ainfo(2)%size_global (1) = ainfo(2)%size_global (1) + 1
1830  ainfo(2)%halo_global (1) = ainfo(2)%halo_global (1) + 1
1831  if ( prc_has_w .and. ( .not. file_history_aggregate ) ) then
1832  ainfo(2)%start_global(1) = ainfo(2)%start_global(1) + 1
1833  else
1834  ainfo(2)%halo_local (1) = ainfo(2)%halo_local (1) + 1
1835  endif
1836  endif
1837 
1838  ! for y
1839  if ( prc_periodic_y .OR. .NOT. file_history_cartesc_boundary ) then
1840  ainfo(3)%size_global (1) = jmax * prc_num_y
1841  ainfo(3)%start_global(1) = js_ing - jhalo
1842  ainfo(3)%halo_global (1) = 0 ! south side
1843  ainfo(3)%halo_global (2) = 0 ! north side
1844  ainfo(3)%halo_local (1) = 0 ! south side
1845  ainfo(3)%halo_local (2) = 0 ! north side
1846  else
1847  ainfo(3)%size_global (1) = jag
1848  ainfo(3)%start_global(1) = jsga
1849  ainfo(3)%halo_global (1) = jhalo ! south side
1850  ainfo(3)%halo_global (2) = jhalo ! north side
1851  ainfo(3)%halo_local (1) = jhalo ! south side
1852  ainfo(3)%halo_local (2) = jhalo ! north side
1853  if ( .not. file_history_aggregate ) then
1854  if( prc_has_s ) ainfo(3)%halo_local(1) = 0
1855  if( prc_has_n ) ainfo(3)%halo_local(2) = 0
1856  end if
1857  endif
1858 
1859  ! for yh
1860  ainfo(4) = ainfo(3)
1861  if ( .NOT. prc_periodic_y .AND. .NOT. file_history_cartesc_boundary ) then
1862  ainfo(4)%size_global (1) = ainfo(4)%size_global (1) + 1
1863  ainfo(4)%halo_global (1) = ainfo(4)%halo_global (1) + 1
1864  if ( prc_has_s .and. ( .not. file_history_aggregate ) ) then
1865  ainfo(4)%start_global(1) = ainfo(4)%start_global(1) + 1
1866  else
1867  ainfo(4)%halo_local (1) = ainfo(4)%halo_local (1) + 1
1868  endif
1869  endif
1870 
1871  if ( file_history_aggregate ) then
1872  ainfo(1)%start_global(1) = 1
1873  ainfo(2)%start_global(1) = 1
1874  ainfo(3)%start_global(1) = 1
1875  ainfo(4)%start_global(1) = 1
1876  end if
1877 
1878  call file_history_set_attribute( "x" , "size_global" , ainfo(1)%size_global (:) )
1879  call file_history_set_attribute( "x" , "start_global", ainfo(1)%start_global(:) )
1880  call file_history_set_attribute( "x" , "halo_global" , ainfo(1)%halo_global (:) )
1881  call file_history_set_attribute( "x" , "halo_local" , ainfo(1)%halo_local (:) )
1882  call file_history_set_attribute( "x" , "periodic" , ainfo(1)%periodic )
1883 
1884  call file_history_set_attribute( "xh", "size_global" , ainfo(2)%size_global (:) )
1885  call file_history_set_attribute( "xh", "start_global", ainfo(2)%start_global(:) )
1886  call file_history_set_attribute( "xh", "halo_global" , ainfo(2)%halo_global (:) )
1887  call file_history_set_attribute( "xh", "halo_local" , ainfo(2)%halo_local (:) )
1888  call file_history_set_attribute( "xh", "periodic" , ainfo(2)%periodic )
1889 
1890  call file_history_set_attribute( "y" , "size_global" , ainfo(3)%size_global (:) )
1891  call file_history_set_attribute( "y" , "start_global", ainfo(3)%start_global(:) )
1892  call file_history_set_attribute( "y" , "halo_global" , ainfo(3)%halo_global (:) )
1893  call file_history_set_attribute( "y" , "halo_local" , ainfo(3)%halo_local (:) )
1894  call file_history_set_attribute( "y" , "periodic" , ainfo(3)%periodic )
1895 
1896  call file_history_set_attribute( "yh", "size_global" , ainfo(4)%size_global (:) )
1897  call file_history_set_attribute( "yh", "start_global", ainfo(4)%start_global(:) )
1898  call file_history_set_attribute( "yh", "halo_global" , ainfo(4)%halo_global (:) )
1899  call file_history_set_attribute( "yh", "halo_local" , ainfo(4)%halo_local (:) )
1900  call file_history_set_attribute( "yh", "periodic" , ainfo(4)%periodic )
1901 
1902  ! map projection info
1903 
1904  if ( minfo%mapping_name /= "" ) then
1905  call file_history_set_attribute( "x" , "standard_name", "projection_x_coordinate" )
1906  call file_history_set_attribute( "xh", "standard_name", "projection_x_coordinate" )
1907  call file_history_set_attribute( "y" , "standard_name", "projection_y_coordinate" )
1908  call file_history_set_attribute( "yh", "standard_name", "projection_y_coordinate" )
1909 
1910  call file_history_set_attribute( minfo%mapping_name, "grid_mapping_name", minfo%mapping_name, add_variable=.true. )
1911 
1912  if ( minfo%false_easting /= undef ) then
1913  call file_history_set_attribute( minfo%mapping_name, & ! [IN]
1914  "false_easting", & ! [IN]
1915  minfo%false_easting ) ! [IN]
1916  endif
1917 
1918  if ( minfo%false_northing /= undef ) then
1919  call file_history_set_attribute( minfo%mapping_name, & ! [IN]
1920  "false_northing", & ! [IN]
1921  minfo%false_northing ) ! [IN]
1922  endif
1923 
1924  if ( minfo%longitude_of_central_meridian /= undef ) then
1925  call file_history_set_attribute( minfo%mapping_name, & ! [IN]
1926  "longitude_of_central_meridian", & ! [IN]
1927  minfo%longitude_of_central_meridian ) ! [IN]
1928  endif
1929 
1930  if ( minfo%longitude_of_projection_origin /= undef ) then
1931  call file_history_set_attribute( minfo%mapping_name, & ! [IN]
1932  "longitude_of_projection_origin", & ! [IN]
1933  minfo%longitude_of_projection_origin ) ! [IN]
1934  endif
1935 
1936  if ( minfo%latitude_of_projection_origin /= undef ) then
1937  call file_history_set_attribute( minfo%mapping_name, & ! [IN]
1938  "latitude_of_projection_origin", & ! [IN]
1939  minfo%latitude_of_projection_origin ) ! [IN]
1940  endif
1941 
1942  if ( minfo%straight_vertical_longitude_from_pole /= undef ) then
1943  call file_history_set_attribute( minfo%mapping_name, & ! [IN]
1944  "straight_vertical_longitude_from_pole", & ! [IN]
1945  minfo%straight_vertical_longitude_from_pole ) ! [IN]
1946  endif
1947 
1948  if ( minfo%standard_parallel(1) /= undef ) then
1949  if ( minfo%standard_parallel(2) /= undef ) then
1950  call file_history_set_attribute( minfo%mapping_name, & ! [IN]
1951  "standard_parallel", & ! [IN]
1952  minfo%standard_parallel(:) ) ! [IN]
1953  else
1954  call file_history_set_attribute( minfo%mapping_name, & ! [IN]
1955  "standard_parallel", & ! [IN]
1956  minfo%standard_parallel(1) ) ! [IN]
1957  endif
1958  endif
1959 
1960  if ( minfo%rotation /= undef ) then
1961  call file_history_set_attribute( minfo%mapping_name, & ! [IN]
1962  "rotation", & ! [IN]
1963  minfo%rotation ) ! [IN]
1964  endif
1965 
1966  endif
1967 
1968  ! area and volume
1969  call file_history_set_attribute( "cell_area", "standard_name", "area" ) ! [IN]
1970  call file_history_set_attribute( "cell_area_uy", "standard_name", "area" ) ! [IN]
1971  call file_history_set_attribute( "cell_area_xv", "standard_name", "area" ) ! [IN]
1972 
1973  call file_history_set_attribute( "cell_area_uyz_x", "standard_name", "area" ) ! [IN]
1974  call file_history_set_attribute( "cell_area_xvz_y", "standard_name", "area" ) ! [IN]
1975  call file_history_set_attribute( "cell_area_uyw_x", "standard_name", "area" ) ! [IN]
1976  call file_history_set_attribute( "cell_area_xvw_y", "standard_name", "area" ) ! [IN]
1977  call file_history_set_attribute( "cell_area_xyz_x", "standard_name", "area" ) ! [IN]
1978  call file_history_set_attribute( "cell_area_uvz_y", "standard_name", "area" ) ! [IN]
1979  call file_history_set_attribute( "cell_area_uvz_x", "standard_name", "area" ) ! [IN]
1980  call file_history_set_attribute( "cell_area_xyz_y", "standard_name", "area" ) ! [IN]
1981 
1982  call file_history_set_attribute( "cell_volume", "standard_name", "volume" ) ! [IN]
1983  call file_history_set_attribute( "cell_volume_xyw", "standard_name", "volume" ) ! [IN]
1984  call file_history_set_attribute( "cell_volume_uyz", "standard_name", "volume" ) ! [IN]
1985  call file_history_set_attribute( "cell_volume_xvz", "standard_name", "volume" ) ! [IN]
1986 
1987  if ( okmax > 0 ) then
1988  call file_history_set_attribute( "cell_volume_xyo", "standard_name", "volume" ) ! [IN]
1989  endif
1990  if ( lkmax > 0 ) then
1991  call file_history_set_attribute( "cell_volume_xyl", "standard_name", "volume" ) ! [IN]
1992  endif
1993  if ( ukmax > 0 ) then
1994  call file_history_set_attribute( "cell_volume_xyu", "standard_name", "volume" ) ! [IN]
1995  endif
1996 
1997  ! SGRID
1998  call file_history_set_attribute( "grid", "cf_role", "grid_topology", add_variable=.true. )
1999  call file_history_set_attribute( "grid", "topology_dimension", (/ 2 /) )
2000  call file_history_set_attribute( "grid", "node_dimensions", "xh yh" )
2001  call file_history_set_attribute( "grid", "face_dimensions", "x: xh (padding: none) y: yh (padding: none)" )
2002  call file_history_set_attribute( "grid", "node_coordinates", "lon_uv lat_uv" )
2003  call file_history_set_attribute( "grid", "face_coordinates", "lon lat" )
2004  call file_history_set_attribute( "grid", "edge1_coordinates", "lon_uy lat_uy" )
2005  call file_history_set_attribute( "grid", "edge2_coordinates", "lon_xv lat_xv" )
2006  call file_history_set_attribute( "grid", "vertical_dimensions", "z: zh (padding: none)" )
2007 
2008  call file_history_set_attribute( "grid_ocean", "cf_role", "grid_topology", add_variable=.true. )
2009  call file_history_set_attribute( "grid_ocean", "topology_dimension", (/ 2 /) )
2010  call file_history_set_attribute( "grid_ocean", "node_dimensions", "xh yh" )
2011  call file_history_set_attribute( "grid_ocean", "face_dimensions", "x: xh (padding: none) y: yh (padding: none)" )
2012  call file_history_set_attribute( "grid_ocean", "node_coordinates", "lon_uv lat_uv" )
2013  call file_history_set_attribute( "grid_ocean", "face_coordinates", "lon lat" )
2014  call file_history_set_attribute( "grid_ocean", "edge1_coordinates", "lon_uy lat_uy" )
2015  call file_history_set_attribute( "grid_ocean", "edge2_coordinates", "lon_xv lat_xv" )
2016  call file_history_set_attribute( "grid_ocean", "vertical_dimensions", "oz: ozh (padding: none)" )
2017 
2018  call file_history_set_attribute( "grid_land", "cf_role", "grid_topology", add_variable=.true. )
2019  call file_history_set_attribute( "grid_land", "topology_dimension", (/ 2 /) )
2020  call file_history_set_attribute( "grid_land", "node_dimensions", "xh yh" )
2021  call file_history_set_attribute( "grid_land", "face_dimensions", "x: xh (padding: none) y: yh (padding: none)" )
2022  call file_history_set_attribute( "grid_land", "node_coordinates", "lon_uv lat_uv" )
2023  call file_history_set_attribute( "grid_land", "face_coordinates", "lon lat" )
2024  call file_history_set_attribute( "grid_land", "edge1_coordinates", "lon_uy lat_uy" )
2025  call file_history_set_attribute( "grid_land", "edge2_coordinates", "lon_xv lat_xv" )
2026  call file_history_set_attribute( "grid_land", "vertical_dimensions", "lz: lzh (padding: none)" )
2027 
2028  call file_history_set_attribute( "grid_urban", "cf_role", "grid_topology", add_variable=.true. )
2029  call file_history_set_attribute( "grid_urban", "topology_dimension", (/ 2 /) )
2030  call file_history_set_attribute( "grid_urban", "node_dimensions", "xh yh" )
2031  call file_history_set_attribute( "grid_urban", "face_dimensions", "x: xh (padding: none) y: yh (padding: none)" )
2032  call file_history_set_attribute( "grid_urban", "node_coordinates", "lon_uv lat_uv" )
2033  call file_history_set_attribute( "grid_urban", "face_coordinates", "lon lat" )
2034  call file_history_set_attribute( "grid_urban", "edge1_coordinates", "lon_uy lat_uy" )
2035  call file_history_set_attribute( "grid_urban", "edge2_coordinates", "lon_xv lat_xv" )
2036  call file_history_set_attribute( "grid_urban", "vertical_dimensions", "uz: uzh (padding: none)" )
2037 
2038  call file_history_set_attribute( "grid_pressure", "cf_role", "grid_topology", add_variable=.true. )
2039  call file_history_set_attribute( "grid_pressure", "topology_dimension", (/ 2 /) )
2040  call file_history_set_attribute( "grid_pressure", "node_dimensions", "xh yh" )
2041  call file_history_set_attribute( "grid_pressure", "face_dimensions", "x: xh (padding: none) y: yh (padding: none)" )
2042  call file_history_set_attribute( "grid_pressure", "node_coordinates", "lon_uv lat_uv" )
2043  call file_history_set_attribute( "grid_pressure", "face_coordinates", "lon lat" )
2044  call file_history_set_attribute( "grid_pressure", "edge1_coordinates", "lon_uy lat_uy" )
2045  call file_history_set_attribute( "grid_pressure", "edge2_coordinates", "lon_xv lat_xv" )
2046  call file_history_set_attribute( "grid_pressure", "vertical_dimensions", "pressure" )
2047 
2048  call file_history_set_attribute( "grid_z", "cf_role", "grid_topology", add_variable=.true. )
2049  call file_history_set_attribute( "grid_z", "topology_dimension", (/ 2 /) )
2050  call file_history_set_attribute( "grid_z", "node_dimensions", "xh yh" )
2051  call file_history_set_attribute( "grid_z", "face_dimensions", "x: xh (padding: none) y: yh (padding: none)" )
2052  call file_history_set_attribute( "grid_z", "node_coordinates", "lon_uv lat_uv" )
2053  call file_history_set_attribute( "grid_z", "face_coordinates", "lon lat" )
2054  call file_history_set_attribute( "grid_z", "edge1_coordinates", "lon_uy lat_uy" )
2055  call file_history_set_attribute( "grid_z", "edge2_coordinates", "lon_xv lat_xv" )
2056  call file_history_set_attribute( "grid_z", "vertical_dimensions", "height_xyw: height (padding: none)" )
2057 
2058  call file_history_set_attribute( "grid_model", "cf_role", "grid_topology", add_variable=.true. )
2059  call file_history_set_attribute( "grid_model", "topology_dimension", (/ 2 /) )
2060  call file_history_set_attribute( "grid_model", "node_dimensions", "FX FY" )
2061  call file_history_set_attribute( "grid_model", "face_dimensions", "CX: FY (padding: none) CY: FY (padding: none)" )
2062  call file_history_set_attribute( "grid_model", "vertical_dimensions", "CZ: FZ (padding: none)" )
2063 
2064  call file_history_set_attribute( "grid_model_global", "cf_role", "grid_topology", add_variable=.true. )
2065  call file_history_set_attribute( "grid_model_global", "topology_dimension", (/ 2 /) )
2066  call file_history_set_attribute( "grid_model_global", "node_dimensions", "FXG FYG" )
2067  call file_history_set_attribute( "grid_model_global", "face_dimensions", "CXG: FYG (padding: none) CYG: FYG (padding: none)" )
2068  call file_history_set_attribute( "grid_model_global", "vertical_dimensions", "CZ: FZ (padding: none)" )
2069 
2070  return
2071  end subroutine file_history_cartesc_set_axes_attributes
2072 
2073 end module scale_file_history_cartesc
scale_atmos_grid_cartesc_index::isb
integer, public isb
Definition: scale_atmos_grid_cartesC_index.F90:63
scale_atmos_grid_cartesc_index::jagb
integer, public jagb
Definition: scale_atmos_grid_cartesC_index.F90:76
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_areazxy_y
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_areazxy_y
virtical area (zxy, normal y) [m2]
Definition: scale_atmos_grid_cartesC_real.F90:75
scale_atmos_grid_cartesc_index::ke
integer, public ke
end point of inner domain: z, local
Definition: scale_atmos_grid_cartesC_index.F90:52
scale_urban_grid_cartesc_real::urban_grid_cartesc_real_vol
real(rp), dimension(:,:,:), allocatable, public urban_grid_cartesc_real_vol
volume of grid cell
Definition: scale_urban_grid_cartesC_real.F90:38
scale_prc::prc_abort
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:342
scale_ocean_grid_cartesc::ocean_grid_cartesc_cz
real(rp), dimension(:), allocatable, public ocean_grid_cartesc_cz
center coordinate [m]: z, local=global
Definition: scale_ocean_grid_cartesC.F90:35
scale_urban_grid_cartesc::urban_grid_cartesc_cdz
real(rp), dimension(:), allocatable, public urban_grid_cartesc_cdz
z-length of control volume [m]
Definition: scale_urban_grid_cartesC.F90:37
scale_atmos_grid_cartesc::atmos_grid_cartesc_cdz
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_cdz
z-length of control volume [m]
Definition: scale_atmos_grid_cartesC.F90:42
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_cz
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_cz
geopotential height [m] (zxy)
Definition: scale_atmos_grid_cartesC_real.F90:38
scale_ocean_grid_cartesc_index::oke
integer, public oke
Definition: scale_ocean_grid_cartesC_index.F90:38
scale_atmos_grid_cartesc_index::ihalo
integer, public ihalo
Definition: scale_atmos_grid_cartesC_index.F90:44
scale_atmos_grid_cartesc_index::jmaxb
integer, public jmaxb
Definition: scale_atmos_grid_cartesC_index.F90:62
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_volwxy
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_volwxy
control volume (wxy) [m3]
Definition: scale_atmos_grid_cartesC_real.F90:84
scale_prc_cartesc::prc_has_s
logical, public prc_has_s
Definition: scale_prc_cartesC.F90:50
scale_land_grid_cartesc::land_grid_cartesc_cdz
real(rp), dimension(:), allocatable, public land_grid_cartesc_cdz
z-length of control volume [m]
Definition: scale_land_grid_cartesC.F90:36
scale_calendar::calendar_get_name
subroutine, public calendar_get_name(name)
Definition: scale_calendar.F90:715
scale_prc_cartesc::prc_periodic_y
logical, public prc_periodic_y
periodic condition or not (Y)?
Definition: scale_prc_cartesC.F90:53
scale_precision
module PRECISION
Definition: scale_precision.F90:14
scale_atmos_grid_cartesc_index::ka
integer, public ka
Definition: scale_atmos_grid_cartesC_index.F90:47
scale_atmos_grid_cartesc_index::isga
integer, public isga
start point of the full domain: cx, global
Definition: scale_atmos_grid_cartesC_index.F90:81
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_lonxv
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_lonxv
longitude at staggered point (xv) [rad,0-2pi]
Definition: scale_atmos_grid_cartesC_real.F90:50
scale_atmos_grid_cartesc_index::isgb
integer, public isgb
start point of the inner domain: x, global
Definition: scale_atmos_grid_cartesC_index.F90:85
scale_file_history_cartesc::file_history_cartesc_set_pres
subroutine, public file_history_cartesc_set_pres(PRES, PRESH, SFC_PRES)
set hydrostatic pressure for pressure coordinate
Definition: scale_file_history_cartesC.F90:260
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_areawuy_x
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_areawuy_x
virtical area (wuy, normal x) [m2]
Definition: scale_atmos_grid_cartesC_real.F90:68
scale_land_grid_cartesc_index::lkmax
integer, public lkmax
Definition: scale_land_grid_cartesC_index.F90:32
scale_atmos_grid_cartesc::atmos_grid_cartesc_cdyg
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_cdyg
center coordinate [m]: y, global
Definition: scale_atmos_grid_cartesC.F90:81
scale_time::time_nowsubsec
real(dp), public time_nowsubsec
subsecond part of current time [sec]
Definition: scale_time.F90:69
scale_prof::prof_rapstart
subroutine, public prof_rapstart(rapname_base, level, disable_barrier)
Start raptime.
Definition: scale_prof.F90:159
scale_atmos_grid_cartesc::atmos_grid_cartesc_fyg
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_fyg
face coordinate [m]: y, global
Definition: scale_atmos_grid_cartesC.F90:78
scale_atmos_grid_cartesc::atmos_grid_cartesc_cdy
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_cdy
y-length of control volume [m]
Definition: scale_atmos_grid_cartesC.F90:61
scale_atmos_grid_cartesc::atmos_grid_cartesc_fx
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_fx
face coordinate [m]: x, local
Definition: scale_atmos_grid_cartesC.F90:57
scale_atmos_grid_cartesc::atmos_grid_cartesc_fdyg
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_fdyg
center coordinate [m]: y, global
Definition: scale_atmos_grid_cartesC.F90:83
scale_file_h::file_rmiss
real(dp), parameter, public file_rmiss
Definition: scale_file_h.F90:49
scale_prc::prc_myrank
integer, public prc_myrank
process num in local communicator
Definition: scale_prc.F90:90
scale_atmos_grid_cartesc::atmos_grid_cartesc_fdy
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_fdy
y-length of grid(j+1) to grid(j) [m]
Definition: scale_atmos_grid_cartesC.F90:63
scale_ocean_grid_cartesc_index::okmax
integer, public okmax
Definition: scale_ocean_grid_cartesC_index.F90:32
scale_atmos_grid_cartesc_index::khalo
integer, parameter, public khalo
Definition: scale_atmos_grid_cartesC_index.F90:43
scale_topography
module TOPOGRAPHY
Definition: scale_topography.F90:11
scale_atmos_grid_cartesc::atmos_grid_cartesc_fbfx
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_fbfx
face buffer factor (0-1): x
Definition: scale_atmos_grid_cartesC.F90:72
scale_calendar
module CALENDAR
Definition: scale_calendar.F90:13
scale_urban_grid_cartesc_index
module urban / grid / icosahedralA / index
Definition: scale_urban_grid_cartesC_index.F90:11
scale_atmos_grid_cartesc_index::imaxg
integer, public imaxg
Definition: scale_atmos_grid_cartesC_index.F90:71
scale_atmos_grid_cartesc_index::je_ing
integer, public je_ing
end point of the inner domain: cy, global
Definition: scale_atmos_grid_cartesC_index.F90:80
scale_file_history::file_history_truncate_3d
procedure(truncate_3d), pointer, public file_history_truncate_3d
Definition: scale_file_history.F90:101
scale_ocean_grid_cartesc
module ocean / grid / cartesianC
Definition: scale_ocean_grid_cartesC.F90:12
scale_io::h_institute
character(len=h_mid), public h_institute
for file header
Definition: scale_io.F90:51
scale_atmos_grid_cartesc_real
module Atmosphere GRID CartesC Real(real space)
Definition: scale_atmos_grid_cartesC_real.F90:11
scale_atmos_grid_cartesc::atmos_grid_cartesc_cdxg
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_cdxg
center coordinate [m]: x, global
Definition: scale_atmos_grid_cartesC.F90:80
scale_urban_grid_cartesc::urban_grid_cartesc_cz
real(rp), dimension(:), allocatable, public urban_grid_cartesc_cz
center coordinate [m]: z, local=global
Definition: scale_urban_grid_cartesC.F90:35
scale_file_history
module file_history
Definition: scale_file_history.F90:15
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_lon
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_lon
longitude [rad,0-2pi]
Definition: scale_atmos_grid_cartesC_real.F90:48
scale_file
module file
Definition: scale_file.F90:15
scale_atmos_grid_cartesc_index::jsgb
integer, public jsgb
start point of the inner domain: y, global
Definition: scale_atmos_grid_cartesC_index.F90:87
scale_atmos_grid_cartesc_index::imax
integer, public imax
Definition: scale_atmos_grid_cartesC_index.F90:37
scale_file_history::file_history_set_nowdate
subroutine, public file_history_set_nowdate(NOWDATE, NOWSUBSEC, NOWSTEP)
set now step
Definition: scale_file_history.F90:1697
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_volzuy
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_volzuy
control volume (zuy) [m3]
Definition: scale_atmos_grid_cartesC_real.F90:85
scale_io::h_source
character(len=h_mid), public h_source
for file header
Definition: scale_io.F90:50
scale_interp_vert::interp_vert_alloc_pres
subroutine, public interp_vert_alloc_pres(Kpres, KA, IA, JA)
Setup.
Definition: scale_interp_vert.F90:345
scale_file_h::file_hshort
integer, parameter, public file_hshort
Definition: scale_file_h.F90:19
scale_atmos_grid_cartesc_index::jmaxg
integer, public jmaxg
Definition: scale_atmos_grid_cartesC_index.F90:72
scale_prc_cartesc::prc_has_n
logical, public prc_has_n
Definition: scale_prc_cartesC.F90:48
scale_land_grid_cartesc_index
module land / grid / cartesianC / index
Definition: scale_land_grid_cartesC_index.F90:11
scale_atmos_grid_cartesc_index::jeb
integer, public jeb
Definition: scale_atmos_grid_cartesC_index.F90:66
scale_file_cartesc::axisattinfo
Definition: scale_file_cartesC.F90:96
scale_prc_cartesc::prc_has_e
logical, public prc_has_e
Definition: scale_prc_cartesC.F90:49
scale_urban_grid_cartesc_index::ukmax
integer, public ukmax
Definition: scale_urban_grid_cartesC_index.F90:32
scale_atmos_grid_cartesc::atmos_grid_cartesc_cbfx
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_cbfx
center buffer factor (0-1): x
Definition: scale_atmos_grid_cartesC.F90:70
scale_prc
module PROCESS
Definition: scale_prc.F90:11
scale_mapprojection
module Map projection
Definition: scale_mapprojection.F90:12
scale_file_history::file_history_set_dim
subroutine, public file_history_set_dim(name, ndims, nzcoords, dims, zcoords, start, count, mapping, area, area_x, area_y, volume, location, grid)
set dimension information
Definition: scale_file_history.F90:1470
scale_atmos_grid_cartesc_index::ie
integer, public ie
end point of inner domain: x, local
Definition: scale_atmos_grid_cartesC_index.F90:54
scale_atmos_grid_cartesc::atmos_grid_cartesc_fbfxg
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_fbfxg
face buffer factor (0-1): x, global
Definition: scale_atmos_grid_cartesC.F90:87
scale_io
module STDIO
Definition: scale_io.F90:10
scale_file_history::file_history_truncate_2d
procedure(truncate_2d), pointer, public file_history_truncate_2d
Definition: scale_file_history.F90:89
scale_atmos_grid_cartesc_index::jag
integer, public jag
Definition: scale_atmos_grid_cartesC_index.F90:74
scale_file_history::file_history_set_axis
subroutine, public file_history_set_axis(name, desc, units, dim, var, bounds, down, gsize, start)
set axis information
Definition: scale_file_history.F90:1579
scale_atmos_grid_cartesc_index::iag
integer, public iag
Definition: scale_atmos_grid_cartesC_index.F90:73
scale_prc_cartesc::prc_periodic_x
logical, public prc_periodic_x
periodic condition or not (X)?
Definition: scale_prc_cartesC.F90:52
scale_atmos_grid_cartesc_index::js_ing
integer, public js_ing
start point of the inner domain: cy, global
Definition: scale_atmos_grid_cartesC_index.F90:79
scale_atmos_grid_cartesc::atmos_grid_cartesc_cbfz
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_cbfz
center buffer factor (0-1): z
Definition: scale_atmos_grid_cartesC.F90:46
scale_ocean_grid_cartesc::ocean_grid_cartesc_fz
real(rp), dimension(:), allocatable, public ocean_grid_cartesc_fz
face coordinate [m]: z, local=global
Definition: scale_ocean_grid_cartesC.F90:36
scale_atmos_grid_cartesc::atmos_grid_cartesc_cxg
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_cxg
center coordinate [m]: x, global
Definition: scale_atmos_grid_cartesC.F90:75
scale_prc::prc_masterrank
integer, parameter, public prc_masterrank
master process in each communicator
Definition: scale_prc.F90:66
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_volzxv
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_volzxv
control volume (zxv) [m3]
Definition: scale_atmos_grid_cartesC_real.F90:86
scale_atmos_grid_cartesc_index
module atmosphere / grid / cartesC index
Definition: scale_atmos_grid_cartesC_index.F90:12
scale_const
module CONSTANT
Definition: scale_const.F90:11
scale_atmos_grid_cartesc_index::ia
integer, public ia
Definition: scale_atmos_grid_cartesC_index.F90:48
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_vol
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_vol
control volume (zxy) [m3]
Definition: scale_atmos_grid_cartesC_real.F90:83
scale_interp_vert::interp_available
logical, public interp_available
topography exists & vertical interpolation has meaning?
Definition: scale_interp_vert.F90:41
scale_atmos_grid_cartesc::atmos_grid_cartesc_fbfy
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_fbfy
face buffer factor (0-1): y
Definition: scale_atmos_grid_cartesC.F90:73
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_latuv
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_latuv
latitude at staggered point (uv) [rad,-pi,pi]
Definition: scale_atmos_grid_cartesC_real.F90:55
scale_landuse::landuse_frac_land
real(rp), dimension(:,:), allocatable, public landuse_frac_land
land fraction
Definition: scale_landuse.F90:54
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_areawxv_y
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_areawxv_y
virtical area (wxv, normal y) [m2]
Definition: scale_atmos_grid_cartesC_real.F90:69
scale_atmos_grid_cartesc_index::ie_ing
integer, public ie_ing
end point of the inner domain: cx, global
Definition: scale_atmos_grid_cartesC_index.F90:78
scale_prc_cartesc::prc_2drank
integer, dimension(:,:), allocatable, public prc_2drank
node index in 2D topology
Definition: scale_prc_cartesC.F90:44
scale_atmos_grid_cartesc_index::kmax
integer, public kmax
Definition: scale_atmos_grid_cartesC_index.F90:36
scale_interp_vert::interp_vert_xi2p
subroutine, public interp_vert_xi2p(Kpres, KA, KS, KE, IA, IS, IE, JA, JS, JE, var, var_P)
Definition: scale_interp_vert.F90:443
scale_urban_grid_cartesc
module urban / grid / cartesianC
Definition: scale_urban_grid_cartesC.F90:12
scale_atmos_grid_cartesc_index::is_ing
integer, public is_ing
start point of the inner domain: cx, global
Definition: scale_atmos_grid_cartesC_index.F90:77
scale_land_grid_cartesc_real
module land / grid / cartesianC / real
Definition: scale_land_grid_cartesC_real.F90:11
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_areazuy_x
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_areazuy_x
virtical area (zuy, normal x) [m2]
Definition: scale_atmos_grid_cartesC_real.F90:66
scale_atmos_grid_cartesc::atmos_grid_cartesc_fy
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_fy
face coordinate [m]: y, local
Definition: scale_atmos_grid_cartesC.F90:58
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_latuy
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_latuy
latitude at staggered point (uy) [rad,-pi,pi]
Definition: scale_atmos_grid_cartesC_real.F90:53
scale_prc_cartesc
module process / cartesC
Definition: scale_prc_cartesC.F90:11
scale_file_history::file_history_truncate_1d
procedure(truncate_1d), pointer, public file_history_truncate_1d
Definition: scale_file_history.F90:77
scale_file_history::file_history_setup
subroutine, public file_history_setup(title, source, institution, time_start, time_interval, time_units, time_since, calendar, default_basename, default_postfix_timelabel, default_zcoord, default_tinterval, default_tunit, default_taverage, default_datatype, myrank)
Setup.
Definition: scale_file_history.F90:342
scale_atmos_grid_cartesc::atmos_grid_cartesc_fbfz
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_fbfz
face buffer factor (0-1): z
Definition: scale_atmos_grid_cartesC.F90:47
scale_time::time_nowstep
integer, public time_nowstep
current step [number]
Definition: scale_time.F90:71
scale_ocean_grid_cartesc::ocean_grid_cartesc_cdz
real(rp), dimension(:), allocatable, public ocean_grid_cartesc_cdz
z-length of control volume [m]
Definition: scale_ocean_grid_cartesC.F90:37
scale_prof
module profiler
Definition: scale_prof.F90:11
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_area
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_area
horizontal area ( xy, normal z) [m2]
Definition: scale_atmos_grid_cartesC_real.F90:65
scale_atmos_grid_cartesc_index::is
integer, public is
start point of inner domain: x, local
Definition: scale_atmos_grid_cartesC_index.F90:53
scale_atmos_grid_cartesc_index::jhalo
integer, public jhalo
Definition: scale_atmos_grid_cartesC_index.F90:45
scale_ocean_grid_cartesc_index
module ocean / grid / cartesianC / index
Definition: scale_ocean_grid_cartesC_index.F90:11
scale_atmos_grid_cartesc::atmos_grid_cartesc_cbfxg
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_cbfxg
center buffer factor (0-1): x, global
Definition: scale_atmos_grid_cartesC.F90:85
scale_file_h
module file_h
Definition: scale_file_h.F90:11
scale_atmos_grid_cartesc_index::ja
integer, public ja
Definition: scale_atmos_grid_cartesC_index.F90:49
scale_time
module TIME
Definition: scale_time.F90:11
scale_atmos_grid_cartesc::atmos_grid_cartesc_fdx
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_fdx
x-length of grid(i+1) to grid(i) [m]
Definition: scale_atmos_grid_cartesC.F90:62
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_areauy
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_areauy
horizontal area ( uy, normal z) [m2]
Definition: scale_atmos_grid_cartesC_real.F90:70
scale_land_grid_cartesc::land_grid_cartesc_cz
real(rp), dimension(:), allocatable, public land_grid_cartesc_cz
center coordinate [m]: z, local=global
Definition: scale_land_grid_cartesC.F90:34
scale_atmos_grid_cartesc::atmos_grid_cartesc_fbfyg
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_fbfyg
face buffer factor (0-1): y, global
Definition: scale_atmos_grid_cartesC.F90:88
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_areazuv_x
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_areazuv_x
virtical area (zuv, normal x) [m2]
Definition: scale_atmos_grid_cartesC_real.F90:74
scale_interp_vert::interp_vert_xi2z
subroutine, public interp_vert_xi2z(KA, KS, KE, IA, IS, IE, JA, JS, JE, Xi, Z, var, var_Z)
Definition: scale_interp_vert.F90:186
scale_interp_vert::interp_vert_setcoef_pres
subroutine, public interp_vert_setcoef_pres(Kpres, KA, KS, KE, IA, IS, IE, JA, JS, JE, PRES, PRESh, SFC_PRES, Paxis)
Definition: scale_interp_vert.F90:374
scale_urban_grid_cartesc_index::uks
integer, public uks
Definition: scale_urban_grid_cartesC_index.F90:37
scale_atmos_grid_cartesc::atmos_grid_cartesc_fxg
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_fxg
face coordinate [m]: x, global
Definition: scale_atmos_grid_cartesC.F90:77
scale_interp_vert::interp_vert_xih2p
subroutine, public interp_vert_xih2p(Kpres, KA, KS, KE, IA, IS, IE, JA, JS, JE, var, var_P)
Definition: scale_interp_vert.F90:482
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_latxv
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_latxv
latitude at staggered point (xv) [rad,-pi,pi]
Definition: scale_atmos_grid_cartesC_real.F90:54
scale_atmos_grid_cartesc_index::ks
integer, public ks
start point of inner domain: z, local
Definition: scale_atmos_grid_cartesC_index.F90:51
scale_atmos_grid_cartesc_index::imaxb
integer, public imaxb
Definition: scale_atmos_grid_cartesC_index.F90:61
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_areazuv_y
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_areazuv_y
virtical area (zuv, normal y) [m2]
Definition: scale_atmos_grid_cartesC_real.F90:72
scale_land_grid_cartesc_real::land_grid_cartesc_real_vol
real(rp), dimension(:,:,:), allocatable, public land_grid_cartesc_real_vol
volume of grid cell
Definition: scale_land_grid_cartesC_real.F90:37
scale_atmos_grid_cartesc_index::iagb
integer, public iagb
Definition: scale_atmos_grid_cartesC_index.F90:75
scale_topography::topography_zsfc
real(rp), dimension(:,:), allocatable, public topography_zsfc
absolute ground height [m]
Definition: scale_topography.F90:38
scale_land_grid_cartesc::land_grid_cartesc_fz
real(rp), dimension(:), allocatable, public land_grid_cartesc_fz
face coordinate [m]: z, local=global
Definition: scale_land_grid_cartesC.F90:35
scale_file::file_get_cftunits
subroutine, public file_get_cftunits(date, tunits)
get unit of time
Definition: scale_file.F90:4832
scale_file_history_cartesc
module file / history_cartesC
Definition: scale_file_history_cartesC.F90:12
scale_ocean_grid_cartesc_real::ocean_grid_cartesc_real_vol
real(rp), dimension(:,:,:), allocatable, public ocean_grid_cartesc_real_vol
volume of grid cell
Definition: scale_ocean_grid_cartesC_real.F90:38
scale_atmos_grid_cartesc::atmos_grid_cartesc_fz
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_fz
face coordinate [m]: z, local
Definition: scale_atmos_grid_cartesC.F90:41
scale_atmos_grid_cartesc::atmos_grid_cartesc_name
character(len=7), parameter, public atmos_grid_cartesc_name
Definition: scale_atmos_grid_cartesC.F90:36
scale_time::time_dtsec
real(dp), public time_dtsec
time interval of model [sec]
Definition: scale_time.F90:33
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_lonuv
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_lonuv
longitude at staggered point (uv) [rad,0-2pi]
Definition: scale_atmos_grid_cartesC_real.F90:51
scale_atmos_grid_cartesc::atmos_grid_cartesc_cyg
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_cyg
center coordinate [m]: y, global
Definition: scale_atmos_grid_cartesC.F90:76
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_areazxy_x
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_areazxy_x
virtical area (zxy, normal x) [m2]
Definition: scale_atmos_grid_cartesC_real.F90:71
scale_atmos_grid_cartesc_index::jmax
integer, public jmax
Definition: scale_atmos_grid_cartesC_index.F90:38
scale_file_history::file_history_aggregate
logical, public file_history_aggregate
Definition: scale_file_history.F90:143
scale_time::time_nowdate
integer, dimension(6), public time_nowdate
current time [YYYY MM DD HH MM SS]
Definition: scale_time.F90:66
scale_atmos_grid_cartesc_index::js
integer, public js
start point of inner domain: y, local
Definition: scale_atmos_grid_cartesC_index.F90:55
scale_prc_cartesc::prc_num_y
integer, public prc_num_y
y length of 2D processor topology
Definition: scale_prc_cartesC.F90:42
scale_atmos_grid_cartesc::atmos_grid_cartesc_cy
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_cy
center coordinate [m]: y, local
Definition: scale_atmos_grid_cartesC.F90:56
scale_atmos_grid_cartesc_index::ieb
integer, public ieb
Definition: scale_atmos_grid_cartesC_index.F90:64
scale_file_history_cartesc::file_history_cartesc_setup
subroutine, public file_history_cartesc_setup
Setup.
Definition: scale_file_history_cartesC.F90:79
scale_prof::prof_rapend
subroutine, public prof_rapend(rapname_base, level, disable_barrier)
Save raptime.
Definition: scale_prof.F90:217
scale_ocean_grid_cartesc_real
module ocean / grid / cartesianC / real
Definition: scale_ocean_grid_cartesC_real.F90:12
scale_const::const_d2r
real(rp), public const_d2r
degree to radian
Definition: scale_const.F90:32
scale_atmos_grid_cartesc::atmos_grid_cartesc_cbfy
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_cbfy
center buffer factor (0-1): y
Definition: scale_atmos_grid_cartesC.F90:71
scale_land_grid_cartesc
module land / grid / cartesianC
Definition: scale_land_grid_cartesC.F90:11
scale_prc_cartesc::prc_num_x
integer, public prc_num_x
x length of 2D processor topology
Definition: scale_prc_cartesC.F90:41
scale_land_grid_cartesc_index::lke
integer, public lke
Definition: scale_land_grid_cartesC_index.F90:41
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_areaxv
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_areaxv
horizontal area ( xv, normal z) [m2]
Definition: scale_atmos_grid_cartesC_real.F90:73
scale_atmos_grid_cartesc::atmos_grid_cartesc_cdx
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_cdx
x-length of control volume [m]
Definition: scale_atmos_grid_cartesC.F90:60
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_lat
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_lat
latitude [rad,-pi,pi]
Definition: scale_atmos_grid_cartesC_real.F90:52
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_fz
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_fz
geopotential height [m] (wxy)
Definition: scale_atmos_grid_cartesC_real.F90:42
scale_atmos_grid_cartesc
module atmosphere / grid / cartesC
Definition: scale_atmos_grid_cartesC.F90:12
scale_landuse
module LANDUSE
Definition: scale_landuse.F90:19
scale_atmos_grid_cartesc::atmos_grid_cartesc_fdz
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_fdz
z-length of grid(i+1) to grid(i) [m]
Definition: scale_atmos_grid_cartesC.F90:43
scale_file_history_cartesc::file_history_cartesc_truncate_1d
subroutine file_history_cartesc_truncate_1d(src, dim_type, zcoord, fill_halo, dst)
truncate 1D data to history buffer
Definition: scale_file_history_cartesC.F90:493
scale_atmos_grid_cartesc::atmos_grid_cartesc_cz
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_cz
center coordinate [m]: z, local
Definition: scale_atmos_grid_cartesC.F90:40
scale_atmos_grid_cartesc::atmos_grid_cartesc_cbfyg
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_cbfyg
center buffer factor (0-1): y, global
Definition: scale_atmos_grid_cartesC.F90:86
scale_urban_grid_cartesc_real
module urban / grid / cartesianC / real
Definition: scale_urban_grid_cartesC_real.F90:12
scale_interp_vert::interp_vert_xih2zh
subroutine, public interp_vert_xih2zh(KA, KS, KE, IA, IS, IE, JA, JS, JE, Xih, Zh, var, var_Z)
Definition: scale_interp_vert.F90:268
scale_const::const_undef
real(rp), public const_undef
Definition: scale_const.F90:41
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_areazxv_y
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_areazxv_y
virtical area (zxv, normal y) [m2]
Definition: scale_atmos_grid_cartesC_real.F90:67
scale_atmos_grid_cartesc_index::jsb
integer, public jsb
Definition: scale_atmos_grid_cartesC_index.F90:65
scale_prc_cartesc::prc_has_w
logical, public prc_has_w
Definition: scale_prc_cartesC.F90:47
scale_interp_vert
module INTERPOLATION vertical
Definition: scale_interp_vert.F90:11
scale_io::io_fid_conf
integer, public io_fid_conf
Config file ID.
Definition: scale_io.F90:56
scale_land_grid_cartesc_index::lks
integer, public lks
Definition: scale_land_grid_cartesC_index.F90:40
scale_time::time_startdaysec
real(dp), public time_startdaysec
second of start time [sec]
Definition: scale_time.F90:75
scale_atmos_grid_cartesc_index::je
integer, public je
end point of inner domain: y, local
Definition: scale_atmos_grid_cartesC_index.F90:56
scale_prc_cartesc::prc_twod
logical, public prc_twod
2D experiment
Definition: scale_prc_cartesC.F90:55
scale_file_cartesc
module file / cartesianC
Definition: scale_file_cartesC.F90:11
scale_mapprojection::mapprojection_mappinginfo
type(mappinginfo), public mapprojection_mappinginfo
Definition: scale_mapprojection.F90:104
scale_atmos_grid_cartesc::atmos_grid_cartesc_cx
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_cx
center coordinate [m]: x, local
Definition: scale_atmos_grid_cartesC.F90:55
scale_atmos_grid_cartesc_index::jsga
integer, public jsga
start point of the full domain: cy, global
Definition: scale_atmos_grid_cartesC_index.F90:83
scale_urban_grid_cartesc_index::uke
integer, public uke
Definition: scale_urban_grid_cartesC_index.F90:38
scale_urban_grid_cartesc::urban_grid_cartesc_fz
real(rp), dimension(:), allocatable, public urban_grid_cartesc_fz
face coordinate [m]: z, local=global
Definition: scale_urban_grid_cartesC.F90:36
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_lonuy
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_lonuy
longitude at staggered point (uy) [rad,0-2pi]
Definition: scale_atmos_grid_cartesC_real.F90:49
scale_atmos_grid_cartesc::atmos_grid_cartesc_fdxg
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_fdxg
center coordinate [m]: x, global
Definition: scale_atmos_grid_cartesC.F90:82
scale_ocean_grid_cartesc_index::oks
integer, public oks
Definition: scale_ocean_grid_cartesC_index.F90:37
scale_file_history_cartesc::file_history_cartesc_truncate_3d
subroutine file_history_cartesc_truncate_3d(src, dim_type, zcoord, fill_halo, dst)
truncate 3D data to history buffer
Definition: scale_file_history_cartesC.F90:636