FMS  2025.04
Flexible Modeling System
diag_integral.F90
1 !***********************************************************************
2 !* Apache License 2.0
3 !*
4 !* This file is part of the GFDL Flexible Modeling System (FMS).
5 !*
6 !* Licensed under the Apache License, Version 2.0 (the "License");
7 !* you may not use this file except in compliance with the License.
8 !* You may obtain a copy of the License at
9 !*
10 !* http://www.apache.org/licenses/LICENSE-2.0
11 !*
12 !* FMS is distributed in the hope that it will be useful, but WITHOUT
13 !* WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied;
14 !* without even the implied warranty of MERCHANTABILITY or FITNESS FOR A
15 !* PARTICULAR PURPOSE. See the License for the specific language
16 !* governing permissions and limitations under the License.
17 !***********************************************************************
18 !> @defgroup diag_integral_mod diag_integral_mod
19 !> @ingroup diag_integral
20 !!
21 !! @author Fei Liu <Fei.Liu@noaa.gov>
22 !!
23 !! @brief This module computes and outputs global and / or hemispheric physics
24 !! integrals.
25 
26 module diag_integral_mod
27 
28 !###############################################################################
29 
30 use platform_mod, only: i8_kind, fms_file_len
31 use time_manager_mod, only: time_type, get_time, set_time, &
33  operator(+), operator(-), &
34  operator(==), operator(>=), &
35  operator(/=)
36 use mpp_mod, only: input_nml_file
37 use fms_mod, only: error_mesg, &
39  fms_init, &
40  mpp_pe, mpp_root_pe,&
41  fatal, write_version_number, &
42  stdlog, string
43 use fms2_io_mod, only: file_exists
44 use constants_mod, only: radius, constants_init
45 use mpp_mod, only: mpp_sum, mpp_init
46 use ensemble_manager_mod, only : get_ensemble_id, get_ensemble_size
47 use platform_mod, only: r4_kind, r8_kind
48 
49 !-------------------------------------------------------------------------------
50 
51 implicit none
52 private
53 
54 !-------------------------------------------------------------------------------
55 ! version number of this module
56 ! Include variable "version" to be written to log file.
57 #include<file_version.h>
58 
59 !-------------------------------------------------------------------------------
60 !------ interfaces ------
61 
62 public &
66 
67 
68 
69 !###############################################################################
70 !> Perform a summation of the named integral field
71 !!
72 !> This interface can be called in any one of three ways:
73 !!
74 !! @code{.f90}
75 !! call sum_diag_integral_field (name, data, is, js)
76 !! call sum_diag_integral_field (name, data, wt, is, js)
77 !! call sum_diag_integral_field (name, data, is, ie, js, je)
78 !! @endcode
79 !!
80 !! in the first option above, `data` may be either
81 !! @code{.f90}
82 !! real, intent(in) :: data(:,:) ![ sum_field_2d ]
83 !! real, intent(in) :: data(:,:,:) ![ sum_field_3d ]
84 !! @endcode
85 !!
86 !! <b> Parameters: </b>
87 !!
88 !! @code{.f90}
89 !! character(len=*), intent(in) :: name
90 !! real, intent(in) :: wt(:,:,:)
91 !! integer, optional, intent(in) :: is, ie, js, je
92 !! @endcode
93 !!
94 !! @param [in] <name> name associated with integral
95 !! @param [in] <data> field of integrands to be summed over
96 !! @param [in] <wt> vertical weighting factor to be applied to integrands
97 !! when summing
98 !! @param [in] <is, ie, js, je> starting/ending i,j indices over which summation
99 !! is to occur
100 !!
101 !! @ingroup diag_integral_mod
103  module procedure sum_field_2d_r4, sum_field_2d_r8
104  module procedure sum_field_2d_hemi_r4, sum_field_2d_hemi_r8
105  module procedure sum_field_3d_r4, sum_field_3d_r8
106  module procedure sum_field_wght_3d_r4, sum_field_wght_3d_r8
107 end interface sum_diag_integral_field
108 
109 !> @addtogroup diag_integral_mod
110 !> @{
111 
112 private &
113 
114 ! from diag_integral_init:
115  set_axis_time, &
116 
117 ! from diag_integral_field_init and sum_diag_integral_field:
118  get_field_index, &
119 
120 ! from diag_integral_output and diag_integral_end:
122 
123 ! from write_field_averages:
125  get_axis_time, &
126 
127 ! from diag_integral_output:
129 
130 ! from sum_diag_integral_field:
132 
133 
134 !-------------------------------------------------------------------------------
135 !------ namelist -------
136 
137 real(r8_kind) :: &
138  output_interval = -1.0_r8_kind !< time interval at which integrals
139  !! are to be output
140 character(len=8) :: &
141  time_units = 'hours' !< time units associated with
142  !! output_interval
143 character(len=FMS_FILE_LEN) :: &
144  file_name = ' ' !< optional integrals output file name
145 logical :: &
146  print_header = .true. !< print a header for the integrals
147  !! file ?
148 integer :: &
149  fields_per_print_line = 4 !< number of fields to write per line
150  !! of output
151 
152 namelist / diag_integral_nml / &
156 
157 !-------------------------------------------------------------------------------
158 !------- private data ------
159 
160 !-------------------------------------------------------------------------------
161 ! variables associated with the determination of when integrals
162 ! are to be written.
163 !-------------------------------------------------------------------------------
164 type (time_type) :: next_alarm_time !< next time at which integrals are to be written
165 type (time_type) :: alarm_interval !< time interval between writing integrals
166 type (time_type) :: zero_time !< time_type variable set to (0,0); used as
167  !! flag to indicate integrals are not being output
168 type (time_type) :: time_init_save !< initial time associated with experiment;
169  !! used as a base for defining time
170 
171 !-------------------------------------------------------------------------------
172 ! variables used in determining weights associated with each
173 ! contribution to the integrand.
174 !-------------------------------------------------------------------------------
175 real(r8_kind), allocatable, dimension(:,:) :: area !< area of each grid box
176 integer :: idim !< x dimension of grid on local processor
177 integer :: jdim !< y dimension of grid on local processor
178 integer :: field_size !< number of columns on global domain
179 real(r8_kind) :: sum_area !< surface area of globe
180 
181 !-------------------------------------------------------------------------------
182 ! variables used to define the integral fields:
183 !-------------------------------------------------------------------------------
184 integer, parameter :: max_len_name = 12 !< maximum length of name associated with integral
185 integer, parameter :: max_num_field = 32 !< maximum number of integrals allowed
186 integer :: num_field = 0 !< number of integrals that have been activated
187 character(len=max_len_name) :: field_name (max_num_field) !< name associated with integral i
188 character(len=16) :: field_format (max_num_field) !< output format for integral i
189 real(r8_kind) :: field_sum (max_num_field) !< integrand for integral i
190 integer :: field_count (max_num_field) !< number of values in integrand i
191 
192 !-------------------------------------------------------------------------------
193 ! variables defining output formats.
194 !-------------------------------------------------------------------------------
195 character(len=160) :: format_text !< format statement for header
196 character(len=160) :: format_data !< format statement for data output
197 integer :: nd !< number of characters in data format statement
198 integer :: nt !< number of characters in text format statement
199 
200 !-------------------------------------------------------------------------------
201 ! miscellaneous variables.
202 !-------------------------------------------------------------------------------
203 integer :: diag_unit = 0 !< unit number for output file
204 logical :: module_is_initialized = .false. !< module is initialized ?
205 
206 
207 
208  contains
209 
210 
211 
212 !###############################################################################
213 !
214 ! PUBLIC SUBROUTINES
215 !
216 !###############################################################################
217 
218 !> @brief diag_integral_init is the constructor for diag_integral_mod.
219 !!
220 !! <b> Template: </b>
221 !!
222 !! @code{.f90}
223 !! call diag_integral_init (Time_init, Time, blon, blat)
224 !! @endcode
225 !!
226 !! <b> Parameters: </b>
227 !!
228 !! @code{.f90}
229 !! type (time_type), intent(in), optional :: Time_init, Time
230 !! real,dimension(:,:), intent(in), optional :: blon, blat, area_in
231 !! @endcode
232 !!
233 !! @param [in] <Time_init> Initial time to start the integral
234 !! @param [in] <Time> current time
235 !! @param [in] <latb> array of model latitudes at cell boundaries [radians]
236 !! @param [in] <lonb> array of model longitudes at cell boundaries [radians]
237 !!
238 subroutine diag_integral_init (Time_init, Time, blon, blat, area_in)
239 
240 type (time_type), intent(in), optional :: time_init !< Initial time to start the integral
241 type (time_type), intent(in), optional :: time !< current time
242 class(*),dimension(:,:), intent(in), optional :: blon !< array of model latitudes at cell boundaries [radians]
243 class(*),dimension(:,:), intent(in), optional :: blat !< array of model longitudes at cell boundaries [radians]
244 class(*),dimension(:,:), intent(in), optional :: area_in
245 
246 !-------------------------------------------------------------------------------
247 ! local variables:
248 ! r2
249 ! rsize
250 ! unit
251 ! io
252 ! ierr
253 ! seconds
254 ! nc
255 ! i,j
256 !-------------------------------------------------------------------------------
257  real(r8_kind) :: rsize
258  integer :: io, ierr, nc, logunit
259  integer :: field_size_local
260  real(r8_kind) :: sum_area_local
261  integer :: ensemble_size(6)
262 !-------------------------------------------------------------------------------
263 ! if routine has already been executed, exit.
264 !-------------------------------------------------------------------------------
265  if (module_is_initialized) return
266 
267 !-------------------------------------------------------------------------------
268 ! verify that modules used by this module that are not called later
269 ! have already been initialized.
270 !-------------------------------------------------------------------------------
271  call fms_init
272  call mpp_init
273  call constants_init
274  call time_manager_init
275 
276 !-------------------------------------------------------------------------------
277 ! if this is the initialization call, proceed. if this was simply
278 ! a verification of previous initialization, return.
279 !-------------------------------------------------------------------------------
280  if (present(time_init) .and. present(time) .and. &
281  present(blon) .and. present(blat) ) then
282 
283 !-------------------------------------------------------------------------------
284 ! read namelist.
285 !-------------------------------------------------------------------------------
286  read (input_nml_file, nml=diag_integral_nml, iostat=io)
287  ierr = check_nml_error(io,'diag_integral_nml')
288 
289 !-------------------------------------------------------------------------------
290 ! write version number and namelist to logfile.
291 !-------------------------------------------------------------------------------
292  call write_version_number ('DIAG_INTEGRAL_MOD', version)
293  logunit = stdlog()
294  if (mpp_pe() == mpp_root_pe() ) &
295  write (logunit, nml=diag_integral_nml)
296 
297 !-------------------------------------------------------------------------------
298 ! save the initial time to time-stamp the integrals which will be
299 ! calculated.
300 !-------------------------------------------------------------------------------
301  time_init_save = time_init
302 
303 !-------------------------------------------------------------------------------
304 ! define the model grid sizes and the total number of columns on
305 ! the processor. sum over all processors and store the global
306 ! number of columns in field_size.
307 !-------------------------------------------------------------------------------
308  idim = size(blon,1) - 1
309  jdim = size(blon,2) - 1
310  field_size_local = idim*jdim
311  rsize = real(field_size_local,r8_kind)
312  call mpp_sum (rsize)
313  field_size = nint(rsize)
314 
315 !-------------------------------------------------------------------------------
316 ! define an array to hold the surface area of each grid column
317 ! so that the integrals may be weighted properly. sum over the
318 ! processor, and then over all processors, storing the total
319 ! global surface area in sum_area.
320 !-------------------------------------------------------------------------------
321  allocate (area(idim,jdim))
322 
323  select type(area_in)
324  type is (real(r4_kind)) ; area = real(area_in,r8_kind)
325  type is (real(r8_kind)) ; area = area_in
326  class default ; call error_mesg('diag_inetgral_mod::diag_integral_init', 'unknown area_in type', fatal)
327  end select
328 
329  sum_area_local = sum(area)
330  sum_area = sum_area_local
331  call mpp_sum (sum_area)
332 
333 !-------------------------------------------------------------------------------
334 ! if integral output is to go to a file, open the file on unit
335 ! diag_unit.
336 !-------------------------------------------------------------------------------
337  if (file_name(1:1) /= ' ' ) then
338  ensemble_size = get_ensemble_size()
339  if (ensemble_size(1) > 1) then
341  endif
342  nc = len_trim(file_name)
343  open(newunit=diag_unit, file=file_name(1:nc), action='write')
344  endif
345 
346 !-------------------------------------------------------------------------------
347 ! define the variables needed to control the time interval of
348 ! output. Zero time is a flag indicating that the alarm is not set,
349 ! i.e., integrals are not desired. otherwise set the next time to
350 ! output integrals to be at the value of nml variable
351 ! output_interval from now.
352 !-------------------------------------------------------------------------------
353  zero_time = set_time(0,0)
354  if (output_interval >= -0.01_r8_kind) then
357  else
359  endif
361 
362 !-------------------------------------------------------------------------------
363 ! deallocate the local array and mark the module as initialized.
364 !-------------------------------------------------------------------------------
365  module_is_initialized = .true.
366  endif ! (present optional arguments)
367 
368 end subroutine diag_integral_init
369 
370 
371 
372 
373 !###############################################################################
374 !> @brief diag_integral_field_init registers and initializes an integral field
375 !!
376 !! <b> Template: </b>
377 !!
378 !! @code{.f90}
379 !! call diag_integral_field_init (name, format)
380 !! @endcode
381 !!
382 !! <b> Parameters: </b>
383 !!
384 !! @code{.f90}
385 !! character(len=*), intent(in) :: name, format
386 !! @endcode
387 !!
388 !! @param [in] <name> Name of the field to be integrated
389 !! @param [in] <format> Output format of the field to be integrated
390 !!
391 subroutine diag_integral_field_init (name, format)
392 
393 character(len=*), intent(in) :: name !< Name of the field to be integrated
394 character(len=*), intent(in) :: format !< Output format of the field to be integrated
395 
396 !-------------------------------------------------------------------------------
397 ! local variables:
398 !-------------------------------------------------------------------------------
399  integer :: field !< index assigned to the current integral
400 
401 !-------------------------------------------------------------------------------
402 ! note: no initialization is required for this interface. all needed
403 ! variables are initialized in the source.
404 !-------------------------------------------------------------------------------
405 
406 !-------------------------------------------------------------------------------
407 ! make sure the integral name is not too long.
408 !-------------------------------------------------------------------------------
409  if (len(name) > max_len_name ) then
410  call error_mesg ('diag_integral_mod', &
411  ' integral name too long', fatal)
412  endif
413 
414 !-------------------------------------------------------------------------------
415 ! check to be sure the integral name has not already been
416 ! initialized.
417 !-------------------------------------------------------------------------------
418  field = get_field_index(name)
419  if (field /= 0) then
420  call error_mesg ('diag_integral_mod', &
421  'integral name already exists', fatal)
422  endif
423 
424 !-------------------------------------------------------------------------------
425 ! prepare to register the integral. make sure that there are not
426 ! more integrals registered than space was provided for; if so, exit.
427 !-------------------------------------------------------------------------------
428  num_field = num_field + 1
429  if (num_field > max_num_field) then
430  call error_mesg ('diag_integral_mod', &
431  'too many fields initialized', fatal)
432  endif
433 
434 !-------------------------------------------------------------------------------
435 ! register the name and output format desired for the given integral.
436 ! initialize its value and the number of grid points that have been
437 ! counted to zero.
438 !-------------------------------------------------------------------------------
439  field_name(num_field) = name
440  field_format(num_field) = format
441  field_sum(num_field) = 0.0_r8_kind
443 
444 end subroutine diag_integral_field_init
445 
446 
447 
448 !###############################################################################
449 !> @brief diag_integral_output determines if this is a timestep on which
450 !! integrals are to be written. if not, it returns; if so, it calls
451 !! write_field_averages.
452 !!
453 !! <b> Template: </b>
454 !!
455 !! @code{.f90}
456 !! call diag_integral_output (Time)
457 !! @endcode
458 !!
459 !! <b> Parameters: </b>
460 !!
461 !! @code{.f90}
462 !! type (time_type), intent(in) :: Time
463 !! @endcode
464 !!
465 !! @param [in] <Time> integral time stamp at the current time
466 !!
467 subroutine diag_integral_output (Time)
468 
469 type (time_type), intent(in) :: time
470 
471 !-------------------------------------------------------------------------------
472 ! be sure module has been initialized.
473 !-------------------------------------------------------------------------------
474  if (.not. module_is_initialized ) then
475  call error_mesg ('diag_integral_mod', &
476  'module has not been initialized', fatal )
477  endif
478 
479 !-------------------------------------------------------------------------------
480 ! see if integral output is desired at this time.
481 !-------------------------------------------------------------------------------
482  if ( diag_integral_alarm(time) ) then
483 !-------------------------------------------------------------------------------
484 ! write the integrals by calling write_field_averages. upon return
485 ! reset the alarm to the next diagnostics time.
486 !-------------------------------------------------------------------------------
487  call write_field_averages (time)
489  endif
490 
491 end subroutine diag_integral_output
492 
493 
494 
495 !###############################################################################
496 !> @brief diag_integral_end is the destructor for diag_integral_mod.
497 !!
498 !! <b> Template: </b>
499 !!
500 !! @code{.f90}
501 !! call diag_integral_end (Time)
502 !! @endcode
503 !!
504 !! <b> Parameters: </b>
505 !!
506 !! @code{.f90}
507 !! type (time_type), intent(in) :: Time
508 !! @endcode
509 !!
510 !! @param [in] <Time> integral time stamp at the current time
511 !!
512 subroutine diag_integral_end (Time)
513 
514 type (time_type), intent(in) :: time
515 
516 !-------------------------------------------------------------------------------
517 ! be sure module has been initialized.
518 !-------------------------------------------------------------------------------
519  if (.not. module_is_initialized ) then
520  call error_mesg ('diag_integral_mod', &
521  'module has not been initialized', fatal )
522  endif
523 
524 !-------------------------------------------------------------------------------
525 ! if the alarm interval was set to Zero_time (meaning no integral
526 ! output during the model run) call write_field_averages to output
527 ! the integrals valid over the entire period of integration.
528 !-------------------------------------------------------------------------------
529  if (alarm_interval == zero_time ) then
530 ! if (Alarm_interval /= Zero_time ) then
531 ! else
532  call write_field_averages (time)
533  endif
534 
535 !-------------------------------------------------------------------------------
536 ! deallocate module variables.
537 !-------------------------------------------------------------------------------
538  deallocate (area)
539  if (diag_unit /= 0) close(diag_unit)
540 
541 !-------------------------------------------------------------------------------
542 ! mark the module as uninitialized.
543 !-------------------------------------------------------------------------------
544  module_is_initialized = .false.
545 
546 end subroutine diag_integral_end
547 
548 
549 
550 !###############################################################################
551 !
552 ! PRIVATE SUBROUTINES
553 !
554 !###############################################################################
555 
556 
557 
558 !###############################################################################
559 !> @brief Function to convert input time to a time_type
560 !!
561 !! <b> Template: </b>
562 !!
563 !! @code{.f90}
564 !! time = set_axis_time (atime, units)
565 !! @endcode
566 !!
567 !! <b> Parameters: </b>
568 !!
569 !! @code{.f90}
570 !! real, intent(in) :: atime
571 !! character(len=*), intent(in) :: units
572 !! type(time_type) :: Time
573 !! @endcode
574 !!
575 !! @param [in] <atime> integral time stamp at the current time
576 !! @param [in] <units> input units, not used
577 !! @param [out] <Time>
578 !!
579 function set_axis_time (atime, units) result (Time)
580 
581 real(r8_kind), intent(in) :: atime !< integral time stamp at the current time
582 character(len=*), intent(in) :: units !< input units, not used
583 type(time_type) :: time
584 
585 !-------------------------------------------------------------------------------
586 ! local variables:
587 !-------------------------------------------------------------------------------
588  integer :: sec !< seconds corresponding to the input
589  !! variable atime
590  integer :: day = 0 !< day component of time_type variable
591 
592 !-------------------------------------------------------------------------------
593 ! convert the input time to seconds, regardless of input units.
594 !-------------------------------------------------------------------------------
595  if (units(1:3) == 'sec') then
596  sec = int(atime + 0.5)
597  else if (units(1:3) == 'min') then
598  sec = int(atime*60. + 0.5)
599  else if (units(1:3) == 'hou') then
600  sec = int(atime*3600. + 0.5)
601  else if (units(1:3) == 'day') then
602  sec = int(atime*86400. + 0.5)
603  else
604  call error_mesg('diag_integral_mod', &
605  'Invalid units sent to set_axis_time', fatal)
606  endif
607 
608 !-------------------------------------------------------------------------------
609 ! convert the time in seconds to a time_type variable.
610 !-------------------------------------------------------------------------------
611  time = set_time(sec, day)
612 
613 end function set_axis_time
614 
615 
616 
617 !###############################################################################
618 !> @brief get_field_index returns returns the index associated with an
619 !! integral name.
620 !!
621 !! <b> Template: </b>
622 !!
623 !! @code{.f90}
624 !! index = get_field_index (name)
625 !! @endcode
626 !!
627 !! <b> Parameters: </b>
628 !!
629 !! @code{.f90}
630 !! character(len=*), intent(in) :: name
631 !! integer :: index
632 !! @endcode
633 !!
634 !! @param [in] <name> Name associated with an integral
635 !! @param [out] <index>
636 !!
637 function get_field_index (name) result (index)
638 
639 character(len=*), intent(in) :: name !< Name associated with an integral
640 integer :: index
641 
642 !-------------------------------------------------------------------------------
643 ! local variables:
644 !-------------------------------------------------------------------------------
645  integer :: nc
646  integer :: i
647 
648  nc = len_trim(name)
649  if (nc > max_len_name) then
650  call error_mesg ('diag_integral_mod', &
651  'name too long', fatal)
652  endif
653 
654 !-------------------------------------------------------------------------------
655 ! search each field name for the current string. when found exit
656 ! with the index. if not found index will be 0 upon return, which
657 ! initiates error condition.
658 !-------------------------------------------------------------------------------
659  index = 0
660  do i = 1, num_field
661  if (name(1:nc) == &
662  field_name(i) (1:len_trim(field_name(i))) ) then
663  index = i
664  exit
665  endif
666  end do
667 
668 end function get_field_index
669 
670 
671 
672 !###############################################################################
673 !> @brief Subroutine to sum multiple fields, average them and then write the
674 !! result to an output file.
675 !!
676 !! <b> Template: </b>
677 !!
678 !! @code{.f90}
679 !! call write_field_averages (Time)
680 !! @endcode
681 !!
682 !! <b> Parameters: </b>
683 !!
684 !! @code{.f90}
685 !! type (time_type), intent(in) :: Time
686 !! @endcode
687 !!
688 !! @param [in] <Time> integral time stamp at the current time
689 !!
690 subroutine write_field_averages (Time)
691 
692 type (time_type), intent(in) :: time !< integral time stamp at the current time
693 
694 !-------------------------------------------------------------------------------
695 ! local variables:
696 ! field_avg
697 ! xtime
698 ! rcount
699 ! nn
700 ! ninc
701 ! nst
702 ! nend
703 ! fields_to_print
704 ! i
705 ! kount
706 !-------------------------------------------------------------------------------
707  real(r8_kind) :: field_avg(max_num_field)
708  real(r8_kind) :: xtime, rcount
709  integer :: nn, ninc, nst, nend, fields_to_print
710  integer :: i, kount
711  integer(i8_kind) :: icount
712  character(len=128) :: xtime_str
713  logical :: use_exp_format
714 
715 !-------------------------------------------------------------------------------
716 ! each header and data format may be different and must be generated
717 ! as needed.
718 !-------------------------------------------------------------------------------
719  fields_to_print = 0
720  do i = 1, num_field
721 
722 !-------------------------------------------------------------------------------
723 ! increment the fields_to_print counter. sum the integrand and the
724 ! number of data points contributing to it over all processors.
725 !-------------------------------------------------------------------------------
726  fields_to_print = fields_to_print + 1
727  rcount = real(field_count(i),r8_kind)
728  call mpp_sum (rcount)
729  call mpp_sum (field_sum(i))
730  icount = int(rcount, i8_kind)
731 
732 !-------------------------------------------------------------------------------
733 ! verify that all the data expected for an integral has been
734 ! obtained.
735 !-------------------------------------------------------------------------------
736  if (icount == 0 ) call error_mesg &
737  ('diag_integral_mod', &
738  'field_count equals zero for field_name ' // &
739  field_name(i)(1:len_trim(field_name(i))), fatal )
740  kount = int(icount/field_size)
741  if ((field_size)*kount /= icount) then
742  print*,"name,pe,kount,field_size,icount,rcount=",trim(field_name(i)),mpp_pe(),kount,field_size,icount,rcount
743  call error_mesg &
744  ('diag_integral_mod', &
745  'field_count not a multiple of field_size', fatal )
746  endif
747 !-------------------------------------------------------------------------------
748 ! define the global integral for field i. reinitialize the point
749 ! and data accumulators.
750 !-------------------------------------------------------------------------------
751  field_avg(fields_to_print) = field_sum(i)/ &
752  (sum_area*real(kount,r8_kind))
753  field_sum(i) = 0.0_r8_kind
754  field_count(i) = 0
755  end do
756 
757 !-------------------------------------------------------------------------------
758 ! only the root pe will write out data.
759 !-------------------------------------------------------------------------------
760  if ( mpp_pe() /= mpp_root_pe() ) return
761 
762 !-------------------------------------------------------------------------------
763 ! define the time associated with the integrals just calculated.
764 !-------------------------------------------------------------------------------
766 
767 !-------------------------------------------------------------------------------
768 ! check if the time value is too long for decimal output
769 !-------------------------------------------------------------------------------
770  xtime_str = trim(string(xtime))
771  use_exp_format = len_trim(xtime_str(1:index(xtime_str, "."))) .ge. 9
772 
773 !-------------------------------------------------------------------------------
774 ! generate the new header and data formats.
775 !-------------------------------------------------------------------------------
776  nst = 1
777  nend = fields_per_print_line
778  ninc = (num_field-1)/fields_per_print_line + 1
779  do nn=1, ninc
780  nst = 1 + (nn-1)*fields_per_print_line
781  nend = min(nn*fields_per_print_line, num_field)
782  if (print_header) call format_text_init (nst, nend)
783  call format_data_init (nst, nend, use_exp_format)
784  if (diag_unit /= 0) then
785  write (diag_unit,format_data(1:nd)) &
786  xtime, (field_avg(i),i=nst,nend)
787  else
788  write (*, format_data(1:nd)) &
789  xtime, (field_avg(i),i=nst,nend)
790  endif
791  end do
792 
793 end subroutine write_field_averages
794 
795 
796 
797 !###############################################################################
798 !> @brief format_text_init generates the header records to be output in the
799 !! integrals file.
800 !!
801 !! <b> Template: </b>
802 !!
803 !! @code{.f90}
804 !! call format_text_init (nst_in, nend_in)
805 !! @endcode
806 !!
807 !! <b> Parameters: </b>
808 !!
809 !! @code{.f90}
810 !! integer, intent(in), optional :: nst_in, nend_in
811 !! @endcode
812 !!
813 !! @param [in] <nst_in, nend_in> starting/ending integral index which will be
814 !! included in this format statement
815 !!
816 subroutine format_text_init (nst_in, nend_in)
817 
818 integer, intent(in), optional :: nst_in !< starting/ending integral index which will be
819  !! included in this format statement
820 integer, intent(in), optional :: nend_in !< starting/ending integral index which will be
821  !! included in this format statement
822 
823 !-------------------------------------------------------------------------------
824 ! local variables:
825 ! i
826 ! nc
827 ! nst
828 ! nend
829 !-------------------------------------------------------------------------------
830  integer :: i, nc, nst, nend
831 
832 !-------------------------------------------------------------------------------
833 ! only the root pe need execute this routine, since only it will
834 ! be outputting integrals.
835 !-------------------------------------------------------------------------------
836  if (mpp_pe() /= mpp_root_pe()) return
837 
838 !-------------------------------------------------------------------------------
839 ! define the starting and ending integral indices that will be
840 ! included in this format statement.
841 !-------------------------------------------------------------------------------
842  if (present (nst_in) ) then
843  nst = nst_in
844  nend = nend_in
845  else
846  nst = 1
847  nend = num_field
848  endif
849 
850 !-------------------------------------------------------------------------------
851 ! define the first 11 characters in the format statement.
852 !-------------------------------------------------------------------------------
853  nt = 11
854  format_text(1:nt) = "('# time"
855 
856 !-------------------------------------------------------------------------------
857 ! generate the rest of the format statement, which will cover
858 ! integral indices nst to nend. if satndard printout is desired,
859 ! cycle through the loop.
860 !-------------------------------------------------------------------------------
861  do i=nst,nend
862  nc = len_trim(field_name(i))
863  format_text(nt+1:nt+nc+5) = ' ' // field_name(i)(1:nc)
864  nt = nt+nc+5
865  end do
866 
867 !-------------------------------------------------------------------------------
868 ! include the end of the format statement.
869 !-------------------------------------------------------------------------------
870  format_text(nt+1:nt+2) = "')"
871  nt = nt+2
872 
873 !-------------------------------------------------------------------------------
874 ! write the format statement to either an output file or to stdout.
875 !-------------------------------------------------------------------------------
876  if (diag_unit /= 0) then
877  write (diag_unit, format_text(1:nt))
878  else
879  write (*, format_text(1:nt))
880  endif
881 
882 end subroutine format_text_init
883 
884 
885 
886 !###############################################################################
887 !> @brief format_text_init generates the format to be output in the
888 !! integrals file.
889 !!
890 !! <b> Template: </b>
891 !!
892 !! @code{.f90}
893 !! call format_data_init (nst_in, nend_in)
894 !! @endcode
895 !!
896 !! <b> Parameters: </b>
897 !!
898 !! @code{.f90}
899 !! integer, intent(in) :: nst_in, nend_in
900 !! @endcode
901 !!
902 !! @param [in] <nst_in, nend_in> starting/ending integral index which will be
903 !! included in this format statement
904 !! @param [in] <use_exp_format> if true, uses exponent notation for the first format code
905 !! to avoid overflow with larger time values
906 !!
907 subroutine format_data_init (nst_in, nend_in, use_exp_format)
908 
909 integer, intent(in) :: nst_in !< starting/ending integral index which will be
910  !! included in this format statement
911 integer, intent(in) :: nend_in !< starting/ending integral index which will be
912  !! included in this format statement
913 logical, intent(in) :: use_exp_format !< if true, uses exponent notation for the first format code
914  !! to avoid overflow with larger time values
915 
916 !-------------------------------------------------------------------------------
917 ! local variables:
918 ! i
919 ! nc
920 ! nst
921 ! nend
922 !-------------------------------------------------------------------------------
923  integer :: i, nc, nst, nend
924 
925 !-------------------------------------------------------------------------------
926 ! define the start of the format, which covers the time stamp of the
927 ! integrals. this section is 9 characters long.
928 !-------------------------------------------------------------------------------
929  nd = 9
930  if( use_exp_format ) then
931  format_data(1:nd) = '(1x,e10.2'
932  else
933  format_data(1:nd) = '(1x,f10.2'
934  endif
935 
936 !-------------------------------------------------------------------------------
937 ! define the indices of the integrals that are to be written by this
938 ! format statement.
939 !-------------------------------------------------------------------------------
940  nst = nst_in
941  nend = nend_in
942 
943 !-------------------------------------------------------------------------------
944 ! complete the data format. use the format defined for the
945 ! particular integral in setting up the format statement.
946 !-------------------------------------------------------------------------------
947  do i=nst,nend
948  nc = len_trim(field_format(i))
949  format_data(nd+1:nd+nc+4) = ',1x,' // field_format(i)(1:nc)
950  nd = nd+nc+4
951  end do
952 
953 !-------------------------------------------------------------------------------
954 ! close the format statement.
955 !-------------------------------------------------------------------------------
956  format_data(nd+1:nd+1) = ')'
957  nd = nd + 1
958 
959 end subroutine format_data_init
960 
961 
962 !###############################################################################
963 !> @brief Function to convert the time_type input variable into units of
964 !! units and returns it in atime.
965 !!
966 !! <b> Template: </b>
967 !!
968 !! @code{.f90}
969 !! atime = get_axis_time (Time, units)
970 !! @endcode
971 !!
972 !! <b> Parameters: </b>
973 !!
974 !! @code{.f90}
975 !! type(time_type), intent(in) :: Time
976 !! character(len=*), intent(in) :: units
977 !! real :: atime
978 !! @endcode
979 !!
980 !! @param [in] <Time> integral time stamp
981 !! @param [in] <units> input units of time_type
982 !! @param [out] <atime>
983 !!
984 !! @return real atime
985 function get_axis_time (Time, units) result (atime)
986 
987 type(time_type), intent(in) :: time !< integral time stamp
988 character(len=*), intent(in) :: units !< input units of time_type
989 real(r8_kind) :: atime
990 
991 !-------------------------------------------------------------------------------
992 ! local variables:
993 !-------------------------------------------------------------------------------
994  integer :: sec, day ! components of time_type variable
995 
996  call get_time (time, sec, day)
997  if (units(1:3) == 'sec') then
998  atime = real(sec,r8_kind) + 86400._r8_kind*real(day,r8_kind)
999  else if (units(1:3) == 'min') then
1000  atime = real(sec,r8_kind)/60._r8_kind + 1440._r8_kind*real(day,r8_kind)
1001  else if (units(1:3) == 'hou') then
1002  atime = real(sec,r8_kind)/3600._r8_kind + 24._r8_kind*real(day,r8_kind)
1003  else if (units(1:3) == 'day') then
1004  atime = real(sec,r8_kind)/86400._r8_kind + real(day,r8_kind)
1005  endif
1006 
1007 end function get_axis_time
1008 
1009 
1010 
1011 !###############################################################################
1012 !> @brief Function to check if it is time to write integrals.
1013 !! if not writing integrals, return.
1014 !!
1015 !! <b> Template: </b>
1016 !!
1017 !! @code{.f90}
1018 !! result = diag_integral_alarm (Time)
1019 !! @endcode
1020 !!
1021 !! <b> Parameters: </b>
1022 !!
1023 !! @code{.f90}
1024 !! type (time_type), intent(in) :: Time
1025 !! logical :: answer
1026 !! @endcode
1027 !!
1028 !! @param [in] <Time> current time
1029 !! @param [out] <answer>
1030 !!
1031 !! @return logical answer
1032 function diag_integral_alarm (Time) result (answer)
1033 
1034 type (time_type), intent(in) :: time !< current time
1035 logical :: answer
1036 
1037  answer = .false.
1038  if (alarm_interval == zero_time) return
1039  if (time >= next_alarm_time) answer = .true.
1040 
1041 end function diag_integral_alarm
1042 
1043 
1044 
1045 !###############################################################################
1046 !> @brief Function to perform a weighted integral in the vertical
1047 !! direction of a 3d data field
1048 !!
1049 !! <b> Template: </b>
1050 !!
1051 !! @code{.f90}
1052 !! data2 = vert_diag_integral (field_data, wt)
1053 !! @endcode
1054 !!
1055 !! <b> Parameters: </b>
1056 !!
1057 !! @code{.f90}
1058 !! real, dimension (:,:,:), intent(in) :: field_data, wt
1059 !! real, dimension (size(field_data,1),size(field_data,2)) :: data2
1060 !! @endcode
1061 !!
1062 !! @param [in] <field_data> integral field data arrays
1063 !! @param [in] <wt> integral field weighting functions
1064 !! @param [out] <data2>
1065 !! @return real array data2
1066 function vert_diag_integral (field_data, wt) result (data2)
1067 real(r8_kind), dimension (:,:,:), intent(in) :: field_data !< integral field data arrays
1068 real(r8_kind), dimension (:,:,:), intent(in) :: wt !< integral field weighting functions
1069 real(r8_kind), dimension (size(field_data,1),size(field_data,2)) :: data2
1070 
1071 !-------------------------------------------------------------------------------
1072 ! local variables:
1073 ! wt2
1074 !-------------------------------------------------------------------------------
1075  real, dimension(size(field_data,1),size(field_data,2)) :: wt2
1076 
1077 !-------------------------------------------------------------------------------
1078  wt2 = sum(wt,3)
1079  if (count(wt2 == 0._r8_kind) > 0) then
1080  call error_mesg ('diag_integral_mod', &
1081  'vert sum of weights equals zero', fatal)
1082  endif
1083  data2 = sum(field_data*wt,3) / wt2
1084 
1085 end function vert_diag_integral
1086 
1087 !> @brief Adds .ens_## to the diag_integral.out file name
1088 !! @return character array updated_file_name
1089 function ensemble_file_name(fname) result(updated_file_name)
1090  character (len=*), intent(inout) :: fname
1091  character (len=FMS_FILE_LEN) :: updated_file_name
1092  integer :: ensemble_id_int
1093  character(len=7) :: ensemble_suffix
1094  character(len=2) :: ensemble_id_char
1095  integer :: i
1096  !> Make sure the file name short enough to handle adding the ensemble number
1097  if (len(trim(fname)) > fms_file_len-7) call error_mesg ('diag_integral_mod :: ensemble_file_name', &
1098  trim(fname)//" is too long and can not support adding ens_XX. Please shorten the "//&
1099  "file_name in the diag_integral_nml", fatal)
1100  !> Get the ensemble ID and convert it to a string
1101  ensemble_id_int = get_ensemble_id()
1102  write(ensemble_id_char,"(I0)") ensemble_id_int
1103  !> Add a 0 if the ensemble is less than 10 (2 digit)
1104  if (ensemble_id_int < 10) then
1105  ensemble_suffix = ".ens_0"//trim(ensemble_id_char)
1106  elseif (ensemble_id_int >= 10 .and. ensemble_id_int < 100) then
1107  ensemble_suffix = ".ens_"//trim(ensemble_id_char)
1108  else
1109  call error_mesg ('diag_integral_mod', &
1110  ' Does not support ensemble sizes over 99.', fatal)
1111  endif
1112  !> Insert the ens_ string in the correct location
1113  !> Loop through to find the last period
1114  do i=len(trim(fname)),2,-1
1115  if (fname(i:i) == ".") then
1116  updated_file_name = fname(1:i-1)//trim(ensemble_suffix)//fname(i:len(fname))
1117  return
1118  endif
1119  enddo
1120  !> Add to the end if there is no period in the file name
1121  updated_file_name = trim(fname)//trim(ensemble_suffix)
1122 end function ensemble_file_name
1123 
1124 #include "diag_integral_r4.fh"
1125 #include "diag_integral_r8.fh"
1126 
1127  end module diag_integral_mod
1128 !> @}
1129 ! close documentation grouping
type(time_type) function, private set_axis_time(atime, units)
Function to convert input time to a time_type.
real(r8_kind) sum_area
surface area of globe
character(len=fms_file_len) file_name
optional integrals output file name
real(r8_kind) function, dimension(size(field_data, 1), size(field_data, 2)), private vert_diag_integral(field_data, wt)
Function to perform a weighted integral in the vertical direction of a 3d data field.
type(time_type) zero_time
time_type variable set to (0,0); used as flag to indicate integrals are not being output
subroutine, public diag_integral_output(Time)
diag_integral_output determines if this is a timestep on which integrals are to be written....
subroutine, public diag_integral_field_init(name, format)
diag_integral_field_init registers and initializes an integral field
integer num_field
number of integrals that have been activated
integer diag_unit
unit number for output file
real(r8_kind) function, private get_axis_time(Time, units)
Function to convert the time_type input variable into units of units and returns it in atime.
integer idim
x dimension of grid on local processor
integer fields_per_print_line
number of fields to write per line
integer nt
number of characters in text format statement
integer jdim
y dimension of grid on local processor
subroutine, private format_data_init(nst_in, nend_in, use_exp_format)
format_text_init generates the format to be output in the integrals file.
character(len=16), dimension(max_num_field) field_format
output format for integral i
real(r8_kind), dimension(:,:), allocatable area
area of each grid box
real(r8_kind) output_interval
time interval at which integrals
integer, parameter max_num_field
maximum number of integrals allowed
character(len=fms_file_len) function ensemble_file_name(fname)
Adds .ens_## to the diag_integral.out file name.
subroutine, private format_text_init(nst_in, nend_in)
format_text_init generates the header records to be output in the integrals file.
type(time_type) time_init_save
initial time associated with experiment; used as a base for defining time
subroutine, public diag_integral_init(Time_init, Time, blon, blat, area_in)
diag_integral_init is the constructor for diag_integral_mod.
integer function, private get_field_index(name)
get_field_index returns returns the index associated with an integral name.
type(time_type) alarm_interval
time interval between writing integrals
integer, parameter max_len_name
maximum length of name associated with integral
character(len=max_len_name), dimension(max_num_field) field_name
name associated with integral i
subroutine, private write_field_averages(Time)
Subroutine to sum multiple fields, average them and then write the result to an output file.
character(len=160) format_data
format statement for data output
character(len=8) time_units
time units associated with
type(time_type) next_alarm_time
next time at which integrals are to be written
logical module_is_initialized
module is initialized ?
integer field_size
number of columns on global domain
character(len=160) format_text
format statement for header
subroutine, public diag_integral_end(Time)
diag_integral_end is the destructor for diag_integral_mod.
real(r8_kind), dimension(max_num_field) field_sum
integrand for integral i
integer, dimension(max_num_field) field_count
number of values in integrand i
integer nd
number of characters in data format statement
logical print_header
print a header for the integrals
logical function, private diag_integral_alarm(Time)
Function to check if it is time to write integrals. if not writing integrals, return.
Perform a summation of the named integral field.
integer function, public get_ensemble_id()
Getter function for ensemble_id.
integer function, dimension(6), public get_ensemble_size()
Returns ensemble size integer array.
integer function, public check_nml_error(IOSTAT, NML_NAME)
Checks the iostat argument that is returned after reading a namelist and determines if the error code...
Definition: fms.F90:523
subroutine, public write_version_number(version, tag, unit)
Prints to the log file (or a specified unit) the version id string and tag name.
Definition: fms.F90:701
subroutine, public fms_init(localcomm, alt_input_nml_path)
Initializes the FMS module and also calls the initialization routines for all modules in the MPP pack...
Definition: fms.F90:286
subroutine, public error_mesg(routine, message, level)
Print notes, warnings and error messages; terminates program for warning and error messages....
Definition: fms.F90:441
subroutine mpp_init(flags, localcomm, test_level, alt_input_nml_path)
Initialize the mpp_mod module. Must be called before any usage.
integer function stdlog()
This function returns the current standard fortran unit numbers for log messages. Log messages,...
Definition: mpp_util.inc:58
integer function mpp_pe()
Returns processor ID.
Definition: mpp_util.inc:406
Reduction operation.
Definition: mpp.F90:596
subroutine, public get_time(Time, seconds, days, ticks, err_msg)
Returns days and seconds ( < 86400 ) corresponding to a time. err_msg should be checked for any error...
subroutine, public time_manager_init()
Initialization routine. Writes the version information to the log file.
Given some number of seconds and days, returns the corresponding time_type.
Type to represent amounts of time. Implemented as seconds and days to allow for larger intervals.