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