FMS
2024.03
Flexible Modeling System
|
fms_diag_reduction_methods_mod contains routines that are meant to be used for error checking and setting up to do the reduction methods More...
Data Types | |
interface | do_time_max |
Does the time_max reduction method. See include/fms_diag_reduction_methods.inc. More... | |
interface | do_time_min |
Does the time_min reduction method. See include/fms_diag_reduction_methods.inc. More... | |
interface | do_time_sum_update |
Sum update updates the buffer for any reductions that involve summation (ie. time_sum, avg, rms, pow) TODO This needs to be extended to integers. More... | |
interface | time_update_done |
Finishes a reduction that involves an average (ie. time_avg, rms, pow) This takes the average at the end of the time step. More... | |
Functions/Subroutines | |
pure character(len=128) function, public | check_indices_order (is_in, ie_in, js_in, je_in) |
Checks improper combinations of is, ie, js, and je. More... | |
do_time_max_r4 | |
do_time_max_r8 | |
do_time_min_r4 | |
do_time_min_r8 | |
do_time_none_r4 | |
do_time_none_r8 | |
do_time_sum_update_r4 | |
do_time_sum_update_r8 | |
logical function, dimension(:,:,:,:), allocatable, public | init_mask (rmask, mask, field) |
Sets the logical mask based on mask or rmask. More... | |
pure real(kind=r8_kind) function, public | set_weight (weight) |
Sets the weight based on the weight passed into send_data (1.0_r8_kind if the weight is not passed in) The weight will be saved as an r8 and converted to r4 as needed. More... | |
sum_update_done_r4 | |
sum_update_done_r8 | |
fms_diag_reduction_methods_mod contains routines that are meant to be used for error checking and setting up to do the reduction methods
interface fms_diag_reduction_methods_mod::do_time_max |
Does the time_max reduction method. See include/fms_diag_reduction_methods.inc.
Definition at line 56 of file fms_diag_reduction_methods.F90.
Private Member Functions | |
do_time_max_r4 | |
do_time_max_r8 | |
interface fms_diag_reduction_methods_mod::do_time_min |
Does the time_min reduction method. See include/fms_diag_reduction_methods.inc.
Definition at line 50 of file fms_diag_reduction_methods.F90.
Private Member Functions | |
do_time_min_r4 | |
do_time_min_r8 | |
interface fms_diag_reduction_methods_mod::do_time_sum_update |
Sum update updates the buffer for any reductions that involve summation (ie. time_sum, avg, rms, pow) TODO This needs to be extended to integers.
Definition at line 63 of file fms_diag_reduction_methods.F90.
Private Member Functions | |
do_time_sum_update_r4 | |
do_time_sum_update_r8 | |
interface fms_diag_reduction_methods_mod::time_update_done |
Finishes a reduction that involves an average (ie. time_avg, rms, pow) This takes the average at the end of the time step.
Definition at line 70 of file fms_diag_reduction_methods.F90.
Private Member Functions | |
sum_update_done_r4 | |
sum_update_done_r8 | |
pure character(len=128) function, public fms_diag_reduction_methods_mod::check_indices_order | ( | integer, intent(in), optional | is_in, |
integer, intent(in), optional | ie_in, | ||
integer, intent(in), optional | js_in, | ||
integer, intent(in), optional | je_in | ||
) |
Checks improper combinations of is, ie, js, and je.
[in] | je_in | Indices passed to fms_diag_accept_data() |
Definition at line 89 of file fms_diag_reduction_methods.F90.
logical function, dimension(:,:,:,:), allocatable, public fms_diag_reduction_methods_mod::init_mask | ( | class(*), dimension(:,:,:,:), intent(in), allocatable | rmask, |
logical, dimension(:,:,:,:), intent(in), allocatable | mask, | ||
class(*), dimension(:,:,:,:), intent(in) | field | ||
) |
Sets the logical mask based on mask or rmask.
[in] | mask | The location of the mask |
[in] | rmask | The masking values |
[in] | field | Field_data |
Definition at line 120 of file fms_diag_reduction_methods.F90.
pure real(kind=r8_kind) function, public fms_diag_reduction_methods_mod::set_weight | ( | class(*), intent(in), optional | weight | ) |
Sets the weight based on the weight passed into send_data (1.0_r8_kind if the weight is not passed in) The weight will be saved as an r8 and converted to r4 as needed.
[in] | weight | The weight use when averaging |
Definition at line 147 of file fms_diag_reduction_methods.F90.