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