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