26 module diag_integral_mod
30 use platform_mod,
only: i8_kind, fms_file_len
33 operator(+),
operator(-), &
34 operator(==),
operator(>=), &
36 use mpp_mod,
only: input_nml_file
43 use fms2_io_mod,
only: file_exists
44 use constants_mod,
only: radius, constants_init
47 use platform_mod,
only: r4_kind, r8_kind
57 #include<file_version.h>
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
140 character(len=8) :: &
143 character(len=FMS_FILE_LEN) :: &
152 namelist / diag_integral_nml / &
175 real(r8_kind),
allocatable,
dimension(:,:) ::
area
240 type (
time_type),
intent(in),
optional :: time_init
241 type (
time_type),
intent(in),
optional :: time
242 class(*),
dimension(:,:),
intent(in),
optional :: blon
243 class(*),
dimension(:,:),
intent(in),
optional :: blat
244 class(*),
dimension(:,:),
intent(in),
optional :: area_in
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)
280 if (
present(time_init) .and.
present(time) .and. &
281 present(blon) .and.
present(blat) )
then
286 read (input_nml_file, nml=diag_integral_nml, iostat=io)
294 if (
mpp_pe() == mpp_root_pe() ) &
295 write (logunit, nml=diag_integral_nml)
308 idim =
size(blon,1) - 1
309 jdim =
size(blon,2) - 1
311 rsize = real(field_size_local,r8_kind)
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)
329 sum_area_local = sum(
area)
339 if (ensemble_size(1) > 1)
then
393 character(len=*),
intent(in) :: name
394 character(len=*),
intent(in) :: format
411 ' integral name too long', fatal)
421 'integral name already exists', fatal)
431 'too many fields initialized', fatal)
476 'module has not been initialized', fatal )
521 'module has not been initialized', fatal )
581 real(r8_kind),
intent(in) :: atime
582 character(len=*),
intent(in) :: units
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)
605 'Invalid units sent to set_axis_time', fatal)
639 character(len=*),
intent(in) :: name
651 'name too long', fatal)
708 real(r8_kind) :: xtime, rcount
709 integer :: nn, ninc, nst, nend, fields_to_print
711 integer(i8_kind) :: icount
712 character(len=128) :: xtime_str
713 logical :: use_exp_format
726 fields_to_print = fields_to_print + 1
730 icount = int(rcount, i8_kind)
737 (
'diag_integral_mod', &
738 'field_count equals zero for field_name ' // &
742 print*,
"name,pe,kount,field_size,icount,rcount=",trim(
field_name(i)),
mpp_pe(),kount,
field_size,icount,rcount
744 (
'diag_integral_mod', &
745 'field_count not a multiple of field_size', fatal )
751 field_avg(fields_to_print) =
field_sum(i)/ &
760 if (
mpp_pe() /= mpp_root_pe() )
return
770 xtime_str = trim(string(xtime))
771 use_exp_format = len_trim(xtime_str(1:index(xtime_str,
"."))) .ge. 9
786 xtime, (field_avg(i),i=nst,nend)
789 xtime, (field_avg(i),i=nst,nend)
818 integer,
intent(in),
optional :: nst_in
820 integer,
intent(in),
optional :: nend_in
830 integer :: i, nc, nst, nend
836 if (
mpp_pe() /= mpp_root_pe())
return
842 if (
present (nst_in) )
then
909 integer,
intent(in) :: nst_in
911 integer,
intent(in) :: nend_in
913 logical,
intent(in) :: use_exp_format
923 integer :: i, nc, nst, nend
930 if( use_exp_format )
then
988 character(len=*),
intent(in) :: units
989 real(r8_kind) :: atime
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)
1067 real(r8_kind),
dimension (:,:,:),
intent(in) :: field_data
1068 real(r8_kind),
dimension (:,:,:),
intent(in) :: wt
1069 real(r8_kind),
dimension (size(field_data,1),size(field_data,2)) :: data2
1075 real,
dimension(size(field_data,1),size(field_data,2)) :: wt2
1079 if (count(wt2 == 0._r8_kind) > 0)
then
1081 'vert sum of weights equals zero', fatal)
1083 data2 = sum(field_data*wt,3) / wt2
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
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)
1102 write(ensemble_id_char,
"(I0)") ensemble_id_int
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)
1110 ' Does not support ensemble sizes over 99.', fatal)
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))
1121 updated_file_name = trim(fname)//trim(ensemble_suffix)
1124 #include "diag_integral_r4.fh"
1125 #include "diag_integral_r8.fh"
1127 end module diag_integral_mod
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...
subroutine, public write_version_number(version, tag, unit)
Prints to the log file (or a specified unit) the version id string and tag name.
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...
subroutine, public error_mesg(routine, message, level)
Print notes, warnings and error messages; terminates program for warning and error messages....
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,...
integer function mpp_pe()
Returns processor ID.
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.