SCALE-RM
scale_file_history.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
10 !-------------------------------------------------------------------------------
11 ! Warning: This file was generated from file/scale_file_history.F90.erb.
12 ! Do not edit this file.
13 !-------------------------------------------------------------------------------
14 #include "scalelib.h"
16  !-----------------------------------------------------------------------------
17  !
18  !++ Used modules
19  !
20  use scale_precision
21  use scale_io
22  use scale_prof
23  use scale_prc, only: &
24  prc_abort
25  !-----------------------------------------------------------------------------
26  implicit none
27  private
28  !-----------------------------------------------------------------------------
29  !
30  !++ Public procedures
31  !
32  public :: file_history_setup
33  public :: file_history_query
34  public :: file_history_reg
35  public :: file_history_put
36  public :: file_history_write
37  public :: file_history_in
38  public :: file_history_set_dim
39  public :: file_history_set_axis
40  public :: file_history_set_associatedcoordinate
41  public :: file_history_set_attribute
42  public :: file_history_set_nowdate
43  public :: file_history_set_disable
44  public :: file_history_finalize
45 
46  interface file_history_query
47  module procedure file_history_query_name
48  module procedure file_history_query_id
49  end interface file_history_query
50 
51  interface file_history_put
52  module procedure file_history_put_0d
53  module procedure file_history_put_1d
54  module procedure file_history_put_2d
55  module procedure file_history_put_3d
56  module procedure file_history_put_4d
57  end interface file_history_put
58 
59  interface file_history_in
60  module procedure file_history_in_0d
61  module procedure file_history_in_1d
62  module procedure file_history_in_2d
63  module procedure file_history_in_3d
64  module procedure file_history_in_4d
65  end interface file_history_in
66 
67  abstract interface
68  subroutine truncate_1d( src, dim_type, zcoord, fill_halo, dsc )
69  import rp, dp
70  real(RP), intent(in) :: src(:)
71  character(len=*), intent(in) :: dim_type
72  character(len=*), intent(in) :: zcoord
73  logical, intent(in) :: fill_halo
74  real(DP), intent(out) :: dsc(:)
75  end subroutine truncate_1d
76  end interface
77  procedure(truncate_1d), pointer :: file_history_truncate_1d => null()
78  public :: file_history_truncate_1d
79  abstract interface
80  subroutine truncate_2d( src, dim_type, zcoord, fill_halo, dsc )
81  import rp, dp
82  real(rp), intent(in) :: src(:,:)
83  character(len=*), intent(in) :: dim_type
84  character(len=*), intent(in) :: zcoord
85  logical, intent(in) :: fill_halo
86  real(dp), intent(out) :: dsc(:)
87  end subroutine truncate_2d
88  end interface
89  procedure(truncate_2d), pointer :: file_history_truncate_2d => null()
90  public :: file_history_truncate_2d
91  abstract interface
92  subroutine truncate_3d( src, dim_type, zcoord, fill_halo, dsc )
93  import rp, dp
94  real(rp), intent(in) :: src(:,:,:)
95  character(len=*), intent(in) :: dim_type
96  character(len=*), intent(in) :: zcoord
97  logical, intent(in) :: fill_halo
98  real(dp), intent(out) :: dsc(:)
99  end subroutine truncate_3d
100  end interface
101  procedure(truncate_3d), pointer :: file_history_truncate_3d => null()
102  public :: file_history_truncate_3d
103  abstract interface
104  subroutine truncate_4d( src, dim_type, zcoord, fill_halo, dsc )
105  import rp, dp
106  real(rp), intent(in) :: src(:,:,:,:)
107  character(len=*), intent(in) :: dim_type
108  character(len=*), intent(in) :: zcoord
109  logical, intent(in) :: fill_halo
110  real(dp), intent(out) :: dsc(:)
111  end subroutine truncate_4d
112  end interface
113  procedure(truncate_4d), pointer :: file_history_truncate_4d => null()
114  public :: file_history_truncate_4d
115 
116 
117  interface file_history_set_associatedcoordinate
121  end interface file_history_set_associatedcoordinate
122 
123  interface file_history_set_attribute
124  module procedure file_history_set_attribute_text
125  module procedure file_history_set_attribute_logical
126  module procedure file_history_set_attribute_int
127  module procedure file_history_set_attribute_int_ary
128  module procedure file_history_set_attribute_float
129  module procedure file_history_set_attribute_float_ary
130  module procedure file_history_set_attribute_double
131  module procedure file_history_set_attribute_double_ary
132  end interface file_history_set_attribute
133 
134 
135  !-----------------------------------------------------------------------------
136  !
137  !++ included parameters
138  !
139  !-----------------------------------------------------------------------------
140  !
141  !++ Public parameters & variables
142  !
143  logical, public :: file_history_aggregate
144  !-----------------------------------------------------------------------------
145  !
146  !++ Private procedures
147  !
148  private :: file_history_create
149  private :: file_history_close
150  private :: file_history_add_variable
151  private :: file_history_write_axes
152  private :: file_history_write_onevar
153  private :: file_history_output_list
154  private :: file_history_check
155 
156  !-----------------------------------------------------------------------------
157  !
158  !++ Private parameters & variables
159  !
160 
161  type request
162  character(len=H_SHORT) :: name
163  character(len=H_SHORT) :: outname
164  character(len=H_LONG) :: basename
165  logical :: postfix_timelabel
166  character(len=H_SHORT) :: zcoord
167  integer :: dstep
168  logical :: taverage
169  integer :: dtype
170  character(len=H_SHORT) :: cell_measures
171  logical :: registered
172  end type request
173 
174  type var_out
175  character(len=H_SHORT) :: name
176  character(len=H_SHORT) :: outname
177  character(len=H_LONG) :: basename
178  logical :: postfix_timelabel
179  character(len=H_SHORT) :: zcoord
180  integer :: zid
181  integer :: dstep
182  logical :: taverage
183  integer :: dtype
184 
185  integer :: fid
186  integer :: vid
187  character(len=H_LONG) :: desc
188  character(len=H_SHORT) :: units
189  character(len=H_SHORT) :: standard_name
190  integer :: dimid
191  character(len=H_SHORT) :: cell_measures
192  integer :: waitstep
193  integer :: laststep_write
194  integer :: laststep_put
195  logical :: flag_clear
196  integer :: size
197  real(dp) :: timesum
198  real(dp), pointer :: varsum(:)
199  logical :: fill_halo
200  end type var_out
201 
202  integer, parameter :: file_history_variant_max = 10
203  type var_in
204  character(len=H_SHORT) :: name
205  integer :: nvariants
206  integer :: variants(file_history_variant_max)
207  end type var_in
208 
209  type dim
210  character(len=H_SHORT) :: name
211  integer :: ndims
212  integer :: nzcoords
213  character(len=H_SHORT), pointer :: dims(:,:)
214  integer , pointer :: start(:,:)
215  integer , pointer :: count(:,:)
216  integer , pointer :: size(:)
217  character(len=H_SHORT), pointer :: zcoords(:)
218  character(len=H_SHORT) :: mapping
219  character(len=H_SHORT) :: area
220  character(len=H_SHORT) :: area_x
221  character(len=H_SHORT) :: area_y
222  character(len=H_SHORT) :: volume
223  character(len=H_SHORT) :: location
224  character(len=H_SHORT) :: grid
225  end type dim
226 
227  type axis
228  character(len=H_SHORT) :: name
229  character(len=H_LONG) :: desc
230  character(len=H_SHORT) :: units
231  character(len=H_SHORT) :: dim
232  integer :: dim_size
233  real(dp), pointer :: var(:)
234  real(dp), pointer :: bounds(:,:)
235  logical :: down
236  integer :: gdim_size ! global dimension size
237  integer :: start ! global array start index
238  end type axis
239 
240  type assoc
241  character(len=H_SHORT) :: name
242  character(len=H_LONG) :: desc
243  character(len=H_SHORT) :: units
244  integer :: ndims
245  character(len=H_SHORT) :: dims(4)
246  integer :: dtype
247  real(dp), pointer :: var(:)
248  integer :: start(4) ! global array start indices
249  integer :: count(4) ! global array request lengths
250  end type assoc
251 
252  integer, parameter :: i_text = 1, i_int = 2, i_float = 3, i_double = 4
253  type attr
254  character(len=H_SHORT) :: varname
255  character(len=H_MID) :: key
256  integer :: type
257  character(len=H_LONG) :: text
258  integer, pointer :: int(:)
259  real(sp), pointer :: float(:)
260  real(dp), pointer :: double(:)
261  logical :: add_variable
262  end type attr
263 
264  ! From upstream side of the library
265  integer :: file_history_myrank
266 
267 
268  real(dp) :: file_history_startdaysec
269  real(dp) :: file_history_dtsec
270  character(len=H_MID) :: file_history_time_since
271 
272  ! From NAMELIST or upstream side of the library
273  character(len=H_MID) :: file_history_title
274  character(len=H_MID) :: file_history_source
275  character(len=H_MID) :: file_history_institution
276 
277  character(len=H_MID) :: file_history_time_units
278  character(len=H_SHORT) :: file_history_calendar
279  logical :: file_history_output_step0 = .false.
280  integer :: file_history_output_wait_step
281  integer :: file_history_output_switch_step
282  integer :: file_history_output_switch_laststep
283  logical :: file_history_error_putmiss = .true.
284 
285  ! working
286  integer, parameter :: file_history_req_max = 1000
287  integer :: file_history_nreqs = 0
288  type(request), allocatable :: file_history_req(:)
289 
290  integer :: file_history_nitems = 0
291  type(var_out), allocatable :: file_history_vars(:)
292 
293  integer :: file_history_nvar_inputs
294  type(var_in), allocatable :: file_history_var_inputs(:)
295 
296  integer, parameter :: file_history_dim_max = 30
297  integer :: file_history_ndims = 0
298  type(dim) :: file_history_dims(file_history_dim_max)
299 
300  integer, parameter :: file_history_axis_max = 100
301  integer :: file_history_naxes = 0
302  type(axis) :: file_history_axes(file_history_axis_max)
303 
304  integer, parameter :: file_history_assoc_max = 40
305  integer :: file_history_nassocs = 0
306  type(assoc) :: file_history_assocs(file_history_assoc_max)
307 
308  integer, parameter :: file_history_attr_max = 200
309  integer :: file_history_nattrs = 0
310  type(attr) :: file_history_attrs(file_history_attr_max)
311 
312  integer :: file_history_nowdate(6)
313  real(dp) :: file_history_nowsubsec
314  integer :: file_history_nowstep
315 
316  character(len=H_MID) :: file_history_options = ''
317 
318 
319  logical :: file_history_disabled = .true.
320 
321  integer :: laststep_write = -1
322  logical :: list_outputed = .false.
323  logical :: debug = .false.
324 
325  !-----------------------------------------------------------------------------
326 contains
327  !-----------------------------------------------------------------------------
329  !-----------------------------------------------------------------------------
330  subroutine file_history_setup( &
331  title, source, institution, &
332  time_start, time_interval, &
333  time_units, time_since, calendar, &
334  default_basename, &
335  default_postfix_timelabel, &
336  default_zcoord, &
337  default_tinterval, &
338  default_tunit, &
339  default_taverage, &
340  default_datatype, &
341  myrank )
342  use scale_file_h, only: &
343  file_real4, &
344  file_real8
345  use scale_file, only: &
347  use scale_calendar, only: &
349  implicit none
350 
351  character(len=*), intent(in) :: title
352  character(len=*), intent(in) :: source
353  character(len=*), intent(in) :: institution
354  real(dp), intent(in) :: time_start
355  real(dp), intent(in) :: time_interval
356 
357  character(len=*), intent(in), optional :: time_units
358  character(len=*), intent(in), optional :: time_since
359  character(len=*), intent(in), optional :: calendar
360  character(len=*), intent(in), optional :: default_basename
361  logical, intent(in), optional :: default_postfix_timelabel
362  character(len=*), intent(in), optional :: default_zcoord
363  real(dp), intent(in), optional :: default_tinterval
364  character(len=*), intent(in), optional :: default_tunit
365  logical, intent(in), optional :: default_taverage
366  character(len=*), intent(in), optional :: default_datatype
367  integer, intent(in), optional :: myrank
368 
369  character(len=H_LONG) :: file_history_default_basename
370  logical :: file_history_default_postfix_timelabel
371  character(len=H_SHORT) :: file_history_default_zcoord
372  real(dp) :: file_history_default_tinterval
373  character(len=H_SHORT) :: file_history_default_tunit
374  logical :: file_history_default_taverage
375  character(len=H_SHORT) :: file_history_default_datatype
378  real(dp) :: file_history_output_wait
379  character(len=H_SHORT) :: file_history_output_wait_tunit
380  real(dp) :: file_history_output_switch_tinterval
381  character(len=H_SHORT) :: file_history_output_switch_tunit
382 
383  namelist / param_file_history / &
384  file_history_title, &
385  file_history_source, &
386  file_history_institution, &
387  file_history_time_units, &
388  file_history_default_basename, &
389  file_history_default_postfix_timelabel, &
390  file_history_default_zcoord, &
391  file_history_default_tinterval, &
392  file_history_default_tunit, &
393  file_history_default_taverage, &
394  file_history_default_datatype, &
395  file_history_output_step0, &
396  file_history_output_wait, &
397  file_history_output_wait_tunit, &
398  file_history_output_switch_tinterval, &
399  file_history_output_switch_tunit, &
400  file_history_error_putmiss, &
402  file_history_options, &
403  debug
404 
405  character(len=H_SHORT) :: name
406  character(len=H_SHORT) :: outname
407  character(len=H_LONG) :: basename
408  logical :: postfix_timelabel
409  character(len=H_SHORT) :: zcoord
410  real(dp) :: tinterval
411  character(len=H_SHORT) :: tunit
412  logical :: taverage
413  character(len=H_SHORT) :: datatype
414 
415  namelist / history_item / &
416  name, &
417  outname, &
418  basename, &
419  postfix_timelabel, &
420  zcoord, &
421  tinterval, &
422  tunit, &
423  taverage, &
424  datatype
425 
426 
427  integer :: reqid
428  real(dp) :: dtsec
429  integer :: dstep
430 
431  integer :: ierr
432  integer :: n, id
433 
434  intrinsic size
435  !---------------------------------------------------------------------------
436 
437  log_newline
438  log_info("FILE_HISTORY_Setup",*) 'Setup'
439 
440  ! setup
441  file_history_myrank = myrank
442 
443  file_history_startdaysec = time_start
444  file_history_dtsec = time_interval
445  if ( present(time_since) ) then
446  file_history_time_since = time_since
447  else
448  file_history_time_since = ''
449  endif
450 
451  if ( present(calendar) ) then
452  file_history_calendar = calendar
453  else
454  file_history_calendar = ""
455  end if
456 
457  file_history_time_units = 'seconds'
458  file_history_default_basename = ''
459  file_history_default_postfix_timelabel = .false.
460  file_history_default_zcoord = ''
461  file_history_default_tinterval = -1.0_dp
462  file_history_default_tunit = 'SEC'
463  file_history_default_taverage = .false.
464  file_history_default_datatype = 'REAL4'
465  file_history_output_wait = 0.0_dp
466  file_history_output_wait_tunit = 'SEC'
467  file_history_output_switch_tinterval = -1.0_dp
468  file_history_output_switch_tunit = 'SEC'
469 
471 
472  !--- read namelist
473  file_history_title = title
474  file_history_source = source
475  file_history_institution = institution
476  if( present(time_units) ) file_history_time_units = time_units
477  if( present(default_basename) ) file_history_default_basename = default_basename
478  if( present(default_postfix_timelabel) ) file_history_default_postfix_timelabel = default_postfix_timelabel
479  if( present(default_zcoord) ) file_history_default_zcoord = default_zcoord
480  if( present(default_tinterval) ) file_history_default_tinterval = default_tinterval
481  if( present(default_tunit) ) file_history_default_tunit = default_tunit
482  if( present(default_taverage) ) file_history_default_taverage = default_taverage
483  if( present(default_datatype) ) file_history_default_datatype = default_datatype
484 
485  !--- read namelist
486  rewind(io_fid_conf)
487  read(io_fid_conf,nml=param_file_history,iostat=ierr)
488  if( ierr < 0 ) then !--- missing
489  log_info("FILE_HISTORY_Setup",*) 'Not found namelist. Default used.'
490  elseif( ierr > 0 ) then !--- fatal error
491  log_error("FILE_HISTORY_Setup",*) 'Not appropriate names in namelist PARAM_FILE_HISTORY. Check!'
492  call prc_abort
493  endif
494  log_nml(param_file_history)
495 
496 
497 
498  if ( file_history_output_wait >= 0.0_dp ) then
499  call calendar_unit2sec( dtsec, file_history_output_wait, file_history_output_wait_tunit )
500  file_history_output_wait_step = int( dtsec / file_history_dtsec )
501  else
502  log_error("FILE_HISTORY_Setup",*) 'FILE_HISTORY_OUTPUT_WAIT must be positive. STOP'
503  call prc_abort
504  endif
505 
506  if ( file_history_output_switch_tinterval >= 0.0_dp ) then
507  call calendar_unit2sec( dtsec, file_history_output_switch_tinterval, file_history_output_switch_tunit )
508  file_history_output_switch_step = int( dtsec / file_history_dtsec )
509  else
510  file_history_output_switch_step = -1
511  endif
512  file_history_output_switch_laststep = 0
513 
514 
515  ! count history request
516  file_history_nreqs = 0
517  if ( io_fid_conf > 0 ) rewind(io_fid_conf)
518  do n = 1, file_history_req_max
519  name = ''
520  outname = 'undefined'
521  basename = file_history_default_basename
522 
523  read(io_fid_conf,nml=history_item,iostat=ierr)
524  if( ierr /= 0 ) exit
525  if( basename == '' .OR. name == '' .OR. outname == '' ) cycle ! invalid HISTORY_ITEM
526 
527  file_history_nreqs = file_history_nreqs + 1
528  enddo
529 
530  if ( file_history_nreqs > file_history_req_max ) then
531  log_error("FILE_HISTORY_Setup",*) 'request of history file is exceed! n >', file_history_req_max
532  call prc_abort
533  elseif( file_history_nreqs == 0 ) then
534  log_info("FILE_HISTORY_Setup",*) 'No history file specified.'
535  return
536  endif
537 
538  allocate( file_history_req(file_history_nreqs) )
539 
540  ! read history request
541  reqid = 0
542  if ( io_fid_conf > 0 ) rewind(io_fid_conf)
543  do n = 1, file_history_req_max
544  ! set default
545  name = ''
546  outname = 'undefined'
547  basename = file_history_default_basename
548  postfix_timelabel = file_history_default_postfix_timelabel
549  zcoord = file_history_default_zcoord
550  tinterval = file_history_default_tinterval
551  tunit = file_history_default_tunit
552  taverage = file_history_default_taverage
553  datatype = file_history_default_datatype
554 
555  read(io_fid_conf,nml=history_item,iostat=ierr)
556  if ( ierr < 0 ) then
557  exit ! no more items
558  elseif( ierr > 0 ) then
559  log_error("FILE_HISTORY_Setup",*) 'Not appropriate names in namelist HISTORY_ITEM. Check!'
560  call prc_abort
561  endif
562  if( basename == '' .OR. name == '' .OR. outname == '' ) cycle ! invalid HISTORY_ITEM
563 
564  log_nml(history_item)
565 
566  ! check duplicated request
567  if ( outname == 'undefined' ) outname = name ! set default name
568  do id = 1, reqid
569  if ( file_history_req(id)%outname == outname ) then
570  log_error("FILE_HISTORY_Setup",*) 'Same name of history output is already registered. Check!', trim(outname)
571  call prc_abort
572  endif
573  enddo
574 
575  reqid = reqid + 1
576 
577  file_history_req(reqid)%name = name
578  file_history_req(reqid)%outname = outname
579  file_history_req(reqid)%basename = basename
580  file_history_req(reqid)%postfix_timelabel = postfix_timelabel
581  if( file_history_output_switch_step >= 0 ) file_history_req(reqid)%postfix_timelabel = .true. ! force true
582  file_history_req(reqid)%zcoord = zcoord
583  file_history_req(reqid)%taverage = taverage
584 
585  call calendar_unit2sec( dtsec, tinterval, tunit )
586  dstep = int( dtsec / file_history_dtsec )
587 
588  if ( dtsec <= 0.d0 ) then
589  log_error("FILE_HISTORY_Setup",*) 'Not appropriate time interval. Check!', trim(name), tinterval, trim(tunit)
590  call prc_abort
591  endif
592 
593  if ( abs(dtsec-real(dstep,kind=dp)*file_history_dtsec) > dtsec*1.e-3_dp ) then
594  log_error("FILE_HISTORY_Setup",*) 'time interval must be a multiple of delta t. (interval,dt)=', dtsec, file_history_dtsec
595  call prc_abort
596  endif
597 
598  file_history_req(reqid)%dstep = dstep
599 
600  if ( datatype == 'REAL4' ) then
601  file_history_req(reqid)%dtype = file_real4
602  elseif( datatype == 'REAL8' ) then
603  file_history_req(reqid)%dtype = file_real8
604  else
605  log_error("FILE_HISTORY_Setup",*) 'Not appropriate DATATYPE. Check!', datatype
606  call prc_abort
607  endif
608 
609  file_history_req(reqid)%registered = .false.
610  enddo
611 
612  log_newline
613  log_info("FILE_HISTORY_Setup",*) 'Number of requested history item : ', file_history_nreqs
614  log_info("FILE_HISTORY_Setup",*) 'Output default data type : ', trim(file_history_default_datatype)
615  log_info("FILE_HISTORY_Setup",*) 'Output value at the initial step? : ', file_history_output_step0
616  if ( file_history_output_wait_step > 0 ) then
617  log_info("FILE_HISTORY_Setup",*) 'Time when the output is suppressed [step] : ', file_history_output_wait_step
618  end if
619  if ( file_history_output_switch_step >= 0 ) then
620  log_info("FILE_HISTORY_Setup",*) 'Interval for switching the file [step] : ', file_history_output_switch_step
621  end if
622  log_info("FILE_HISTORY_Setup",*) 'Check if requested item is not registered? : ', file_history_error_putmiss
623 
624  file_history_nitems = 0
625  allocate( file_history_vars(file_history_nreqs) )
626 
627  file_history_nvar_inputs = 0
628  allocate( file_history_var_inputs(file_history_nreqs) )
629 
630  file_history_truncate_1d => file_history_truncate_1d_default
631  file_history_truncate_2d => file_history_truncate_2d_default
632  file_history_truncate_3d => file_history_truncate_3d_default
633  file_history_truncate_4d => file_history_truncate_4d_default
634 
635  file_history_disabled = .false.
636 
637  return
638  end subroutine file_history_setup
639 
640  !-----------------------------------------------------------------------------
642  !-----------------------------------------------------------------------------
643  subroutine file_history_reg( &
644  name, desc, unit, &
645  itemid, &
646  standard_name, &
647  ndims, dim_type, &
648  cell_measures, &
649  fill_halo )
650  implicit none
651 
652  character(len=*), intent(in) :: name
653  character(len=*), intent(in) :: desc
654  character(len=*), intent(in) :: unit
655 
656  integer, intent(out) :: itemid
657 
658  character(len=*), intent(in), optional :: standard_name
659  integer, intent(in), optional :: ndims
660  character(len=*), intent(in), optional :: dim_type
661  character(len=*), intent(in), optional :: cell_measures
662  logical, intent(in), optional :: fill_halo
663 
664  character(len=H_SHORT) :: standard_name_
665  character(len=H_SHORT) :: cell_measures_
666  integer :: dimid, iid
667  integer :: n
668  !---------------------------------------------------------------------------
669 
670  itemid = -1
671  if ( file_history_nreqs == 0 ) return
672 
673  itemid = file_history_find_id( name )
674  if ( itemid > 0 ) return ! already registered
675 
676  call prof_rapstart('FILE_HISTORY_OUT', 2)
677 
678  if ( len_trim(name) >= h_short ) then
679  log_error("FILE_HISTORY_reg",'(1x,A,I2,A,A)') 'Length of history name should be <= ', h_short-1 ,' chars. name=', trim(name)
680  call prc_abort
681  endif
682 
683  ! standard_name
684  if ( present(standard_name) ) then
685  standard_name_ = standard_name
686  else
687  standard_name_ = ""
688  end if
689 
690  ! get dimension id
691  if ( file_history_ndims < 1 ) then
692  log_error("FILE_HISTORY_reg",*) 'at least one dim_type must be registerd with FILE_HISTORY_set_dim. name=', trim(name)
693  call prc_abort
694  end if
695  dimid = -1
696  if ( present(dim_type) ) then
697  do n = 1, file_history_ndims
698  if ( file_history_dims(n)%name == dim_type ) then
699  dimid = n
700  exit
701  end if
702  end do
703  if ( dimid == -1 ) then
704  log_error("FILE_HISTORY_reg",*) 'dim_type must be registerd with FILE_HISTORY_set_dim: ', trim(dim_type) ,' name=', trim(name)
705  call prc_abort
706  end if
707  else if ( present(ndims) ) then
708  do n = 1, file_history_ndims
709  if ( file_history_dims(n)%ndims == ndims ) then
710  dimid = n
711  exit
712  end if
713  end do
714  if ( dimid == -1 ) then
715  log_error("FILE_HISTORY_reg",'(a,i1,a)') 'dim_type of ', ndims, 'D must be registerd with FILE_HISTORY_set_dim. name=', trim(name)
716  call prc_abort
717  end if
718  else
719  ! ndims = 3 is assumed as default
720  do n = 1, file_history_ndims
721  if ( file_history_dims(n)%ndims == 3 ) then
722  dimid = n
723  exit
724  end if
725  end do
726  if ( dimid == -1 ) then
727  log_error("FILE_HISTORY_reg",'(a,i1,a)') 'dim_type or ndims must be specified. name=', trim(name)
728  call prc_abort
729  end if
730  end if
731 
732  if ( present(cell_measures) ) then
733  select case ( cell_measures )
734  case ( "area" )
735  if ( file_history_dims(dimid)%area == "" ) then
736  log_error("FILE_HISTORY_reg",*) 'area is not supported for cell_measures. name=', trim(name)
737  call prc_abort
738  end if
739  case ( "area_z" )
740  if ( file_history_dims(dimid)%area == "" ) then
741  log_error("FILE_HISTORY_reg",*) 'area_z is not supported for cell_measures. name=', trim(name)
742  call prc_abort
743  end if
744  case ( "area_x" )
745  if ( file_history_dims(dimid)%area_x == "" ) then
746  log_error("FILE_HISTORY_reg",*) 'area_x is not supported for cell_measures. name=', trim(name)
747  call prc_abort
748  end if
749  case ( "area_y" )
750  if ( file_history_dims(dimid)%area_y == "" ) then
751  log_error("FILE_HISTORY_reg",*) 'area_y is not supported for cell_measures. name=', trim(name)
752  call prc_abort
753  end if
754  case ( "volume" )
755  if ( file_history_dims(dimid)%volume == "" ) then
756  log_error("FILE_HISTORY_reg",*) 'volume is not supported for cell_measures. name=', trim(name)
757  call prc_abort
758  end if
759  case default
760  log_error("FILE_HISTORY_reg",*) 'cell_measures must be "area" or "volume". name=', trim(name)
761  call prc_abort
762  end select
763  cell_measures_ = cell_measures
764  else if ( file_history_dims(dimid)%ndims == 2 ) then
765  cell_measures_ = "area"
766  else if ( file_history_dims(dimid)%ndims == 3 ) then
767  cell_measures_ = "volume"
768  else
769  cell_measures_ = ""
770  end if
771 
772  if ( file_history_dims(dimid)%nzcoords > 1 ) then
773 
774  itemid = -1
775  do n = 1, file_history_dims(dimid)%nzcoords
776  if ( file_history_dims(dimid)%zcoords(n) == "model" ) then
777  call file_history_add_variable( name, desc, unit, standard_name_, & ! (in)
778  dimid, & ! (in)
779  file_history_dims(dimid)%zcoords(n), & ! (in)
780  iid, & ! (out)
781  cell_measures = cell_measures_, & ! (in)
782  fill_halo = fill_halo ) ! (in)
783  else
784  call file_history_add_variable( name, desc, unit, standard_name_, & ! (in)
785  dimid, & ! (in)
786  file_history_dims(dimid)%zcoords(n), & ! (in)
787  iid, & ! (out)
788  fill_halo = fill_halo ) ! (in)
789  end if
790  if ( iid > 0 ) itemid = iid
791  end do
792 
793  else
794 
795  call file_history_add_variable( name, desc, unit, standard_name_, & ! (in)
796  dimid, & ! (in)
797  "model", & ! (in)
798  itemid, & ! (out)
799  cell_measures = cell_measures_, & ! (in)
800  fill_halo = fill_halo ) ! (in)
801 
802  end if
803 
804  call prof_rapend('FILE_HISTORY_OUT', 2)
805 
806  return
807  end subroutine file_history_reg
808  !-----------------------------------------------------------------------------
809  ! interface FILE_HISTORY_Put
810  !-----------------------------------------------------------------------------
811  subroutine file_history_put_0d( &
812  itemid, &
813  var )
814  use scale_file_h, only: &
815  rmiss => file_rmiss
816  use scale_const, only: &
817  undef => const_undef, &
818  eps => const_eps
819  implicit none
820 
821  integer, intent(in) :: itemid
822  real(rp), intent(in) :: var
823 
824  integer :: dimid
825  real(dp), allocatable :: buffer(:)
826  real(dp) :: dt
827  integer :: idx
828  logical :: do_put
829 
830  integer :: i, id
831 
832  intrinsic shape
833  !---------------------------------------------------------------------------
834 
835  if ( file_history_disabled ) return
836  if ( itemid < 0 ) return
837 
838  call file_history_query( itemid, do_put )
839  if ( .not. do_put ) return
840 
841  call prof_rapstart('FILE_HISTORY_OUT', 2)
842 
843  do i = 1, file_history_var_inputs(itemid)%nvariants
844  id = file_history_var_inputs(itemid)%variants(i)
845 
846  dt = ( file_history_nowstep - file_history_vars(id)%laststep_put ) * file_history_dtsec
847 
848  if ( dt < eps .AND. ( .NOT. file_history_vars(id)%taverage ) ) then
849  log_error("FILE_HISTORY_Put_0D",*) 'variable was put two times before output!: ', &
850  trim(file_history_vars(id)%name), file_history_nowstep, file_history_vars(id)%laststep_put
851  call prc_abort
852  endif
853 
854  if ( file_history_vars(id)%flag_clear ) then ! time to purge
855  file_history_vars(id)%timesum = 0.0_dp
856  if ( file_history_vars(id)%taverage ) file_history_vars(id)%varsum(:) = 0.0_dp
857  endif
858 
859  dimid = file_history_vars(id)%dimid
860  if ( file_history_vars(id)%taverage ) then
861  if ( file_history_vars(id)%varsum(1) /= rmiss ) then
862  if ( var /= undef ) then
863  file_history_vars(id)%varsum(1) = file_history_vars(id)%varsum(1) + var * dt
864  else
865  file_history_vars(id)%varsum(1) = rmiss
866  end if
867  end if
868  file_history_vars(id)%timesum = file_history_vars(id)%timesum + dt
869  else
870  file_history_vars(id)%varsum(1) = var
871  file_history_vars(id)%timesum = 0.0_dp
872  endif
873 
874  file_history_vars(id)%laststep_put = file_history_nowstep
875  file_history_vars(id)%flag_clear = .false.
876 
877  end do ! variants
878 
879  call prof_rapend('FILE_HISTORY_OUT', 2)
880 
881  return
882  end subroutine file_history_put_0d
883 
884  !-----------------------------------------------------------------------------
886  !-----------------------------------------------------------------------------
887  subroutine file_history_in_0d( &
888  var, &
889  name, desc, unit, &
890  standard_name, &
891  dim_type )
892  implicit none
893 
894  real(RP), intent(in) :: var
895  character(len=*), intent(in) :: name
896  character(len=*), intent(in) :: desc
897  character(len=*), intent(in) :: unit
898 
899  character(len=*), intent(in), optional :: standard_name
900  character(len=*), intent(in), optional :: dim_type
901 
902  logical, parameter :: fill_halo = .false.
903 
904  integer, parameter :: ndim = 0
905  integer :: itemid
906  logical :: do_put
907  !---------------------------------------------------------------------------
908 
909  if ( file_history_disabled ) return
910 
911  ! Check whether the item has been already registered
912  call file_history_reg( name, desc, unit, & ! [IN]
913  itemid, & ! [OUT]
914  standard_name=standard_name, & ! [IN]
915  ndims=ndim, & ! [IN]
916  dim_type=dim_type, & ! [IN]
917  fill_halo=fill_halo ) ! [IN]
918 
919  if ( itemid < 0 ) return
920 
921  ! Check whether it is time to input the item
922  call file_history_query( itemid, do_put ) ! [IN], [OUT]
923 
924  if ( do_put ) call file_history_put( itemid, var )
925 
926  return
927  end subroutine file_history_in_0d
928 
929  !-----------------------------------------------------------------------------
930  subroutine file_history_put_1d( &
931  itemid, &
932  var )
933  use scale_file_h, only: &
934  rmiss => file_rmiss
935  use scale_const, only: &
936  undef => const_undef, &
937  eps => const_eps
938  implicit none
939 
940  integer, intent(in) :: itemid
941  real(RP), intent(in) :: var(:)
942 
943  integer :: dimid
944  real(DP), allocatable :: buffer(:)
945  real(DP) :: dt
946  integer :: idx
947  logical :: do_put
948 
949  integer :: i, id
950 
951  intrinsic shape
952  !---------------------------------------------------------------------------
953 
954  if ( file_history_disabled ) return
955  if ( itemid < 0 ) return
956 
957  call file_history_query( itemid, do_put )
958  if ( .not. do_put ) return
959 
960  call prof_rapstart('FILE_HISTORY_OUT', 2)
961 
962  do i = 1, file_history_var_inputs(itemid)%nvariants
963  id = file_history_var_inputs(itemid)%variants(i)
964 
965  dt = ( file_history_nowstep - file_history_vars(id)%laststep_put ) * file_history_dtsec
966 
967  if ( dt < eps .AND. ( .NOT. file_history_vars(id)%taverage ) ) then
968  log_error("FILE_HISTORY_Put_1D",*) 'variable was put two times before output!: ', &
969  trim(file_history_vars(id)%name), file_history_nowstep, file_history_vars(id)%laststep_put
970  call prc_abort
971  endif
972 
973  if ( file_history_vars(id)%flag_clear ) then ! time to purge
974  file_history_vars(id)%timesum = 0.0_dp
975  if ( file_history_vars(id)%taverage ) file_history_vars(id)%varsum(:) = 0.0_dp
976  endif
977 
978  dimid = file_history_vars(id)%dimid
979  if ( file_history_vars(id)%taverage ) then
980  allocate( buffer( file_history_vars(id)%size ) )
981  call file_history_truncate_1d( var(:), & ! (in)
982  file_history_dims(dimid)%name, & ! (in)
983  file_history_vars(id)%zcoord, & ! (in)
984  file_history_vars(id)%fill_halo, & ! (in)
985  buffer(:) ) ! (out)
986  do idx = 1, file_history_vars(id)%size
987  if ( file_history_vars(id)%varsum(idx) /= rmiss ) then
988  if ( buffer(idx) /= undef ) then
989  file_history_vars(id)%varsum(idx) = file_history_vars(id)%varsum(idx) + buffer(idx) * dt
990  else
991  file_history_vars(id)%varsum(idx) = rmiss
992  end if
993  end if
994  enddo
995  deallocate( buffer )
996  file_history_vars(id)%timesum = file_history_vars(id)%timesum + dt
997  else
998  call file_history_truncate_1d( var(:), & ! (in)
999  file_history_dims(dimid)%name, & ! (in)
1000  file_history_vars(id)%zcoord, & ! (in)
1001  file_history_vars(id)%fill_halo, & ! (in)
1002  file_history_vars(id)%varsum(:) ) ! (out)
1003  file_history_vars(id)%timesum = 0.0_dp
1004  endif
1005 
1006  file_history_vars(id)%laststep_put = file_history_nowstep
1007  file_history_vars(id)%flag_clear = .false.
1008 
1009  end do ! variants
1010 
1011  call prof_rapend('FILE_HISTORY_OUT', 2)
1012 
1013  return
1014  end subroutine file_history_put_1d
1015 
1016  !-----------------------------------------------------------------------------
1018  !-----------------------------------------------------------------------------
1019  subroutine file_history_in_1d( &
1020  var, &
1021  name, desc, unit, &
1022  standard_name, &
1023  dim_type )
1024  implicit none
1025 
1026  real(RP), intent(in) :: var(:)
1027  character(len=*), intent(in) :: name
1028  character(len=*), intent(in) :: desc
1029  character(len=*), intent(in) :: unit
1030 
1031  character(len=*), intent(in), optional :: standard_name
1032  character(len=*), intent(in), optional :: dim_type
1033 
1034  logical, parameter :: fill_halo = .false.
1035 
1036  integer, parameter :: ndim = 1
1037  integer :: itemid
1038  logical :: do_put
1039  !---------------------------------------------------------------------------
1040 
1041  if ( file_history_disabled ) return
1042 
1043  ! Check whether the item has been already registered
1044  call file_history_reg( name, desc, unit, & ! [IN]
1045  itemid, & ! [OUT]
1046  standard_name=standard_name, & ! [IN]
1047  ndims=ndim, & ! [IN]
1048  dim_type=dim_type, & ! [IN]
1049  fill_halo=fill_halo ) ! [IN]
1050 
1051  if ( itemid < 0 ) return
1052 
1053  ! Check whether it is time to input the item
1054  call file_history_query( itemid, do_put ) ! [IN], [OUT]
1055 
1056  if ( do_put ) call file_history_put( itemid, var(:) )
1057 
1058  return
1059  end subroutine file_history_in_1d
1060 
1061  !-----------------------------------------------------------------------------
1062  subroutine file_history_put_2d( &
1063  itemid, &
1064  var )
1065  use scale_file_h, only: &
1066  rmiss => file_rmiss
1067  use scale_const, only: &
1068  undef => const_undef, &
1069  eps => const_eps
1070  implicit none
1071 
1072  integer, intent(in) :: itemid
1073  real(RP), intent(in) :: var(:,:)
1074 
1075  integer :: dimid
1076  real(DP), allocatable :: buffer(:)
1077  real(DP) :: dt
1078  integer :: idx
1079  logical :: do_put
1080 
1081  integer :: i, id
1082 
1083  intrinsic shape
1084  !---------------------------------------------------------------------------
1085 
1086  if ( file_history_disabled ) return
1087  if ( itemid < 0 ) return
1088 
1089  call file_history_query( itemid, do_put )
1090  if ( .not. do_put ) return
1091 
1092  call prof_rapstart('FILE_HISTORY_OUT', 2)
1093 
1094  do i = 1, file_history_var_inputs(itemid)%nvariants
1095  id = file_history_var_inputs(itemid)%variants(i)
1096 
1097  dt = ( file_history_nowstep - file_history_vars(id)%laststep_put ) * file_history_dtsec
1098 
1099  if ( dt < eps .AND. ( .NOT. file_history_vars(id)%taverage ) ) then
1100  log_error("FILE_HISTORY_Put_2D",*) 'variable was put two times before output!: ', &
1101  trim(file_history_vars(id)%name), file_history_nowstep, file_history_vars(id)%laststep_put
1102  call prc_abort
1103  endif
1104 
1105  if ( file_history_vars(id)%flag_clear ) then ! time to purge
1106  file_history_vars(id)%timesum = 0.0_dp
1107  if ( file_history_vars(id)%taverage ) file_history_vars(id)%varsum(:) = 0.0_dp
1108  endif
1109 
1110  dimid = file_history_vars(id)%dimid
1111  if ( file_history_vars(id)%taverage ) then
1112  allocate( buffer( file_history_vars(id)%size ) )
1113  call file_history_truncate_2d( var(:,:), & ! (in)
1114  file_history_dims(dimid)%name, & ! (in)
1115  file_history_vars(id)%zcoord, & ! (in)
1116  file_history_vars(id)%fill_halo, & ! (in)
1117  buffer(:) ) ! (out)
1118  do idx = 1, file_history_vars(id)%size
1119  if ( file_history_vars(id)%varsum(idx) /= rmiss ) then
1120  if ( buffer(idx) /= undef ) then
1121  file_history_vars(id)%varsum(idx) = file_history_vars(id)%varsum(idx) + buffer(idx) * dt
1122  else
1123  file_history_vars(id)%varsum(idx) = rmiss
1124  end if
1125  end if
1126  enddo
1127  deallocate( buffer )
1128  file_history_vars(id)%timesum = file_history_vars(id)%timesum + dt
1129  else
1130  call file_history_truncate_2d( var(:,:), & ! (in)
1131  file_history_dims(dimid)%name, & ! (in)
1132  file_history_vars(id)%zcoord, & ! (in)
1133  file_history_vars(id)%fill_halo, & ! (in)
1134  file_history_vars(id)%varsum(:) ) ! (out)
1135  file_history_vars(id)%timesum = 0.0_dp
1136  endif
1137 
1138  file_history_vars(id)%laststep_put = file_history_nowstep
1139  file_history_vars(id)%flag_clear = .false.
1140 
1141  end do ! variants
1142 
1143  call prof_rapend('FILE_HISTORY_OUT', 2)
1144 
1145  return
1146  end subroutine file_history_put_2d
1147 
1148  !-----------------------------------------------------------------------------
1150  !-----------------------------------------------------------------------------
1151  subroutine file_history_in_2d( &
1152  var, &
1153  name, desc, unit, &
1154  standard_name, &
1155  dim_type, &
1156  fill_halo )
1157  implicit none
1158 
1159  real(RP), intent(in) :: var(:,:)
1160  character(len=*), intent(in) :: name
1161  character(len=*), intent(in) :: desc
1162  character(len=*), intent(in) :: unit
1163 
1164  character(len=*), intent(in), optional :: standard_name
1165  character(len=*), intent(in), optional :: dim_type
1166  logical, intent(in), optional :: fill_halo
1167 
1168  integer, parameter :: ndim = 2
1169  integer :: itemid
1170  logical :: do_put
1171  !---------------------------------------------------------------------------
1172 
1173  if ( file_history_disabled ) return
1174 
1175  ! Check whether the item has been already registered
1176  call file_history_reg( name, desc, unit, & ! [IN]
1177  itemid, & ! [OUT]
1178  standard_name=standard_name, & ! [IN]
1179  ndims=ndim, & ! [IN]
1180  dim_type=dim_type, & ! [IN]
1181  fill_halo=fill_halo ) ! [IN]
1182 
1183  if ( itemid < 0 ) return
1184 
1185  ! Check whether it is time to input the item
1186  call file_history_query( itemid, do_put ) ! [IN], [OUT]
1187 
1188  if ( do_put ) call file_history_put( itemid, var(:,:) )
1189 
1190  return
1191  end subroutine file_history_in_2d
1192 
1193  !-----------------------------------------------------------------------------
1194  subroutine file_history_put_3d( &
1195  itemid, &
1196  var )
1197  use scale_file_h, only: &
1198  rmiss => file_rmiss
1199  use scale_const, only: &
1200  undef => const_undef, &
1201  eps => const_eps
1202  implicit none
1203 
1204  integer, intent(in) :: itemid
1205  real(RP), intent(in) :: var(:,:,:)
1206 
1207  integer :: dimid
1208  real(DP), allocatable :: buffer(:)
1209  real(DP) :: dt
1210  integer :: idx
1211  logical :: do_put
1212 
1213  integer :: i, id
1214 
1215  intrinsic shape
1216  !---------------------------------------------------------------------------
1217 
1218  if ( file_history_disabled ) return
1219  if ( itemid < 0 ) return
1220 
1221  call file_history_query( itemid, do_put )
1222  if ( .not. do_put ) return
1223 
1224  call prof_rapstart('FILE_HISTORY_OUT', 2)
1225 
1226  do i = 1, file_history_var_inputs(itemid)%nvariants
1227  id = file_history_var_inputs(itemid)%variants(i)
1228 
1229  dt = ( file_history_nowstep - file_history_vars(id)%laststep_put ) * file_history_dtsec
1230 
1231  if ( dt < eps .AND. ( .NOT. file_history_vars(id)%taverage ) ) then
1232  log_error("FILE_HISTORY_Put_3D",*) 'variable was put two times before output!: ', &
1233  trim(file_history_vars(id)%name), file_history_nowstep, file_history_vars(id)%laststep_put
1234  call prc_abort
1235  endif
1236 
1237  if ( file_history_vars(id)%flag_clear ) then ! time to purge
1238  file_history_vars(id)%timesum = 0.0_dp
1239  if ( file_history_vars(id)%taverage ) file_history_vars(id)%varsum(:) = 0.0_dp
1240  endif
1241 
1242  dimid = file_history_vars(id)%dimid
1243  if ( file_history_vars(id)%taverage ) then
1244  allocate( buffer( file_history_vars(id)%size ) )
1245  call file_history_truncate_3d( var(:,:,:), & ! (in)
1246  file_history_dims(dimid)%name, & ! (in)
1247  file_history_vars(id)%zcoord, & ! (in)
1248  file_history_vars(id)%fill_halo, & ! (in)
1249  buffer(:) ) ! (out)
1250  do idx = 1, file_history_vars(id)%size
1251  if ( file_history_vars(id)%varsum(idx) /= rmiss ) then
1252  if ( buffer(idx) /= undef ) then
1253  file_history_vars(id)%varsum(idx) = file_history_vars(id)%varsum(idx) + buffer(idx) * dt
1254  else
1255  file_history_vars(id)%varsum(idx) = rmiss
1256  end if
1257  end if
1258  enddo
1259  deallocate( buffer )
1260  file_history_vars(id)%timesum = file_history_vars(id)%timesum + dt
1261  else
1262  call file_history_truncate_3d( var(:,:,:), & ! (in)
1263  file_history_dims(dimid)%name, & ! (in)
1264  file_history_vars(id)%zcoord, & ! (in)
1265  file_history_vars(id)%fill_halo, & ! (in)
1266  file_history_vars(id)%varsum(:) ) ! (out)
1267  file_history_vars(id)%timesum = 0.0_dp
1268  endif
1269 
1270  file_history_vars(id)%laststep_put = file_history_nowstep
1271  file_history_vars(id)%flag_clear = .false.
1272 
1273  end do ! variants
1274 
1275  call prof_rapend('FILE_HISTORY_OUT', 2)
1276 
1277  return
1278  end subroutine file_history_put_3d
1279 
1280  !-----------------------------------------------------------------------------
1282  !-----------------------------------------------------------------------------
1283  subroutine file_history_in_3d( &
1284  var, &
1285  name, desc, unit, &
1286  standard_name, &
1287  dim_type, &
1288  fill_halo )
1289  implicit none
1290 
1291  real(RP), intent(in) :: var(:,:,:)
1292  character(len=*), intent(in) :: name
1293  character(len=*), intent(in) :: desc
1294  character(len=*), intent(in) :: unit
1295 
1296  character(len=*), intent(in), optional :: standard_name
1297  character(len=*), intent(in), optional :: dim_type
1298  logical, intent(in), optional :: fill_halo
1299 
1300  integer, parameter :: ndim = 3
1301  integer :: itemid
1302  logical :: do_put
1303  !---------------------------------------------------------------------------
1304 
1305  if ( file_history_disabled ) return
1306 
1307  ! Check whether the item has been already registered
1308  call file_history_reg( name, desc, unit, & ! [IN]
1309  itemid, & ! [OUT]
1310  standard_name=standard_name, & ! [IN]
1311  ndims=ndim, & ! [IN]
1312  dim_type=dim_type, & ! [IN]
1313  fill_halo=fill_halo ) ! [IN]
1314 
1315  if ( itemid < 0 ) return
1316 
1317  ! Check whether it is time to input the item
1318  call file_history_query( itemid, do_put ) ! [IN], [OUT]
1319 
1320  if ( do_put ) call file_history_put( itemid, var(:,:,:) )
1321 
1322  return
1323  end subroutine file_history_in_3d
1324 
1325  !-----------------------------------------------------------------------------
1326  subroutine file_history_put_4d( &
1327  itemid, &
1328  var )
1329  use scale_file_h, only: &
1330  rmiss => file_rmiss
1331  use scale_const, only: &
1332  undef => const_undef, &
1333  eps => const_eps
1334  implicit none
1335 
1336  integer, intent(in) :: itemid
1337  real(RP), intent(in) :: var(:,:,:,:)
1338 
1339  integer :: dimid
1340  real(DP), allocatable :: buffer(:)
1341  real(DP) :: dt
1342  integer :: idx
1343  logical :: do_put
1344 
1345  integer :: i, id
1346 
1347  intrinsic shape
1348  !---------------------------------------------------------------------------
1349 
1350  if ( file_history_disabled ) return
1351  if ( itemid < 0 ) return
1352 
1353  call file_history_query( itemid, do_put )
1354  if ( .not. do_put ) return
1355 
1356  call prof_rapstart('FILE_HISTORY_OUT', 2)
1357 
1358  do i = 1, file_history_var_inputs(itemid)%nvariants
1359  id = file_history_var_inputs(itemid)%variants(i)
1360 
1361  dt = ( file_history_nowstep - file_history_vars(id)%laststep_put ) * file_history_dtsec
1362 
1363  if ( dt < eps .AND. ( .NOT. file_history_vars(id)%taverage ) ) then
1364  log_error("FILE_HISTORY_Put_4D",*) 'variable was put two times before output!: ', &
1365  trim(file_history_vars(id)%name), file_history_nowstep, file_history_vars(id)%laststep_put
1366  call prc_abort
1367  endif
1368 
1369  if ( file_history_vars(id)%flag_clear ) then ! time to purge
1370  file_history_vars(id)%timesum = 0.0_dp
1371  if ( file_history_vars(id)%taverage ) file_history_vars(id)%varsum(:) = 0.0_dp
1372  endif
1373 
1374  dimid = file_history_vars(id)%dimid
1375  if ( file_history_vars(id)%taverage ) then
1376  allocate( buffer( file_history_vars(id)%size ) )
1377  call file_history_truncate_4d( var(:,:,:,:), & ! (in)
1378  file_history_dims(dimid)%name, & ! (in)
1379  file_history_vars(id)%zcoord, & ! (in)
1380  file_history_vars(id)%fill_halo, & ! (in)
1381  buffer(:) ) ! (out)
1382  do idx = 1, file_history_vars(id)%size
1383  if ( file_history_vars(id)%varsum(idx) /= rmiss ) then
1384  if ( buffer(idx) /= undef ) then
1385  file_history_vars(id)%varsum(idx) = file_history_vars(id)%varsum(idx) + buffer(idx) * dt
1386  else
1387  file_history_vars(id)%varsum(idx) = rmiss
1388  end if
1389  end if
1390  enddo
1391  deallocate( buffer )
1392  file_history_vars(id)%timesum = file_history_vars(id)%timesum + dt
1393  else
1394  call file_history_truncate_4d( var(:,:,:,:), & ! (in)
1395  file_history_dims(dimid)%name, & ! (in)
1396  file_history_vars(id)%zcoord, & ! (in)
1397  file_history_vars(id)%fill_halo, & ! (in)
1398  file_history_vars(id)%varsum(:) ) ! (out)
1399  file_history_vars(id)%timesum = 0.0_dp
1400  endif
1401 
1402  file_history_vars(id)%laststep_put = file_history_nowstep
1403  file_history_vars(id)%flag_clear = .false.
1404 
1405  end do ! variants
1406 
1407  call prof_rapend('FILE_HISTORY_OUT', 2)
1408 
1409  return
1410  end subroutine file_history_put_4d
1411 
1412  !-----------------------------------------------------------------------------
1414  !-----------------------------------------------------------------------------
1415  subroutine file_history_in_4d( &
1416  var, &
1417  name, desc, unit, &
1418  standard_name, &
1419  dim_type, &
1420  fill_halo )
1421  implicit none
1422 
1423  real(RP), intent(in) :: var(:,:,:,:)
1424  character(len=*), intent(in) :: name
1425  character(len=*), intent(in) :: desc
1426  character(len=*), intent(in) :: unit
1427 
1428  character(len=*), intent(in), optional :: standard_name
1429  character(len=*), intent(in), optional :: dim_type
1430  logical, intent(in), optional :: fill_halo
1431 
1432  integer, parameter :: ndim = 4
1433  integer :: itemid
1434  logical :: do_put
1435  !---------------------------------------------------------------------------
1436 
1437  if ( file_history_disabled ) return
1438 
1439  ! Check whether the item has been already registered
1440  call file_history_reg( name, desc, unit, & ! [IN]
1441  itemid, & ! [OUT]
1442  standard_name=standard_name, & ! [IN]
1443  ndims=ndim, & ! [IN]
1444  dim_type=dim_type, & ! [IN]
1445  fill_halo=fill_halo ) ! [IN]
1446 
1447  if ( itemid < 0 ) return
1448 
1449  ! Check whether it is time to input the item
1450  call file_history_query( itemid, do_put ) ! [IN], [OUT]
1451 
1452  if ( do_put ) call file_history_put( itemid, var(:,:,:,:) )
1453 
1454  return
1455  end subroutine file_history_in_4d
1456 
1457 
1458  !-----------------------------------------------------------------------------
1460  !-----------------------------------------------------------------------------
1461  subroutine file_history_set_dim( &
1462  name, &
1463  ndims, nzcoords, &
1464  dims, zcoords, &
1465  start, count, &
1466  mapping, &
1467  area, area_x, area_y, &
1468  volume, &
1469  location, grid )
1470  implicit none
1471 
1472  character(len=*), intent(in) :: name
1473  integer, intent(in) :: ndims
1474  integer, intent(in) :: nzcoords
1475  character(len=*), intent(in) :: dims(ndims,nzcoords)
1476  character(len=*), intent(in) :: zcoords(nzcoords)
1477  integer, intent(in) :: start(ndims,nzcoords)
1478  integer, intent(in) :: count(ndims,nzcoords)
1479 
1480  character(len=*), intent(in), optional :: mapping
1481  character(len=*), intent(in), optional :: area
1482  character(len=*), intent(in), optional :: area_x
1483  character(len=*), intent(in), optional :: area_y
1484  character(len=*), intent(in), optional :: volume
1485  character(len=*), intent(in), optional :: location
1486  character(len=*), intent(in), optional :: grid
1487 
1488  integer :: id
1489  integer :: size, n, m
1490 
1491  if ( file_history_ndims >= file_history_dim_max ) then
1492  log_error("FILE_HISTORY_Set_Dim",*) 'number of dimension exceed max limit: ', file_history_dim_max
1493  call prc_abort
1494  end if
1495  file_history_ndims = file_history_ndims + 1
1496  id = file_history_ndims
1497 
1498  allocate( file_history_dims(id)%dims (max(ndims,1),nzcoords) )
1499  allocate( file_history_dims(id)%start(max(ndims,1),nzcoords) )
1500  allocate( file_history_dims(id)%count(max(ndims,1),nzcoords) )
1501  allocate( file_history_dims(id)%zcoords(nzcoords) )
1502  allocate( file_history_dims(id)%size(nzcoords) )
1503 
1504  file_history_dims(id)%name = name
1505  file_history_dims(id)%ndims = ndims
1506  file_history_dims(id)%nzcoords = nzcoords
1507  file_history_dims(id)%zcoords(:) = zcoords(:)
1508  if ( ndims > 0 ) then
1509  file_history_dims(id)%dims(:,:) = dims(:,:)
1510  file_history_dims(id)%start(:,:) = start(:,:)
1511  file_history_dims(id)%count(:,:) = count(:,:)
1512  else ! 0D
1513  file_history_dims(id)%dims(1,1) = ""
1514  file_history_dims(id)%start(1,1) = 1
1515  file_history_dims(id)%count(1,1) = 1
1516  end if
1517 
1518  do m = 1, nzcoords
1519  size = 1
1520  do n = 1, ndims
1521  size = size * count(n,m)
1522  end do
1523  file_history_dims(id)%size(m) = size
1524  end do
1525 
1526  if ( present(mapping) ) then
1527  file_history_dims(id)%mapping = mapping
1528  else
1529  file_history_dims(id)%mapping = ""
1530  end if
1531 
1532  if ( present(area) ) then
1533  file_history_dims(id)%area = area
1534  else
1535  file_history_dims(id)%area = ""
1536  end if
1537  if ( present(area_x) ) then
1538  file_history_dims(id)%area_x = area_x
1539  else
1540  file_history_dims(id)%area_x = ""
1541  end if
1542  if ( present(area_y) ) then
1543  file_history_dims(id)%area_y = area_y
1544  else
1545  file_history_dims(id)%area_y = ""
1546  end if
1547  if ( present(volume) ) then
1548  file_history_dims(id)%volume = volume
1549  else
1550  file_history_dims(id)%volume = ""
1551  end if
1552 
1553  if ( present(location) ) then
1554  file_history_dims(id)%location = location
1555  if ( present(grid) ) then
1556  file_history_dims(id)%grid = "grid_"//trim(grid)
1557  else
1558  file_history_dims(id)%grid = "grid"
1559  end if
1560  else
1561  file_history_dims(id)%location = ""
1562  file_history_dims(id)%grid = ""
1563  end if
1564 
1565  return
1566  end subroutine file_history_set_dim
1567 
1568  !-----------------------------------------------------------------------------
1570  !-----------------------------------------------------------------------------
1571  subroutine file_history_set_axis( &
1572  name, desc, units, &
1573  dim, &
1574  var, &
1575  bounds, &
1576  down, &
1577  gsize, &
1578  start )
1579  implicit none
1580 
1581  character(len=*), intent(in) :: name
1582  character(len=*), intent(in) :: desc
1583  character(len=*), intent(in) :: units
1584  character(len=*), intent(in) :: dim
1585  real(rp), intent(in) :: var(:)
1586 
1587  real(rp), intent(in), optional :: bounds(:,:)
1588  logical, intent(in), optional :: down
1589  integer, intent(in), optional :: gsize ! global dim size
1590  integer, intent(in), optional :: start ! global subarray start indices
1591 
1592  integer :: dim_size
1593  integer :: id
1594 
1595  intrinsic size
1596  !---------------------------------------------------------------------------
1597 
1598  dim_size = size(var)
1599 
1600  if ( file_history_naxes >= file_history_axis_max ) then
1601  log_error("FILE_HISTORY_Set_Axis",*) 'Number of axis exceeds the limit.'
1602  call prc_abort
1603  endif
1604 
1605  file_history_naxes = file_history_naxes + 1
1606  id = file_history_naxes
1607 
1608  allocate( file_history_axes(id)%var(dim_size) )
1609 
1610  file_history_axes(id)%name = name
1611  file_history_axes(id)%desc = desc
1612  file_history_axes(id)%units = units
1613  file_history_axes(id)%dim = dim
1614  file_history_axes(id)%dim_size = dim_size
1615  file_history_axes(id)%var(:) = var(:)
1616 
1617  if ( present(down) ) then
1618  file_history_axes(id)%down = down
1619  else
1620  file_history_axes(id)%down = .false.
1621  endif
1622  if ( present(gsize) ) then ! global dimension size
1623  file_history_axes(id)%gdim_size = gsize
1624  else
1625  file_history_axes(id)%gdim_size = -1
1626  end if
1627  if ( present(start) ) then ! global subarray starting indices
1628  file_history_axes(id)%start = start
1629  else
1630  file_history_axes(id)%start = 1
1631  end if
1632 
1633  if ( present(bounds) ) then
1634  allocate( file_history_axes(id)%bounds(2,dim_size) )
1635  file_history_axes(id)%bounds(:,:) = bounds(:,:)
1636  end if
1637 
1638  return
1639  end subroutine file_history_set_axis
1640 
1641  !-----------------------------------------------------------------------------
1642  subroutine file_history_write
1643  use scale_file, only: &
1644  file_enddef, &
1645  file_flush
1646  implicit none
1647 
1648  integer :: fid, prev_fid
1649  integer :: id
1650  !---------------------------------------------------------------------------
1651 
1652  if ( file_history_disabled ) return
1653 
1654  call prof_rapstart('FILE_HISTORY_OUT', 2)
1655 
1656  ! Write registered history variables to history file
1657  do id = 1, file_history_nitems
1658  call file_history_write_onevar( id, file_history_nowstep ) ! [IN]
1659  enddo
1660 
1661  ! when using PnetCDF, the above FILE_HISTORY_Write() only posts write requests
1662  ! Now we need to commit the requests to the file
1663  prev_fid = -1
1664  do id = 1, file_history_nitems
1665  fid = file_history_vars(id)%fid
1666  if ( fid > 0 .AND. fid /= prev_fid ) then
1667  call file_flush( fid )
1668  prev_fid = fid
1669  endif
1670  enddo
1671 
1672  ! check time to switching output file
1673  if ( file_history_output_switch_step >= 0 &
1674  .AND. file_history_nowstep-file_history_output_switch_laststep > file_history_output_switch_step ) then
1675 
1676  call file_history_close
1677 
1678  log_info("FILE_HISTORY_Write",*) 'FILE_HISTORY file is switched.'
1679 
1680  do id = 1, file_history_nitems
1681  file_history_vars(id)%fid = -1 ! reset
1682  file_history_vars(id)%vid = -1 ! reset
1683  enddo
1684 
1685  file_history_output_switch_laststep = file_history_nowstep - 1
1686  endif
1687 
1688  call prof_rapend('FILE_HISTORY_OUT', 2)
1689 
1690  return
1691  end subroutine file_history_write
1692 
1693  !-----------------------------------------------------------------------------
1695  !-----------------------------------------------------------------------------
1696  subroutine file_history_set_nowdate( NOWDATE, NOWSUBSEC, NOWSTEP )
1697  integer, intent(in) :: nowdate(:)
1698  real(dp), intent(in) :: nowsubsec
1699  integer, intent(in) :: nowstep
1700 
1701  file_history_nowdate(:) = nowdate(:)
1702  file_history_nowsubsec = nowsubsec
1703  file_history_nowstep = nowstep
1704 
1705  return
1706  end subroutine file_history_set_nowdate
1707 
1708  !-----------------------------------------------------------------------------
1710  !-----------------------------------------------------------------------------
1711  subroutine file_history_set_disable( switch )
1712  implicit none
1713 
1714  logical, intent(in) :: switch
1715  !---------------------------------------------------------------------------
1716 
1717  file_history_disabled = switch
1718 
1719  return
1720  end subroutine file_history_set_disable
1721 
1722  !-----------------------------------------------------------------------------
1724  !-----------------------------------------------------------------------------
1725  subroutine file_history_finalize
1727  call file_history_close
1728 
1729  return
1730  end subroutine file_history_finalize
1731 
1732 
1733  !! private procedures
1734 
1735  !-----------------------------------------------------------------------------
1736  subroutine file_history_check( &
1737  name, zcoord, &
1738  itemid )
1739  implicit none
1740 
1741  character(len=*), intent(in) :: name
1742  character(len=*), intent(in) :: zcoord
1743  integer, intent(out) :: itemid
1744 
1745  integer :: id, i
1746  !---------------------------------------------------------------------------
1747 
1748  !--- search existing item
1749  do itemid = 1, file_history_nvar_inputs
1750  if ( name == file_history_var_inputs(itemid)%name ) then ! match name
1751  do i = 1, file_history_var_inputs(itemid)%nvariants
1752  id = file_history_var_inputs(itemid)%variants(i)
1753  !--- check zcoord
1754  if ( file_history_vars(id)%zcoord == zcoord ) return
1755  end do
1756  end if
1757  end do
1758  itemid = -1
1759 
1760  return
1761  end subroutine file_history_check
1762 
1763  !-----------------------------------------------------------------------------
1764  subroutine file_history_add_variable( &
1765  name, desc, units, &
1766  standard_name, &
1767  dimid, &
1768  zcoord, &
1769  itemid, &
1770  cell_measures, &
1771  fill_halo )
1772  use scale_file_h, only: &
1774  implicit none
1775  character(len=*), intent(in) :: name
1776  character(len=*), intent(in) :: desc
1777  character(len=*), intent(in) :: units
1778  character(len=*), intent(in) :: standard_name
1779  integer, intent(in) :: dimid
1780  character(len=*), intent(in) :: zcoord
1781  integer, intent(out) :: itemid
1782  character(len=*), intent(in), optional :: cell_measures
1783  logical, intent(in), optional :: fill_halo
1784 
1785  integer :: reqid, zid, id
1786  logical :: existed
1787  integer :: n, m
1788 
1789  intrinsic size
1790  !---------------------------------------------------------------------------
1791 
1792  call file_history_check( name, zcoord, & ! (in)
1793  itemid )
1794 
1795  if ( itemid > 0 ) return
1796 
1797  do reqid = 1, file_history_nreqs
1798 
1799  if ( file_history_req(reqid)%registered ) cycle
1800  if ( name /= file_history_req(reqid)%name ) cycle
1801 
1802  if ( file_history_dims(dimid)%nzcoords == 1 .or. &
1803  zcoord == file_history_req(reqid)%zcoord ) then
1804 
1805  file_history_req(reqid)%registered = .true.
1806 
1807  file_history_nitems = file_history_nitems + 1
1808  id = file_history_nitems
1809 
1810  file_history_vars(id)%name = file_history_req(reqid)%name
1811  file_history_vars(id)%outname = file_history_req(reqid)%outname
1812  file_history_vars(id)%basename = file_history_req(reqid)%basename
1813  file_history_vars(id)%postfix_timelabel = file_history_req(reqid)%postfix_timelabel
1814  file_history_vars(id)%zcoord = zcoord
1815  file_history_vars(id)%dstep = file_history_req(reqid)%dstep
1816  file_history_vars(id)%taverage = file_history_req(reqid)%taverage
1817  file_history_vars(id)%dtype = file_history_req(reqid)%dtype
1818 
1819  file_history_vars(id)%zid = -1
1820  do zid = 1, file_history_dims(dimid)%nzcoords
1821  if ( file_history_dims(dimid)%zcoords(zid) == file_history_vars(id)%zcoord ) then
1822  file_history_vars(id)%zid = zid
1823  exit
1824  end if
1825  end do
1826  if ( zid < 0 ) then
1827  log_error("FILE_HISTORY_Add_Variable",*) 'z-coordinate ', trim(file_history_vars(id)%zcoord), ' is not found for dimension ', trim(file_history_dims(dimid)%name)
1828  call prc_abort
1829  end if
1830 
1831  file_history_vars(id)%fid = -1
1832  file_history_vars(id)%vid = -1
1833  file_history_vars(id)%desc = desc
1834  file_history_vars(id)%units = units
1835  file_history_vars(id)%standard_name = standard_name
1836  file_history_vars(id)%dimid = dimid
1837  if ( present(cell_measures) ) then
1838  file_history_vars(id)%cell_measures = cell_measures
1839  else
1840  file_history_vars(id)%cell_measures = ""
1841  end if
1842  if ( present(fill_halo) ) then
1843  file_history_vars(id)%fill_halo = fill_halo
1844  else
1845  file_history_vars(id)%fill_halo = .false.
1846  end if
1847 
1848  file_history_vars(id)%waitstep = file_history_output_wait_step
1849  if ( file_history_output_step0 .AND. file_history_nowstep == 1 ) then
1850  file_history_vars(id)%laststep_write = 1 - file_history_vars(id)%dstep
1851  else
1852  file_history_vars(id)%laststep_write = 1
1853  endif
1854  file_history_vars(id)%laststep_put = file_history_vars(id)%laststep_write
1855  file_history_vars(id)%flag_clear = .true.
1856  file_history_vars(id)%size = file_history_dims(dimid)%size(zid)
1857  allocate( file_history_vars(id)%varsum( file_history_vars(id)%size ) )
1858 
1859  file_history_vars(id)%timesum = 0.0_dp
1860 
1861  if ( debug ) then
1862  log_info("FILE_HISTORY_Add_Variable",*) '[HISTORY] Item registration No.= ', id
1863  log_info_cont(*) 'Item name : ', trim(file_history_vars(id)%name)
1864  log_info_cont(*) 'Output name : ', trim(file_history_vars(id)%outname)
1865  log_info_cont(*) 'Description : ', trim(file_history_vars(id)%desc)
1866  log_info_cont(*) 'Unit : ', trim(file_history_vars(id)%units)
1867  log_info_cont(*) 'Basename of output file : ', trim(file_history_vars(id)%basename)
1868  log_info_cont(*) 'Add timelabel to the filename? : ', file_history_vars(id)%postfix_timelabel
1869  log_info_cont(*) 'Zcoord : ', trim(file_history_vars(id)%zcoord)
1870  log_info_cont(*) 'Interval [step] : ', file_history_vars(id)%dstep
1871  log_info_cont(*) 'Time Average? : ', file_history_vars(id)%taverage
1872  log_info_cont(*) 'Datatype : ', trim(file_dtypelist(file_history_vars(id)%dtype))
1873  log_info_cont(*) 'axis name : ', ( trim(file_history_dims(dimid)%dims(n,zid))//" ", n=1, file_history_dims(dimid)%ndims )
1874  endif
1875 
1876  existed = .false.
1877  do m = 1, file_history_nvar_inputs
1878  if ( file_history_var_inputs(m)%name == name ) then
1879  file_history_var_inputs(m)%nvariants = file_history_var_inputs(m)%nvariants + 1
1880  if ( file_history_var_inputs(m)%nvariants > file_history_variant_max ) then
1881  log_error("FILE_HISTORY_Add_Variable",*) 'Number of variant for ', trim(name), ' excees limit!'
1882  call prc_abort
1883  end if
1884  file_history_var_inputs(m)%variants(file_history_var_inputs(m)%nvariants) = id
1885  itemid = m
1886  existed = .true.
1887  exit
1888  end if
1889  end do
1890  if ( .not. existed ) then
1891  file_history_nvar_inputs = file_history_nvar_inputs + 1
1892  itemid = file_history_nvar_inputs
1893  file_history_var_inputs(itemid)%name = name
1894  file_history_var_inputs(itemid)%nvariants = 1
1895  file_history_var_inputs(itemid)%variants(1) = id
1896  end if
1897 
1898  endif ! match item?
1899 
1900  enddo
1901 
1902  return
1903  end subroutine file_history_add_variable
1904 
1905  !-----------------------------------------------------------------------------
1906  subroutine file_history_create( &
1907  id, &
1908  options )
1909  use scale_file_h, only: &
1910  file_real8, &
1911  file_real4
1912  use scale_file, only: &
1913  file_create, &
1914  file_flush, &
1915  file_redef, &
1916  file_set_option, &
1917  file_def_axis, &
1921  file_add_variable, &
1922  file_set_attribute
1923  use scale_time, only: &
1925  implicit none
1926 
1927  integer, intent(in) :: id
1928 
1929  character(len=*), intent(in) :: options ! 'filetype1:key1=val1&filetype2:key2=val2&...'
1930 
1931  integer :: fid
1932  character(len=H_MID) :: tunits
1933  character(len=H_LONG) :: basename_mod
1934  logical :: fileexisted
1935  integer(8) :: array_size
1936  integer :: dim_size
1937  integer :: dtype
1938  integer :: dimid, zid
1939  integer :: ndims
1940  character(len=H_SHORT) :: dims(3)
1941  real(dp) :: dtsec
1942 
1943  character(len=H_MID) :: timelabel
1944 
1945  integer :: ic, ie, is, lo
1946  integer :: m
1947  !---------------------------------------------------------------------------
1948 
1949  fid = file_history_vars(id)%fid
1950 
1951  if ( fid >= 0 ) return ! file already exists
1952 
1953  if ( file_history_time_since == '' ) then
1954  tunits = trim(file_history_time_units)
1955  else
1956  tunits = trim(file_history_time_units)//' since '//trim(file_history_time_since)
1957  endif
1958 
1959  if ( file_history_vars(id)%postfix_timelabel ) then
1960  call time_time2label( file_history_nowdate, file_history_nowsubsec, & ! [IN]
1961  timelabel ) ! [OUT]
1962  basename_mod = trim(file_history_vars(id)%basename)//'_'//trim(timelabel)
1963  else
1964  basename_mod = trim(file_history_vars(id)%basename)
1965  endif
1966 
1967  call file_create( basename_mod, & ! [IN]
1968  file_history_title, & ! [IN]
1969  file_history_source, & ! [IN]
1970  file_history_institution, & ! [IN]
1971  fid, fileexisted, & ! [OUT]
1972  rankid = file_history_myrank, & ! [IN]
1973  aggregate = file_history_aggregate, & ! [IN]
1974  time_units = tunits, & ! [IN]
1975  calendar = file_history_calendar ) ! [IN]
1976 
1977  if ( .not. fileexisted ) then
1978 
1979  ! write options
1980  ic = -1 ! index of ':'
1981  ie = -1 ! index of '='
1982  is = 1 ! start index
1983  lo = len_trim(options)
1984  if ( lo > 0 ) then
1985  do m = 1, lo+1
1986  if ( m == lo+1 .OR. options(m:m) == '&' ) then
1987  if ( ic == -1 .OR. ie == -1 ) then
1988  log_error("FILE_HISTORY_Create",*)'option is invalid: ', trim(options)
1989  call prc_abort
1990  endif
1991  call file_set_option( fid, options(is:ic-1), options(ic+1:ie-1), options(ie+1:m -1) ) ! [IN]
1992  ic = -1
1993  ie = -1
1994  is = m+1
1995  elseif( options(m:m) == ':' ) then
1996  ic = m
1997  elseif( options(m:m) == '=' ) then
1998  ie = m
1999  endif
2000  enddo
2001  endif
2002 
2003  if ( rp == dp ) then
2004  dtype = file_real8
2005  else
2006  dtype = file_real4
2007  end if
2008 
2009  ! define registered history axis variables in the newly created file
2010  ! actual writing axis variables are deferred to FILE_HISTORY_WriteAxes
2011  do m = 1, file_history_naxes
2012  if ( file_history_aggregate ) then ! for shared-file I/O, define axis in its global size
2013  dim_size = file_history_axes(m)%gdim_size ! axis global size
2014  if ( dim_size < 1 ) then
2015  log_error("FILE_HISTORY_Create",*) 'gsize is not set by FILE_HISTORY_Set_Axis'
2016  log_error_cont(*) 'It is necessary for aggregate file'
2017  call prc_abort
2018  end if
2019  else
2020  dim_size = file_history_axes(m)%dim_size
2021  endif
2022  call file_def_axis( fid, & ! [IN]
2023  file_history_axes(m)%name, & ! [IN]
2024  file_history_axes(m)%desc, & ! [IN]
2025  file_history_axes(m)%units, & ! [IN]
2026  file_history_axes(m)%dim, & ! [IN]
2027  dtype, dim_size, & ! [IN]
2028  bounds=associated(file_history_axes(m)%bounds) ) ! [IN]
2029  if ( file_history_axes(m)%down ) then
2030  call file_set_attribute( fid, file_history_axes(m)%name, 'positive', 'down' ) ! [IN]
2031  endif
2032  enddo
2033 
2034  ! define registered history associated coordinate variables in the newly created file
2035  ! actual writing coordinate variables are deferred to FILE_HISTORY_WriteAxes
2036  do m = 1, file_history_nassocs
2037  ndims = file_history_assocs(m)%ndims
2038  call file_def_associatedcoordinate( fid, & ! [IN]
2039  file_history_assocs(m)%name, & ! [IN]
2040  file_history_assocs(m)%desc, & ! [IN]
2041  file_history_assocs(m)%units, & ! [IN]
2042  file_history_assocs(m)%dims(1:ndims), & ! [IN]
2043  file_history_assocs(m)%dtype ) ! [IN]
2044  enddo
2045 
2046  ! attributes
2047  do m = 1, file_history_nattrs
2048 
2049  if ( file_history_attrs(m)%add_variable ) then
2050  ! associated variable
2051  call file_add_associatedvariable( fid, file_history_attrs(m)%varname )
2052  end if
2053 
2054  select case ( file_history_attrs(m)%type )
2055  case ( i_text )
2056  call file_set_attribute( fid, & ! [IN]
2057  file_history_attrs(m)%varname, & ! [IN]
2058  file_history_attrs(m)%key, & ! [IN]
2059  file_history_attrs(m)%text ) ! [IN]
2060  case ( i_int )
2061  call file_set_attribute( fid, & ! [IN]
2062  file_history_attrs(m)%varname, & ! [IN]
2063  file_history_attrs(m)%key, & ! [IN]
2064  file_history_attrs(m)%int(:) ) ! [IN]
2065  case ( i_float )
2066  call file_set_attribute( fid, & ! [IN]
2067  file_history_attrs(m)%varname, & ! [IN]
2068  file_history_attrs(m)%key, & ! [IN]
2069  file_history_attrs(m)%float(:) ) ! [IN]
2070  case ( i_double )
2071  call file_set_attribute( fid, & ! [IN]
2072  file_history_attrs(m)%varname, & ! [IN]
2073  file_history_attrs(m)%key, & ! [IN]
2074  file_history_attrs(m)%double(:) ) ! [IN]
2075  end select
2076 
2077  end do
2078 
2079  else
2080  call file_flush( fid )
2081  call file_redef( fid )
2082  end if
2083 
2084  do m = 1, file_history_nitems
2085  if ( file_history_vars(id)%basename == file_history_vars(m)%basename .and. &
2086  file_history_vars(m)%fid < 0 ) then
2087  ! Add new variable
2088  file_history_vars(m)%fid = fid
2089  dtsec = real(file_history_vars(m)%dstep,kind=dp) * file_history_dtsec
2090  dimid = file_history_vars(m)%dimid
2091  zid = file_history_vars(m)%zid
2092  ndims = file_history_dims(dimid)%ndims
2093  dims(1:ndims) = file_history_dims(dimid)%dims(1:ndims,zid)
2094  call file_add_variable( file_history_vars(m)%fid, & ! [IN]
2095  file_history_vars(m)%outname, & ! [IN]
2096  file_history_vars(m)%desc, & ! [IN]
2097  file_history_vars(m)%units, & ! [IN]
2098  file_history_vars(m)%standard_name, & ! [IN]
2099  dims(1:ndims), & ! [IN]
2100  file_history_vars(m)%dtype, & ! [IN]
2101  dtsec, & ! [IN]
2102  file_history_vars(m)%vid, & ! [OUT]
2103  time_avg=file_history_vars(m)%taverage ) ! [IN]
2104  if ( file_history_dims(dimid)%mapping /= "" ) then
2105  call file_set_attribute( file_history_vars(m)%fid, file_history_vars(m)%outname, & ! [IN]
2106  'grid_mapping', file_history_dims(dimid)%mapping ) ! [IN]
2107  endif
2108 
2109  select case( file_history_vars(m)%cell_measures )
2110  case ( "area", "area_z" )
2111  if ( file_history_dims(dimid)%area /= "" ) then
2112  call file_set_attribute( file_history_vars(m)%fid, file_history_vars(m)%outname, & ! [IN]
2113  'cell_measures', "area: "//trim(file_history_dims(dimid)%area) ) ! [IN]
2114  end if
2115  case ( "area_x" )
2116  if ( file_history_dims(dimid)%area_x /= "" ) then
2117  call file_set_attribute( file_history_vars(m)%fid, file_history_vars(m)%outname, & ! [IN]
2118  'cell_measures', "area: "//trim(file_history_dims(dimid)%area_x) ) ! [IN]
2119  end if
2120  case ( "area_y" )
2121  if ( file_history_dims(dimid)%area_x /= "" ) then
2122  call file_set_attribute( file_history_vars(m)%fid, file_history_vars(m)%outname, & ! [IN]
2123  'cell_measures', "area: "//trim(file_history_dims(dimid)%area_y) ) ! [IN]
2124  end if
2125  case ( "volume" )
2126  if ( file_history_dims(dimid)%area_x /= "" ) then
2127  call file_set_attribute( file_history_vars(m)%fid, file_history_vars(m)%outname, & ! [IN]
2128  'cell_measures', "volume: "//trim(file_history_dims(dimid)%volume) ) ! [IN]
2129  end if
2130  end select
2131 
2132  if ( file_history_dims(dimid)%location /= "" ) then
2133  if ( file_history_vars(m)%zcoord == "model" ) then
2134  call file_set_attribute( file_history_vars(m)%fid, file_history_vars(m)%outname, & ! [IN]
2135  'grid', file_history_dims(dimid)%grid ) ! [IN]
2136  else
2137  call file_set_attribute( file_history_vars(m)%fid, file_history_vars(m)%outname, & ! [IN]
2138  'grid', trim(file_history_dims(dimid)%grid)//'_'//trim(file_history_vars(id)%zcoord) ) ! [IN]
2139  end if
2140  call file_set_attribute( file_history_vars(m)%fid, file_history_vars(m)%outname, & ! [IN]
2141  'location', file_history_dims(dimid)%location ) ! [IN]
2142  end if
2143  end if
2144  end do
2145 
2146 
2147  ! allows PnetCDF to allocate an internal buffer
2148  ! to aggregate write requests for history variables
2149  array_size = 0
2150  do m = 1, file_history_nitems
2151  if ( file_history_vars(m)%fid == file_history_vars(id)%fid ) then
2152  array_size = array_size + file_history_vars(m)%size
2153  end if
2154  end do
2155  call file_attach_buffer( file_history_vars(id)%fid, array_size * dp )
2156 
2157 
2158  if ( .not. fileexisted ) call file_history_write_axes(id) ! [IN]
2159 
2160  return
2161  end subroutine file_history_create
2162 
2163  subroutine file_history_close
2164  use scale_file, only: &
2166  file_close
2167  implicit none
2168 
2169  integer :: fid, prev_fid
2170  integer :: id
2171  !---------------------------------------------------------------------------
2172 
2173  prev_fid = -1
2174  do id = 1, file_history_nitems
2175  fid = file_history_vars(id)%fid
2176  file_history_vars(id)%fid = -1
2177  if ( fid > 0 .AND. fid /= prev_fid ) then
2178  call file_detach_buffer( fid ) ! Release the internal buffer previously allowed to be used by PnetCDF
2179  call file_close( fid )
2180  prev_fid = fid
2181  endif
2182  enddo
2183 
2184  return
2185  end subroutine file_history_close
2186 
2187  !-----------------------------------------------------------------------------
2188  ! interface FILE_HISTORY_SetAssociatedCoordinate
2189  !-----------------------------------------------------------------------------
2191  name, &
2192  desc, &
2193  units, &
2194  dims, &
2195  var, &
2196  datatype, &
2197  start )
2198  use scale_file_h, only: &
2199  file_real4, &
2200  file_real8
2201  implicit none
2202 
2203  character(len=*), intent(in) :: name
2204  character(len=*), intent(in) :: desc
2205  character(len=*), intent(in) :: units
2206  character(len=*), intent(in) :: dims(:)
2207  real(RP), intent(in) :: var(:)
2208  character(len=*), intent(in), optional :: datatype
2209  integer, intent(in), optional :: start(:)
2210 
2211  integer :: dtype
2212  integer :: dim_size
2213  integer :: id
2214 
2215  intrinsic size, shape, reshape
2216  !---------------------------------------------------------------------------
2217 
2218  if ( present(datatype) ) then
2219  if ( datatype == 'REAL4' ) then
2220  dtype = file_real4
2221  elseif( datatype == 'REAL8' ) then
2222  dtype = file_real8
2223  else
2224  log_error("FILE_HISTORY_Set_AssociatedCoordinate_1D",*) 'Not appropriate datatype. Check!', datatype
2225  call prc_abort
2226  endif
2227  else if ( rp == sp ) then
2228  dtype = file_real4
2229  else
2230  dtype = file_real8
2231  endif
2232 
2233  dim_size = size(var)
2234 
2235  if ( file_history_nassocs < file_history_assoc_max ) then
2236  file_history_nassocs = file_history_nassocs + 1
2237  id = file_history_nassocs
2238 
2239  allocate( file_history_assocs(id)%var(dim_size) )
2240 
2241  file_history_assocs(id)%name = name
2242  file_history_assocs(id)%desc = desc
2243  file_history_assocs(id)%units = units
2244  file_history_assocs(id)%ndims = 1
2245  file_history_assocs(id)%dims(:) = ''
2246  file_history_assocs(id)%dims(1:1) = dims(1:1)
2247  file_history_assocs(id)%dtype = dtype
2248  file_history_assocs(id)%var(:) = real(reshape( var, (/ dim_size /) ),kind=dp)
2249 
2250  ! start and count are used for parallel I/O to a single shared file
2251  ! since var is reshaped into 1D array, we need to preserve its original shape in count
2252  file_history_assocs(id)%count(1:1) = shape(var)
2253  if ( present(start) ) then
2254  file_history_assocs(id)%start(1:1) = start(1:1)
2255  else
2256  file_history_assocs(id)%start = (/ 1, 1, 1, 1 /)
2257  end if
2258  else
2259  log_error("FILE_HISTORY_Set_AssociatedCoordinate_1D",*) 'Number of associate coordinates exceeds the limit.'
2260  call prc_abort
2261  endif
2262 
2263  return
2265 
2266  !-----------------------------------------------------------------------------
2268  name, &
2269  desc, &
2270  units, &
2271  dims, &
2272  var, &
2273  datatype, &
2274  start )
2275  use scale_file_h, only: &
2276  file_real4, &
2277  file_real8
2278  implicit none
2279 
2280  character(len=*), intent(in) :: name
2281  character(len=*), intent(in) :: desc
2282  character(len=*), intent(in) :: units
2283  character(len=*), intent(in) :: dims(:)
2284  real(RP), intent(in) :: var(:,:)
2285  character(len=*), intent(in), optional :: datatype
2286  integer, intent(in), optional :: start(:)
2287 
2288  integer :: dtype
2289  integer :: dim_size
2290  integer :: id
2291 
2292  intrinsic size, shape, reshape
2293  !---------------------------------------------------------------------------
2294 
2295  if ( present(datatype) ) then
2296  if ( datatype == 'REAL4' ) then
2297  dtype = file_real4
2298  elseif( datatype == 'REAL8' ) then
2299  dtype = file_real8
2300  else
2301  log_error("FILE_HISTORY_Set_AssociatedCoordinate_2D",*) 'Not appropriate datatype. Check!', datatype
2302  call prc_abort
2303  endif
2304  else if ( rp == sp ) then
2305  dtype = file_real4
2306  else
2307  dtype = file_real8
2308  endif
2309 
2310  dim_size = size(var)
2311 
2312  if ( file_history_nassocs < file_history_assoc_max ) then
2313  file_history_nassocs = file_history_nassocs + 1
2314  id = file_history_nassocs
2315 
2316  allocate( file_history_assocs(id)%var(dim_size) )
2317 
2318  file_history_assocs(id)%name = name
2319  file_history_assocs(id)%desc = desc
2320  file_history_assocs(id)%units = units
2321  file_history_assocs(id)%ndims = 2
2322  file_history_assocs(id)%dims(:) = ''
2323  file_history_assocs(id)%dims(1:2) = dims(1:2)
2324  file_history_assocs(id)%dtype = dtype
2325  file_history_assocs(id)%var(:) = real(reshape( var, (/ dim_size /) ),kind=dp)
2326 
2327  ! start and count are used for parallel I/O to a single shared file
2328  ! since var is reshaped into 1D array, we need to preserve its original shape in count
2329  file_history_assocs(id)%count(1:2) = shape(var)
2330  if ( present(start) ) then
2331  file_history_assocs(id)%start(1:2) = start(1:2)
2332  else
2333  file_history_assocs(id)%start = (/ 1, 1, 1, 1 /)
2334  end if
2335  else
2336  log_error("FILE_HISTORY_Set_AssociatedCoordinate_2D",*) 'Number of associate coordinates exceeds the limit.'
2337  call prc_abort
2338  endif
2339 
2340  return
2342 
2343  !-----------------------------------------------------------------------------
2345  name, &
2346  desc, &
2347  units, &
2348  dims, &
2349  var, &
2350  datatype, &
2351  start )
2352  use scale_file_h, only: &
2353  file_real4, &
2354  file_real8
2355  implicit none
2356 
2357  character(len=*), intent(in) :: name
2358  character(len=*), intent(in) :: desc
2359  character(len=*), intent(in) :: units
2360  character(len=*), intent(in) :: dims(:)
2361  real(RP), intent(in) :: var(:,:,:)
2362  character(len=*), intent(in), optional :: datatype
2363  integer, intent(in), optional :: start(:)
2364 
2365  integer :: dtype
2366  integer :: dim_size
2367  integer :: id
2368 
2369  intrinsic size, shape, reshape
2370  !---------------------------------------------------------------------------
2371 
2372  if ( present(datatype) ) then
2373  if ( datatype == 'REAL4' ) then
2374  dtype = file_real4
2375  elseif( datatype == 'REAL8' ) then
2376  dtype = file_real8
2377  else
2378  log_error("FILE_HISTORY_Set_AssociatedCoordinate_3D",*) 'Not appropriate datatype. Check!', datatype
2379  call prc_abort
2380  endif
2381  else if ( rp == sp ) then
2382  dtype = file_real4
2383  else
2384  dtype = file_real8
2385  endif
2386 
2387  dim_size = size(var)
2388 
2389  if ( file_history_nassocs < file_history_assoc_max ) then
2390  file_history_nassocs = file_history_nassocs + 1
2391  id = file_history_nassocs
2392 
2393  allocate( file_history_assocs(id)%var(dim_size) )
2394 
2395  file_history_assocs(id)%name = name
2396  file_history_assocs(id)%desc = desc
2397  file_history_assocs(id)%units = units
2398  file_history_assocs(id)%ndims = 3
2399  file_history_assocs(id)%dims(:) = ''
2400  file_history_assocs(id)%dims(1:3) = dims(1:3)
2401  file_history_assocs(id)%dtype = dtype
2402  file_history_assocs(id)%var(:) = real(reshape( var, (/ dim_size /) ),kind=dp)
2403 
2404  ! start and count are used for parallel I/O to a single shared file
2405  ! since var is reshaped into 1D array, we need to preserve its original shape in count
2406  file_history_assocs(id)%count(1:3) = shape(var)
2407  if ( present(start) ) then
2408  file_history_assocs(id)%start(1:3) = start(1:3)
2409  else
2410  file_history_assocs(id)%start = (/ 1, 1, 1, 1 /)
2411  end if
2412  else
2413  log_error("FILE_HISTORY_Set_AssociatedCoordinate_3D",*) 'Number of associate coordinates exceeds the limit.'
2414  call prc_abort
2415  endif
2416 
2417  return
2419 
2420  !-----------------------------------------------------------------------------
2421  ! interface FILE_HISTORY_Set_Attribute
2422  !-----------------------------------------------------------------------------
2423  subroutine file_history_set_attribute_text( &
2424  varname, &
2425  key, val, &
2426  add_variable )
2427  use scale_prc, only: &
2428  prc_abort
2429  use scale_file, only: &
2430  file_set_attribute
2431  implicit none
2432  character(len=*), intent(in) :: varname
2433  character(len=*), intent(in) :: key
2434  character(len=*), intent(in) :: val
2435  logical, intent(in), optional :: add_variable
2436 
2437  integer :: id
2438  !---------------------------------------------------------------------------
2439 
2440  file_history_nattrs = file_history_nattrs + 1
2441  if ( file_history_nattrs > file_history_attr_max ) then
2442  log_error("FILE_HISTORY_Set_Attribute_Text",*) 'number of attributes exceeds the limit'
2443  call prc_abort
2444  end if
2445 
2446  id = file_history_nattrs
2447 
2448  file_history_attrs(id)%varname = varname
2449  file_history_attrs(id)%key = key
2450  file_history_attrs(id)%text = val
2451  file_history_attrs(id)%type = i_text
2452 
2453  if ( present(add_variable) ) then
2454  file_history_attrs(id)%add_variable = add_variable
2455  else
2456  file_history_attrs(id)%add_variable = .false.
2457  end if
2458 
2459  return
2460  end subroutine file_history_set_attribute_text
2461 
2463  varname, &
2464  key, val, &
2465  add_variable )
2466  use scale_file, only: &
2467  file_set_attribute
2468  implicit none
2469  character(len=*), intent(in) :: varname
2470  character(len=*), intent(in) :: key
2471  logical, intent(in) :: val
2472  logical, intent(in), optional :: add_variable
2473 
2474  character(len=5) :: buf
2475  !---------------------------------------------------------------------------
2476 
2477  if ( val ) then
2478  buf = "true"
2479  else
2480  buf = "false"
2481  end if
2482 
2483  call file_history_set_attribute_text( varname, key, buf, add_variable=add_variable )
2484 
2485  return
2486  end subroutine file_history_set_attribute_logical
2487 
2488  !-----------------------------------------------------------------------------
2490  varname, &
2491  key, val, &
2492  add_variable )
2493  use scale_prc, only: &
2494  prc_abort
2495  use scale_file, only: &
2496  file_set_attribute
2497  implicit none
2498  character(len=*), intent(in) :: varname
2499  character(len=*), intent(in) :: key
2500  integer, intent(in) :: val(:)
2501  logical, intent(in), optional :: add_variable
2502 
2503  integer :: id
2504 
2505  intrinsic size
2506  !---------------------------------------------------------------------------
2507 
2508  file_history_nattrs = file_history_nattrs + 1
2509  if ( file_history_nattrs > file_history_attr_max ) then
2510  log_error("FILE_HISTORY_Set_Attribute_Int",*) 'number of attributes exceeds the limit'
2511  call prc_abort
2512  end if
2513 
2514  id = file_history_nattrs
2515 
2516  allocate( file_history_attrs(id)%int( size(val) ) )
2517 
2518  file_history_attrs(id)%varname = varname
2519  file_history_attrs(id)%key = key
2520  file_history_attrs(id)%int(:) = val(:)
2521  file_history_attrs(id)%type = i_int
2522 
2523  if ( present(add_variable) ) then
2524  file_history_attrs(id)%add_variable = add_variable
2525  else
2526  file_history_attrs(id)%add_variable = .false.
2527  end if
2528 
2529  return
2530  end subroutine file_history_set_attribute_int_ary
2531 
2532  !-----------------------------------------------------------------------------
2533  subroutine file_history_set_attribute_int( &
2534  varname, &
2535  key, val, &
2536  add_variable )
2537  implicit none
2538  character(len=*), intent(in) :: varname
2539  character(len=*), intent(in) :: key
2540  integer, intent(in) :: val
2541  logical, intent(in), optional :: add_variable
2542 
2543  integer :: ary(1)
2544  !---------------------------------------------------------------------------
2545 
2546  ary(1) = val
2547  call file_history_set_attribute_int_ary( varname, & ! (in)
2548  key, ary(:), & ! (in)
2549  add_variable=add_variable ) ! (in)
2550 
2551  return
2552  end subroutine file_history_set_attribute_int
2553 
2554  !-----------------------------------------------------------------------------
2555  subroutine file_history_set_attribute_float_ary( &
2556  varname, &
2557  key, val, &
2558  add_variable )
2559  use scale_prc, only: &
2560  prc_abort
2561  use scale_file, only: &
2562  file_set_attribute
2563  implicit none
2564  character(len=*), intent(in) :: varname
2565  character(len=*), intent(in) :: key
2566  real(SP), intent(in) :: val(:)
2567  logical, intent(in), optional :: add_variable
2568 
2569  integer :: id
2570 
2571  intrinsic size
2572  !---------------------------------------------------------------------------
2573 
2574  file_history_nattrs = file_history_nattrs + 1
2575  if ( file_history_nattrs > file_history_attr_max ) then
2576  log_error("FILE_HISTORY_Set_Attribute_Float",*) 'number of attributes exceeds the limit'
2577  call prc_abort
2578  end if
2579 
2580  id = file_history_nattrs
2581 
2582  allocate( file_history_attrs(id)%float( size(val) ) )
2583 
2584  file_history_attrs(id)%varname = varname
2585  file_history_attrs(id)%key = key
2586  file_history_attrs(id)%float(:) = val(:)
2587  file_history_attrs(id)%type = i_float
2588 
2589  if ( present(add_variable) ) then
2590  file_history_attrs(id)%add_variable = add_variable
2591  else
2592  file_history_attrs(id)%add_variable = .false.
2593  end if
2594 
2595  return
2596  end subroutine file_history_set_attribute_float_ary
2597 
2598  subroutine file_history_set_attribute_float( &
2599  varname, &
2600  key, val, &
2601  add_variable )
2602  implicit none
2603  character(len=*), intent(in) :: varname
2604  character(len=*), intent(in) :: key
2605  real(SP), intent(in) :: val
2606  logical, intent(in), optional :: add_variable
2607 
2608  real(SP) :: ary(1)
2609 
2610  ary(1) = val
2611  call file_history_set_attribute_float_ary( varname, & ! (in)
2612  key, ary(:), & ! (in)
2613  add_variable ) ! (in)
2614 
2615  return
2616  end subroutine file_history_set_attribute_float
2617  !-----------------------------------------------------------------------------
2618  subroutine file_history_set_attribute_double_ary( &
2619  varname, &
2620  key, val, &
2621  add_variable )
2622  use scale_prc, only: &
2623  prc_abort
2624  use scale_file, only: &
2625  file_set_attribute
2626  implicit none
2627  character(len=*), intent(in) :: varname
2628  character(len=*), intent(in) :: key
2629  real(DP), intent(in) :: val(:)
2630  logical, intent(in), optional :: add_variable
2631 
2632  integer :: id
2633 
2634  intrinsic size
2635  !---------------------------------------------------------------------------
2636 
2637  file_history_nattrs = file_history_nattrs + 1
2638  if ( file_history_nattrs > file_history_attr_max ) then
2639  log_error("FILE_HISTORY_Set_Attribute_Double",*) 'number of attributes exceeds the limit'
2640  call prc_abort
2641  end if
2642 
2643  id = file_history_nattrs
2644 
2645  allocate( file_history_attrs(id)%double( size(val) ) )
2646 
2647  file_history_attrs(id)%varname = varname
2648  file_history_attrs(id)%key = key
2649  file_history_attrs(id)%double(:) = val(:)
2650  file_history_attrs(id)%type = i_double
2651 
2652  if ( present(add_variable) ) then
2653  file_history_attrs(id)%add_variable = add_variable
2654  else
2655  file_history_attrs(id)%add_variable = .false.
2656  end if
2657 
2658  return
2659  end subroutine file_history_set_attribute_double_ary
2660 
2661  subroutine file_history_set_attribute_double( &
2662  varname, &
2663  key, val, &
2664  add_variable )
2665  implicit none
2666  character(len=*), intent(in) :: varname
2667  character(len=*), intent(in) :: key
2668  real(DP), intent(in) :: val
2669  logical, intent(in), optional :: add_variable
2670 
2671  real(DP) :: ary(1)
2672 
2673  ary(1) = val
2674  call file_history_set_attribute_double_ary( varname, & ! (in)
2675  key, ary(:), & ! (in)
2676  add_variable ) ! (in)
2677 
2678  return
2679  end subroutine file_history_set_attribute_double
2680 
2681  !-----------------------------------------------------------------------------
2682  ! interface FILE_HOSTORY_Query
2683  !-----------------------------------------------------------------------------
2684  subroutine file_history_query_id( &
2685  itemid, &
2686  answer )
2687  integer, intent(in) :: itemid
2688  logical, intent(out) :: answer
2689 
2690  integer :: id, i
2691 
2692  answer = .false.
2693  if ( file_history_disabled ) return
2694  if ( itemid < 0 ) return
2695 
2696  do i = 1, file_history_var_inputs(itemid)%nvariants
2697  id = file_history_var_inputs(itemid)%variants(i)
2698  if ( file_history_vars(id)%taverage ) then
2699  answer = .true.
2700  return
2701  else if ( file_history_nowstep >= file_history_vars(id)%laststep_write + file_history_vars(id)%dstep ) then
2702  answer = .true.
2703  return
2704  endif
2705  end do
2706 
2707  return
2708  end subroutine file_history_query_id
2709 
2710  subroutine file_history_query_name( &
2711  name, &
2712  answer )
2713  implicit none
2714 
2715  character(len=*), intent(in) :: name
2716 
2717  logical, intent(out) :: answer
2718 
2719  integer :: itemid
2720  !---------------------------------------------------------------------------
2721 
2722  answer = .false.
2723  if ( file_history_disabled ) return
2724 
2725  do itemid = 1, file_history_nvar_inputs
2726  if ( file_history_var_inputs(itemid)%name == name ) then
2727  call file_history_query_id( itemid, answer ) ! [IN], [OUT]
2728  return
2729  end if
2730  end do
2731 
2732  return
2733  end subroutine file_history_query_name
2734 
2735  !-----------------------------------------------------------------------------
2736  subroutine file_history_write_axes(id)
2737  use scale_file, only: &
2738  file_enddef, &
2739  file_flush, &
2740  file_write_axis, &
2741  file_write_associatedcoordinate
2742  use scale_prc, only: &
2743  prc_abort
2744  implicit none
2745  integer, intent(in) :: id
2746 
2747  integer :: start(1)
2748  integer :: m, fid
2749  !---------------------------------------------------------------------------
2750 
2751  if ( id < 0 ) return
2752 
2753  fid = file_history_vars(id)%fid
2754  call file_enddef( fid )
2755 
2756  ! write registered history variables to file
2757  do m = 1, file_history_naxes
2758  if ( file_history_axes(m)%start > 0 ) then
2759  start(1) = file_history_axes(m)%start
2760 
2761  call file_write_axis( fid, & ! [IN]
2762  file_history_axes(m)%name, & ! [IN]
2763  file_history_axes(m)%var, & ! [IN]
2764  start ) ! [IN]
2765 
2766  if ( associated(file_history_axes(m)%bounds) ) then
2767  call file_write_associatedcoordinate( fid, & ! [IN]
2768  trim(file_history_axes(m)%name)//'_bnds', & ! [IN]
2769  file_history_axes(m)%bounds(:,:), & ! [IN]
2770  (/ 1, start(1) /) ) ! [IN]
2771  end if
2772 
2773  end if
2774  end do
2775 
2776  do m = 1, file_history_nassocs
2777  call file_write_associatedcoordinate( fid, & ! [IN]
2778  file_history_assocs(m)%name, & ! [IN]
2779  file_history_assocs(m)%var, & ! [IN]
2780  file_history_assocs(m)%start, & ! [IN]
2781  file_history_assocs(m)%count, & ! [IN]
2782  file_history_assocs(m)%ndims ) ! [IN]
2783  enddo
2784 
2785  ! for PnetCDF I/O, flush all pending nonblocking write requests
2786  call file_flush( fid )
2787 
2788  return
2789  end subroutine file_history_write_axes
2790 
2791  !-----------------------------------------------------------------------------
2792  subroutine file_history_write_onevar( &
2793  id, &
2794  step_now )
2795  use scale_file_h, only: &
2796  rmiss => file_rmiss
2797  use scale_calendar, only: &
2799  use scale_file, only: &
2800  file_write
2801  implicit none
2802 
2803  integer, intent(in) :: id
2804  integer, intent(in) :: step_now
2805 
2806  integer :: dimid, zid
2807  real(DP) :: time_str, time_end
2808  real(DP) :: sec_str, sec_end
2809  integer :: i
2810  !---------------------------------------------------------------------------
2811 
2812  if( file_history_nreqs == 0 ) return
2813 
2814  if ( step_now < file_history_vars(id)%laststep_write + file_history_vars(id)%dstep ) then
2815  return
2816  endif
2817 
2818  if ( file_history_vars(id)%flag_clear ) then
2819  if ( file_history_output_step0 .AND. file_history_nowstep == 1 ) then
2820  do i = 1, file_history_vars(id)%size
2821  file_history_vars(id)%varsum(i) = rmiss
2822  end do
2823  else if ( file_history_error_putmiss ) then
2824  log_error("FILE_HISTORY_Write_OneVar",*) 'The time interval of history output ', trim(file_history_vars(id)%name), &
2825  ' and the time interval of its related scheme are inconsistent.'
2826  log_error_cont(*) 'Please check the namelist PARAM_TIME, PARAM_FILE_HISTORY, and HISTORY_ITEM.'
2827  log_error_cont(*) 'Please set FILE_HISTORY_ERROR_PUTMISS in the namelist PARAM_FILE_HISTORY to .false.', &
2828  ' when you want to disable this check.'
2829  log_error_cont(*) 'The time interval of history output ', trim(file_history_vars(id)%name), &
2830  ' and the time interval of its related scheme are inconsistent.', &
2831  ' Please see detail in log file.'
2832  call prc_abort
2833  else
2834  log_warn("FILE_HISTORY_Write_OneVar",*) 'Output value is not updated in this step.', &
2835  ' NAME = ', trim(file_history_vars(id)%name), &
2836  ', OUTNAME = ', trim(file_history_vars(id)%outname)
2837  endif
2838  endif
2839 
2840  if ( .NOT. file_history_vars(id)%flag_clear .AND. file_history_vars(id)%taverage ) then
2841  do i = 1, file_history_vars(id)%size
2842  if ( file_history_vars(id)%varsum(i) /= rmiss ) then
2843  file_history_vars(id)%varsum(i) = file_history_vars(id)%varsum(i) / file_history_vars(id)%timesum
2844  end if
2845  end do
2846  endif
2847 
2848  call file_history_output_list
2849 
2850  if ( step_now > file_history_vars(id)%waitstep ) then
2851  if ( laststep_write < step_now ) then ! log only once in this step
2852  log_progress(*) 'output history'
2853  endif
2854 
2855  ! Note this subroutine must be called after all FILE_HISTORY_reg calls are completed
2856  ! Write registered history axes to history file
2857  call file_history_create( id, options = file_history_options ) ! [IN]
2858 
2859  sec_str = file_history_startdaysec + real(file_history_vars(id)%laststep_write-1,kind=dp) * file_history_dtsec
2860  sec_end = file_history_startdaysec + real(step_now -1,kind=dp) * file_history_dtsec
2861 
2862  ! convert time units [sec]->[sec,min,hour,day,month,year]
2863  call calendar_sec2unit( time_str, sec_str, file_history_time_units )
2864  call calendar_sec2unit( time_end, sec_end, file_history_time_units )
2865 
2866  dimid = file_history_vars(id)%dimid
2867  zid = file_history_vars(id)%zid
2868  if ( file_history_dims(dimid)%count(1,zid) > 0 ) then
2869 
2870  ! for one-file-per-process I/O method, count(1) == 1 always
2871  ! for one file shared by all processes, count(1) >= 0,
2872  ! being 0 indicates a 1D history variable, which will only be written by the
2873  ! south-most processes in parallel, or a z axis to be written by rank 0 only
2874 
2875  call file_write( file_history_vars(id)%vid, & ! [IN]
2876  file_history_vars(id)%varsum(:), & ! [IN]
2877  time_str, & ! [IN]
2878  time_end, & ! [IN]
2879  ndims=file_history_dims(dimid)%ndims, & ! ndims before reshape
2880  count=file_history_dims(dimid)%count(:,zid), & ! global subarray lengths
2881  start=file_history_dims(dimid)%start(:,zid) ) ! global subarray start indices
2882  end if
2883  else
2884  if ( laststep_write < step_now ) then
2885  log_progress(*) 'history output is suppressed'
2886  endif
2887  endif
2888 
2889  file_history_vars(id)%laststep_write = step_now
2890  file_history_vars(id)%flag_clear = .true.
2891 
2892  laststep_write = step_now ! remember for multiple call in the same step
2893 
2894  return
2895  end subroutine file_history_write_onevar
2896 
2897  !-----------------------------------------------------------------------------
2898  subroutine file_history_output_list
2899  implicit none
2900 
2901  real(DP) :: dtsec
2902  integer :: id
2903  !---------------------------------------------------------------------------
2904 
2905  if ( list_outputed ) then
2906  return
2907  endif
2908 
2909  if ( file_history_nitems /= file_history_nreqs ) then
2910 
2911  if ( .not. ( file_history_output_step0 .and. file_history_nowstep == 1 ) ) then
2912 
2913  log_info("FILE_HISTORY_Output_List",*) '[HISTORY] All of requested variable by the namelist HISTORY_ITEM did not find.'
2914  do id = 1, file_history_nreqs
2915  log_info("FILE_HISTORY_Output_List",'(A,A24,A,L1)') 'NAME : ', file_history_req(id)%name, &
2916  ', registered? : ', file_history_req(id)%registered
2917  enddo
2918  log_info("FILE_HISTORY_Output_List",*) 'Please set FILE_HISTORY_ERROR_PUTMISS in the namelist PARAM_FILE_HISTORY to .false.', &
2919  ' when you want to disable this check.'
2920 
2921  if ( file_history_error_putmiss ) then
2922  log_error("FILE_HISTORY_Output_List",*) 'Requested variables by the namelist HISTORY_ITEM did not find. Please see detail in log file.'
2923  call prc_abort
2924  endif
2925  end if
2926 
2927  return
2928  end if
2929 
2930  log_info("FILE_HISTORY_Output_List",*) '[HISTORY] Output item list '
2931  log_info_cont('(1x,A,I4)') 'Number of history item :', file_history_nreqs
2932  log_info_cont(*) 'ITEM :OUTNAME ', &
2933  ': size:interval[sec]: step:timeavg?:zcoord'
2934  log_info_cont(*) '=================================================', &
2935  '================================================='
2936 
2937 
2938  do id = 1, file_history_nitems
2939  dtsec = real(file_history_vars(id)%dstep,kind=dp) * file_history_dtsec
2940 
2941  log_info_cont('(1x,A24,1x,A24,1x,I8,1x,F13.3,1x,I8,1x,L8,1x,A8)') &
2942  file_history_vars(id)%name, &
2943  file_history_vars(id)%outname, &
2944  file_history_vars(id)%size, &
2945  dtsec, &
2946  file_history_vars(id)%dstep, &
2947  file_history_vars(id)%taverage, &
2948  file_history_vars(id)%zcoord
2949  enddo
2950 
2951  log_info_cont(*) '=================================================', &
2952  '================================================='
2953 
2954  list_outputed = .true.
2955 
2956  return
2957  end subroutine file_history_output_list
2958 
2959  function file_history_find_id( name )
2960  character(len=*), intent(in) :: name
2961  integer :: FILE_HISTORY_find_id
2962 
2963  integer :: itemid
2964 
2965  do itemid = 1, file_history_nvar_inputs
2966  if ( file_history_var_inputs(itemid)%name == name ) then
2967  file_history_find_id = itemid
2968  return
2969  end if
2970  end do
2971 
2972  file_history_find_id = -1
2973 
2974  return
2975  end function file_history_find_id
2976 
2977  function file_history_get_size( &
2978  dims, ndims )
2979  character(len=*), intent(in) :: dims(:)
2980  integer, intent(in) :: ndims
2981  integer :: FILE_HISTORY_get_size
2982 
2983  integer :: len
2984  integer :: n, i
2985 
2986  file_history_get_size = 1
2987  do n = 1, ndims
2988  len = -1
2989  do i = 1, file_history_naxes
2990  if ( file_history_axes(i)%name == dims(n) ) then
2991  len = file_history_axes(i)%dim_size
2992  exit
2993  end if
2994  end do
2995  if ( len < 0 ) then
2996  log_error("FILE_HISTORY_get_size",*) 'dimension name is not found: ', dims(n)
2997  call prc_abort
2998  end if
2999  file_history_get_size = file_history_get_size * len
3000  end do
3001 
3002  return
3003  end function file_history_get_size
3004 
3005  subroutine file_history_truncate_1d_default( &
3006  src, &
3007  dim_type, zcoord, fill_halo, &
3008  dsc )
3009  real(RP), intent(in) :: src(:)
3010  character(len=*), intent(in) :: dim_type
3011  character(len=*), intent(in) :: zcoord
3012  logical, intent(in) :: fill_halo
3013  real(DP), intent(out) :: dsc(:)
3014 
3015  dsc(:) = src(:)
3016 
3017  return
3018  end subroutine file_history_truncate_1d_default
3019  subroutine file_history_truncate_2d_default( &
3020  src, &
3021  dim_type, zcoord, fill_halo, &
3022  dsc )
3023  real(RP), intent(in) :: src(:,:)
3024  character(len=*), intent(in) :: dim_type
3025  character(len=*), intent(in) :: zcoord
3026  logical, intent(in) :: fill_halo
3027  real(DP), intent(out) :: dsc(:)
3028 
3029  integer :: i, j
3030  integer :: idx
3031 
3032  intrinsic size
3033 
3034  idx = 1
3035  do j = 1, size(src,2)
3036  do i = 1, size(src,1)
3037  dsc(idx) = src(i, j)
3038  idx = idx + 1
3039  end do
3040  end do
3041 
3042  return
3043  end subroutine file_history_truncate_2d_default
3044  subroutine file_history_truncate_3d_default( &
3045  src, &
3046  dim_type, zcoord, fill_halo, &
3047  dsc )
3048  real(RP), intent(in) :: src(:,:,:)
3049  character(len=*), intent(in) :: dim_type
3050  character(len=*), intent(in) :: zcoord
3051  logical, intent(in) :: fill_halo
3052  real(DP), intent(out) :: dsc(:)
3053 
3054  integer :: k, i, j
3055  integer :: idx
3056 
3057  intrinsic size
3058 
3059  idx = 1
3060  do j = 1, size(src,3)
3061  do i = 1, size(src,2)
3062  do k = 1, size(src,1)
3063  dsc(idx) = src(k, i, j)
3064  idx = idx + 1
3065  end do
3066  end do
3067  end do
3068 
3069  return
3070  end subroutine file_history_truncate_3d_default
3071  subroutine file_history_truncate_4d_default( &
3072  src, &
3073  dim_type, zcoord, fill_halo, &
3074  dsc )
3075  real(RP), intent(in) :: src(:,:,:,:)
3076  character(len=*), intent(in) :: dim_type
3077  character(len=*), intent(in) :: zcoord
3078  logical, intent(in) :: fill_halo
3079  real(DP), intent(out) :: dsc(:)
3080 
3081  integer :: l, k, i, j
3082  integer :: idx
3083 
3084  intrinsic size
3085 
3086  idx = 1
3087  do j = 1, size(src,4)
3088  do i = 1, size(src,3)
3089  do k = 1, size(src,2)
3090  do l = 1, size(src,1)
3091  dsc(idx) = src(l, k, i, j)
3092  idx = idx + 1
3093  end do
3094  end do
3095  end do
3096  end do
3097 
3098  return
3099  end subroutine file_history_truncate_4d_default
3100 
3101 end module scale_file_history
3102 
3103 
3104 
3105 !--
3106 ! vi:set readonly sw=4 ts=8
3107 !
3108 !Local Variables:
3109 !mode: f90
3110 !buffer-read-only:t
3111 !End:
3112 !
3113 !++
scale_file_history::file_history_set_attribute_double
subroutine file_history_set_attribute_double(varname, key, val, add_variable)
Definition: scale_file_history.F90:2665
scale_file_history::file_history_set_attribute_int_ary
subroutine file_history_set_attribute_int_ary(varname, key, val, add_variable)
Definition: scale_file_history.F90:2493
scale_precision::sp
integer, parameter, public sp
Definition: scale_precision.F90:31
scale_prc::prc_abort
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:342
scale_file_history::file_history_set_attribute_logical
subroutine file_history_set_attribute_logical(varname, key, val, add_variable)
Definition: scale_file_history.F90:2466
scale_file_history::file_history_in_3d
subroutine file_history_in_3d(var, name, desc, unit, standard_name, dim_type, fill_halo)
Wrapper routine of FILE_HISTORY_reg + FILE_HISTORY_put.
Definition: scale_file_history.F90:1289
scale_file_history::file_history_finalize
subroutine, public file_history_finalize
finalization
Definition: scale_file_history.F90:1726
scale_file_history::i
logical, public i
Definition: scale_file_history.F90:143
scale_file_h::file_real4
integer, parameter, public file_real4
Definition: scale_file_h.F90:24
scale_file::file_enddef
subroutine, public file_enddef(fid)
Definition: scale_file.F90:4585
scale_precision
module PRECISION
Definition: scale_precision.F90:14
scale_file_history::aggregate
logical, public aggregate
Definition: scale_file_history.F90:143
scale_file::file_attach_buffer
subroutine, public file_attach_buffer(fid, buf_amount)
Definition: scale_file.F90:4643
scale_file::file_flush
subroutine, public file_flush(fid)
Definition: scale_file.F90:4709
scale_file_h::file_dtypelist
character(len=file_hshort), dimension(0:4), public file_dtypelist
Definition: scale_file_h.F90:75
scale_file_history::file_history_set_disable
subroutine, public file_history_set_disable(switch)
set switch to turn on/off history
Definition: scale_file_history.F90:1712
scale_const::const_eps
real(rp), public const_eps
small number
Definition: scale_const.F90:33
scale_prof::prof_rapstart
subroutine, public prof_rapstart(rapname_base, level, disable_barrier)
Start raptime.
Definition: scale_prof.F90:159
scale_calendar::calendar_sec2unit
subroutine, public calendar_sec2unit(value, second, unit)
Convert several second to specified unit.
Definition: scale_calendar.F90:458
scale_file_h::file_rmiss
real(dp), parameter, public file_rmiss
Definition: scale_file_h.F90:49
scale_calendar
module CALENDAR
Definition: scale_calendar.F90:13
scale_file_history::file_history_truncate_3d
procedure(truncate_3d), pointer, public file_history_truncate_3d
Definition: scale_file_history.F90:101
scale_file_history::file_history_set_associatedcoordinate_3d
subroutine file_history_set_associatedcoordinate_3d(name, desc, units, dims, var, datatype, start)
Definition: scale_file_history.F90:2352
scale_file_history
module file_history
Definition: scale_file_history.F90:15
float
typedef float(real32_t)
scale_file
module file
Definition: scale_file.F90:15
scale_file_history::file_history_set_attribute_float
subroutine file_history_set_attribute_float(varname, key, val, add_variable)
Definition: scale_file_history.F90:2602
scale_file_history::file_history_set_nowdate
subroutine, public file_history_set_nowdate(NOWDATE, NOWSUBSEC, NOWSTEP)
set now step
Definition: scale_file_history.F90:1697
scale_file::file_def_axis
subroutine, public file_def_axis(fid, name, desc, units, dim_name, dtype, dim_size, bounds)
Definition: scale_file.F90:667
scale_prc
module PROCESS
Definition: scale_prc.F90:11
scale_precision::rp
integer, parameter, public rp
Definition: scale_precision.F90:41
scale_file_history::file_history_set_dim
subroutine, public file_history_set_dim(name, ndims, nzcoords, dims, zcoords, start, count, mapping, area, area_x, area_y, volume, location, grid)
set dimension information
Definition: scale_file_history.F90:1470
scale_file_history::file_history_set_attribute_int
subroutine file_history_set_attribute_int(varname, key, val, add_variable)
Definition: scale_file_history.F90:2537
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_file_history::file_history_set_axis
subroutine, public file_history_set_axis(name, desc, units, dim, var, bounds, down, gsize, start)
set axis information
Definition: scale_file_history.F90:1579
scale_file::file_close
subroutine, public file_close(fid, abort)
Definition: scale_file.F90:4738
scale_file::file_def_associatedcoordinate
subroutine, public file_def_associatedcoordinate(fid, name, desc, units, dim_names, dtype)
Definition: scale_file.F90:1037
scale_file_history::switch
logical, public switch
Definition: scale_file_history.F90:143
scale_const
module CONSTANT
Definition: scale_const.F90:11
scale_file_history::file_history_in_4d
subroutine file_history_in_4d(var, name, desc, unit, standard_name, dim_type, fill_halo)
Wrapper routine of FILE_HISTORY_reg + FILE_HISTORY_put.
Definition: scale_file_history.F90:1421
scale_file_history::file_history_truncate_1d
procedure(truncate_1d), pointer, public file_history_truncate_1d
Definition: scale_file_history.F90:77
scale_file_history::file_history_setup
subroutine, public file_history_setup(title, source, institution, time_start, time_interval, time_units, time_since, calendar, default_basename, default_postfix_timelabel, default_zcoord, default_tinterval, default_tunit, default_taverage, default_datatype, myrank)
Setup.
Definition: scale_file_history.F90:342
scale_io::h_short
integer, parameter, public h_short
Character length (short=16)
Definition: scale_io.F90:44
scale_file_history::file_history_set_attribute_text
subroutine file_history_set_attribute_text(varname, key, val, add_variable)
Definition: scale_file_history.F90:2427
scale_prof
module profiler
Definition: scale_prof.F90:11
scale_precision::dp
integer, parameter, public dp
Definition: scale_precision.F90:32
scale_file_h
module file_h
Definition: scale_file_h.F90:11
scale_time
module TIME
Definition: scale_time.F90:11
scale_file::file_add_associatedvariable
subroutine, public file_add_associatedvariable(fid, vname, existed)
Definition: scale_file.F90:422
scale_file_history::file_history_in_2d
subroutine file_history_in_2d(var, name, desc, unit, standard_name, dim_type, fill_halo)
Wrapper routine of FILE_HISTORY_reg + FILE_HISTORY_put.
Definition: scale_file_history.F90:1157
scale_file_history::file_history_set_associatedcoordinate_1d
subroutine file_history_set_associatedcoordinate_1d(name, desc, units, dims, var, datatype, start)
Definition: scale_file_history.F90:2198
scale_file::file_aggregate
logical, public file_aggregate
Definition: scale_file.F90:182
scale_file_history::file_history_write
subroutine, public file_history_write
Definition: scale_file_history.F90:1643
scale_file_history::file_history_aggregate
logical, public file_history_aggregate
Definition: scale_file_history.F90:143
scale_file::file_redef
subroutine, public file_redef(fid)
Definition: scale_file.F90:4613
scale_file_history::file_history_reg
subroutine, public file_history_reg(name, desc, unit, itemid, standard_name, ndims, dim_type, cell_measures, fill_halo)
Register/Append variable to history file.
Definition: scale_file_history.F90:650
double
typedef double(real64_t)
scale_prof::prof_rapend
subroutine, public prof_rapend(rapname_base, level, disable_barrier)
Save raptime.
Definition: scale_prof.F90:217
scale_file_history::file_history_in_0d
subroutine file_history_in_0d(var, name, desc, unit, standard_name, dim_type)
Wrapper routine of FILE_HISTORY_reg + FILE_HISTORY_put.
Definition: scale_file_history.F90:892
scale_file::file_detach_buffer
subroutine, public file_detach_buffer(fid)
Definition: scale_file.F90:4677
scale_file_history::file_history_set_associatedcoordinate_2d
subroutine file_history_set_associatedcoordinate_2d(name, desc, units, dims, var, datatype, start)
Definition: scale_file_history.F90:2275
scale_time::time_time2label
subroutine, public time_time2label(date, subsec, timelabel)
generate time label
Definition: scale_time.F90:107
scale_const::const_undef
real(rp), public const_undef
Definition: scale_const.F90:41
scale_file_history::file_history_in_1d
subroutine file_history_in_1d(var, name, desc, unit, standard_name, dim_type)
Wrapper routine of FILE_HISTORY_reg + FILE_HISTORY_put.
Definition: scale_file_history.F90:1024
scale_io::io_fid_conf
integer, public io_fid_conf
Config file ID.
Definition: scale_io.F90:56
scale_file_h::file_real8
integer, parameter, public file_real8
Definition: scale_file_h.F90:25
scale_file::file_create
subroutine, public file_create(basename, title, source, institution, fid, existed, rankid, single, aggregate, time_units, calendar, append)
create file fid is >= 1
Definition: scale_file.F90:267
scale_file::file_set_option
subroutine, public file_set_option(fid, filetype, key, val)
Definition: scale_file.F90:456
scale_file_history::file_history_truncate_4d
procedure(truncate_4d), pointer, public file_history_truncate_4d
Definition: scale_file_history.F90:113
scale_calendar::calendar_unit2sec
subroutine, public calendar_unit2sec(second, value, unit)
Convert several units to second.
Definition: scale_calendar.F90:424