FMS  2024.03
Flexible Modeling System
diag_integral_mod

This module computes and outputs global and / or hemispheric physics integrals. More...

Data Types

interface  sum_diag_integral_field
 Perform a summation of the named integral field. More...
 

Functions/Subroutines

logical function, private diag_integral_alarm (Time)
 Function to check if it is time to write integrals. if not writing integrals, return. More...
 
subroutine, public diag_integral_end (Time)
 diag_integral_end is the destructor for diag_integral_mod. More...
 
subroutine, public diag_integral_field_init (name, format)
 diag_integral_field_init registers and intializes an integral field More...
 
subroutine, public diag_integral_init (Time_init, Time, blon, blat, area_in)
 diag_integral_init is the constructor for diag_integral_mod. More...
 
subroutine, public diag_integral_output (Time)
 diag_integral_output determines if this is a timestep on which integrals are to be written. if not, it returns; if so, it calls write_field_averages. More...
 
character(len=fms_file_len) function ensemble_file_name (fname)
 Adds .ens_## to the diag_integral.out file name. More...
 
subroutine, private format_data_init (nst_in, nend_in)
 format_text_init generates the format to be output in the integrals file. More...
 
subroutine, private format_text_init (nst_in, nend_in)
 format_text_init generates the header records to be output in the integrals file. More...
 
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. More...
 
integer function, private get_field_index (name)
 get_field_index returns returns the index associated with an integral name. More...
 
type(time_type) function, private set_axis_time (atime, units)
 Function to convert input time to a time_type. More...
 
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. More...
 
subroutine, private write_field_averages (Time)
 Subroutine to sum multiple fields, average them and then write the result to an output file. More...
 

Variables

