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