type(time_typealarm_interval
 time interval between writing integrals
 
real(r8_kind), dimension(:,:), allocatable area
 area of each grid box
 
integer diag_unit = 0
 unit number for output file
 
logical do_format_data = .true.
 a data format needs to be generated ?
 
integer, dimension(max_num_field) field_count
 number of values in integrand i
 
character(len=16), dimension(max_num_field) field_format
 output format for integral i
 
character(len=max_len_name), dimension(max_num_field) field_name
 name associated with integral i
 
integer field_size
 number of columns on global domain
 
real(r8_kind), dimension(max_num_field) field_sum
 integrand for integral i
 
integer fields_per_print_line = 4
 number of fields to write per line
 
character(len=fms_file_len) file_name = ' '
 optional integrals output file name
 
character(len=160) format_data
 format statement for data output
 
character(len=160) format_text
 format statement for header
 
integer idim
 x dimension of grid on local processor
 
integer jdim
 y dimension of grid on local processor
 
integer, parameter max_len_name = 12
 maximum length of name associated with integral
 
integer, parameter max_num_field = 32
 maximum number of integrals allowed
 
logical module_is_initialized = .false.
 module is initialized ?
 
integer nd
 number of characters in data format statement
 
type(time_typenext_alarm_time
 next time at which integrals are to be written
 
integer nt
 number of characters in text format statement
 
integer num_field = 0
 number of integrals that have been activated
 
real(r8_kind) output_interval = -1.0_r8_kind
 time interval at which integrals
 
logical print_header = .true.
 print a header for the integrals
 
real(r8_kind) sum_area
 surface area of globe
 
type(time_typetime_init_save
 initial time associated with experiment; used as a base for defining time
 
character(len=8) time_units = 'hours'
 time units associated with
 
type(time_typezero_time
 time_type variable set to (0,0); used as flag to indicate integrals are not being output
 

Detailed Description

This module computes and outputs global and / or hemispheric physics integrals.

Author
Fei Liu Fei.L.nosp@m.iu@n.nosp@m.oaa.g.nosp@m.ov

Data Type Documentation

◆ diag_integral_mod::sum_diag_integral_field

interface diag_integral_mod::sum_diag_integral_field

Perform a summation of the named integral field.

This interface can be called in any one of three ways:

call sum_diag_integral_field (name, data, is, js)
call sum_diag_integral_field (name, data, wt, is, js)
call sum_diag_integral_field (name, data, is, ie, js, je)

in the first option above, data may be either

real, intent(in) :: data(:,:) ![ sum_field_2d ]
real, intent(in) :: data(:,:,:) ![ sum_field_3d ]

Parameters:

character(len=*), intent(in) :: name
real, intent(in) :: wt(:,:,:)
integer, optional, intent(in) :: is, ie, js, je
Parameters
[in]<name>name associated with integral
[in]<data>field of integrands to be summed over
[in]<wt>vertical weighting factor to be applied to integrands when summing
[in]<is,ie,js,je>starting/ending i,j indices over which summation is to occur

Definition at line 103 of file diag_integral.F90.

Private Member Functions

 sum_field_2d_hemi_r4
 
 sum_field_2d_hemi_r8
 
 sum_field_2d_r4
 
 sum_field_2d_r8
 
 sum_field_3d_r4
 
 sum_field_3d_r8
 
 sum_field_wght_3d_r4
 
 sum_field_wght_3d_r8
 

Function/Subroutine Documentation

◆ diag_integral_alarm()

logical function, private diag_integral_mod::diag_integral_alarm ( type (time_type), intent(in)  Time)
private

Function to check if it is time to write integrals. if not writing integrals, return.

Template:

result = diag_integral_alarm(time)

Parameters:

type (time_type), intent(in) :: Time
logical :: answer
Parameters
[in]<Time>current time
[out]<answer>
Returns
logical answer
Parameters
[in]timecurrent time

Definition at line 1023 of file diag_integral.F90.

◆ diag_integral_end()

subroutine, public diag_integral_mod::diag_integral_end ( type (time_type), intent(in)  Time)

diag_integral_end is the destructor for diag_integral_mod.

Template:

call diag_integral_end (time)

Parameters:

type (time_type), intent(in) :: Time
Parameters
[in]<Time>integral time stamp at the current time

Definition at line 514 of file diag_integral.F90.

◆ diag_integral_field_init()

subroutine, public diag_integral_mod::diag_integral_field_init ( character(len=*), intent(in)  name,
character(len=*), intent(in)  format 
)

diag_integral_field_init registers and intializes an integral field

Template:

call diag_integral_field_init (name, format)

Parameters:

character(len=*), intent(in) :: name, format
Parameters
[in]<name>Name of the field to be integrated
[in]<format>Output format of the field to be integrated
[in]nameName of the field to be integrated
[in]formatOutput format of the field to be integrated

Definition at line 393 of file diag_integral.F90.

◆ diag_integral_init()

subroutine, public diag_integral_mod::diag_integral_init ( type (time_type), intent(in), optional  Time_init,
type (time_type), intent(in), optional  Time,
class(*), dimension(:,:), intent(in), optional  blon,
class(*), dimension(:,:), intent(in), optional  blat,
class(*), dimension(:,:), intent(in), optional  area_in 
)

diag_integral_init is the constructor for diag_integral_mod.

Template:

call diag_integral_init (time_init, time, blon, blat)

Parameters:

type (time_type), intent(in), optional :: Time_init, Time
real,dimension(:,:), intent(in), optional :: blon, blat, area_in
Parameters
[in]<Time_init>Initial time to start the integral
[in]<Time>current time
[in]<latb>array of model latitudes at cell boundaries [radians]
[in]<lonb>array of model longitudes at cell boundaries [radians]
[in]time_initInitial time to start the integral
[in]timecurrent time
[in]blonarray of model latitudes at cell boundaries [radians]
[in]blatarray of model longitudes at cell boundaries [radians]

Definition at line 240 of file diag_integral.F90.

◆ diag_integral_output()

subroutine, public diag_integral_mod::diag_integral_output ( type (time_type), intent(in)  Time)

diag_integral_output determines if this is a timestep on which integrals are to be written. if not, it returns; if so, it calls write_field_averages.

Template:

call diag_integral_output (time)

Parameters:

type (time_type), intent(in) :: Time
Parameters
[in]<Time>integral time stamp at the current time

Definition at line 469 of file diag_integral.F90.

◆ ensemble_file_name()

character (len=fms_file_len) function diag_integral_mod::ensemble_file_name ( character (len=*), intent(inout)  fname)
private

Adds .ens_## to the diag_integral.out file name.

Returns
character array updated_file_name

Make sure the file name short enough to handle adding the ensemble number

Get the ensemble ID and convert it to a string

Add a 0 if the ensemble is less than 10 (2 digit)

Insert the ens_ string in the correct location Loop through to find the last period

Add to the end if there is no period in the file name

Definition at line 1080 of file diag_integral.F90.

◆ format_data_init()

subroutine, private diag_integral_mod::format_data_init ( integer, intent(in), optional  nst_in,
integer, intent(in), optional  nend_in 
)
private

format_text_init generates the format to be output in the integrals file.

Template:

call format_data_init (nst_in, nend_in)

Parameters:

integer, intent(in), optional :: nst_in, nend_in
Parameters
[in]<nst_in,nend_in>starting/ending integral index which will be included in this format statement
[in]nst_instarting/ending integral index which will be included in this format statement
[in]nend_instarting/ending integral index which will be included in this format statement

Definition at line 899 of file diag_integral.F90.

◆ format_text_init()

subroutine, private diag_integral_mod::format_text_init ( integer, intent(in), optional  nst_in,
integer, intent(in), optional  nend_in 
)
private

format_text_init generates the header records to be output in the integrals file.

Template:

call format_text_init (nst_in, nend_in)

Parameters:

integer, intent(in), optional :: nst_in, nend_in
Parameters
[in]<nst_in,nend_in>starting/ending integral index which will be included in this format statement
[in]nst_instarting/ending integral index which will be included in this format statement
[in]nend_instarting/ending integral index which will be included in this format statement

Definition at line 810 of file diag_integral.F90.

◆ get_axis_time()

real(r8_kind) function, private diag_integral_mod::get_axis_time ( type(time_type), intent(in)  Time,
character(len=*), intent(in)  units 
)
private

Function to convert the time_type input variable into units of units and returns it in atime.

Template:

atime = get_axis_time(time, units)

Parameters:

type(time_type), intent(in) :: Time
character(len=*), intent(in) :: units
real :: atime
Parameters
[in]<Time>integral time stamp
[in]<units>input units of time_type
[out]<atime>
Returns
real atime
Parameters
[in]timeintegral time stamp
[in]unitsinput units of time_type

Definition at line 976 of file diag_integral.F90.

◆ get_field_index()

integer function, private diag_integral_mod::get_field_index ( character(len=*), intent(in)  name)
private

get_field_index returns returns the index associated with an integral name.

Template:

index = get_field_index(name)

Parameters:

character(len=*), intent(in) :: name
integer :: index
Parameters
[in]<name>Name associated with an integral
[out]<index>
[in]nameName associated with an integral

Definition at line 639 of file diag_integral.F90.

◆ set_axis_time()

type(time_type) function, private diag_integral_mod::set_axis_time ( real(r8_kind), intent(in)  atime,
character(len=*), intent(in)  units 
)
private

Function to convert input time to a time_type.

Template:

time = set_axis_time(atime, units)

Parameters:

real, intent(in) :: atime
character(len=*), intent(in) :: units
type(time_type) :: Time
Parameters
[in]<atime>integral time stamp at the current time
[in]<units>input units, not used
[out]<Time>
[in]atimeintegral time stamp at the current time
[in]unitsinput units, not used

Definition at line 581 of file diag_integral.F90.

◆ vert_diag_integral()

real(r8_kind) function, dimension (size(field_data,1),size(field_data,2)), private diag_integral_mod::vert_diag_integral ( real(r8_kind), dimension (:,:,:), intent(in)  field_data,
real(r8_kind), dimension (:,:,:), intent(in)  wt 
)
private

Function to perform a weighted integral in the vertical direction of a 3d data field.

Template:

data2 = vert_diag_integral(field_data, wt)

Parameters:

real, dimension (:,:,:), intent(in) :: field_data, wt
real, dimension (size(field_data,1),size(field_data,2)) :: data2
Parameters
[in]<field_data>integral field data arrays
[in]<wt>integral field weighting functions
[out]<data2>
Returns
real array data2
Parameters
[in]field_dataintegral field data arrays
[in]wtintegral field weighting functions

Definition at line 1057 of file diag_integral.F90.

◆ write_field_averages()

subroutine, private diag_integral_mod::write_field_averages ( type (time_type), intent(in)  Time)
private

Subroutine to sum multiple fields, average them and then write the result to an output file.

Template:

call write_field_averages (time)

Parameters:

type (time_type), intent(in) :: Time
Parameters
[in]<Time>integral time stamp at the current time
[in]timeintegral time stamp at the current time

Definition at line 692 of file diag_integral.F90.