-
Notifications
You must be signed in to change notification settings - Fork 138
New issue
Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.
By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.
Already on GitHub? Sign in to your account
fix diag_send_complete loops and add get_file_ids #1407
Merged
Merged
Changes from 1 commit
Commits
Show all changes
2 commits
Select commit
Hold shift + click to select a range
File filter
Filter by extension
Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
There are no files selected for viewing
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
Original file line number | Diff line number | Diff line change |
---|---|---|
|
@@ -46,6 +46,7 @@ module fms_diag_object_mod | |
use omp_lib | ||
#endif | ||
use mpp_domains_mod, only: domain1d, domain2d, domainUG, null_domain2d | ||
use fms_string_utils_mod, only: string | ||
use platform_mod | ||
implicit none | ||
private | ||
|
@@ -648,43 +649,58 @@ subroutine fms_diag_send_complete(this, time_step) | |
class(*), pointer :: input_data_buffer(:,:,:,:) | ||
character(len=128) :: error_string | ||
type(fmsDiagIbounds_type) :: bounds | ||
integer, dimension(:), allocatable :: file_ids !< Array of file IDs for a field | ||
logical, parameter :: DEBUG_SC = .true. !< turn on output for debugging | ||
|
||
!< Update the current model time by adding the time_step | ||
this%current_model_time = this%current_model_time + time_step | ||
|
||
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! | ||
!! In the future, this may be parallelized for offloading | ||
file_loop: do ifile = 1, size(this%FMS_diag_files) | ||
diag_file => this%FMS_diag_files(ifile) | ||
field_outer_if: if (size(diag_file%FMS_diag_file%get_field_ids()) .ge. 1) then | ||
allocate (file_field_ids(size(diag_file%FMS_diag_file%get_field_ids() ))) | ||
file_field_ids = diag_file%FMS_diag_file%get_field_ids() | ||
field_loop: do ifield = 1, size(file_field_ids) | ||
! loop through each field | ||
field_loop: do ifield = 1, size(this%FMS_diag_fields) | ||
diag_field => this%FMS_diag_fields(ifield) | ||
if(.not. diag_field%is_registered()) cycle | ||
if(DEBUG_SC) call mpp_error(NOTE, "fms_diag_send_complete:: var: "//diag_field%get_varname()) | ||
! get files the field is in | ||
allocate (file_ids(size(diag_field%get_file_ids() ))) | ||
file_ids = diag_field%get_file_ids() | ||
math = diag_field%get_math_needs_to_be_done() | ||
! if doing math loop through each file for given field | ||
if (size(file_ids) .ge. 1 .and. math) then | ||
file_loop: do ifile = 1, size(file_ids) | ||
diag_file => this%FMS_diag_files(ifile) | ||
! if the file is not allocated go away | ||
if(.not. allocated(diag_file%FMS_diag_file)) then | ||
if(DEBUG_SC) call mpp_error(NOTE, "file id:"//string(ifile)//" not allocated for field:"//diag_field%get_varname()) | ||
cycle | ||
endif | ||
! If the field is not registered go away | ||
if (.not. diag_file%FMS_diag_file%is_field_registered(ifield)) cycle | ||
|
||
diag_field => this%FMS_diag_fields(file_field_ids(ifield)) | ||
!> Check if math needs to be done | ||
math = diag_field%get_math_needs_to_be_done() | ||
calling_math: if (math) then | ||
if (.not. diag_file%FMS_diag_file%is_field_registered(ifield)) then | ||
if(DEBUG_SC) call mpp_error(NOTE, "file id:"//string(ifile)//" not registered for field:"//diag_field%get_varname()) | ||
cycle | ||
endif | ||
! Check if buffer alloc'd | ||
has_input_buff: if (diag_field%has_input_data_buffer()) then | ||
input_data_buffer => diag_field%get_data_buffer() | ||
! reset bounds, allocate output buffer, and update it with reduction | ||
call bounds%reset_bounds_from_array_4D(input_data_buffer) | ||
call this%allocate_diag_field_output_buffers(input_data_buffer, file_field_ids(ifield)) | ||
error_string = this%fms_diag_do_reduction(input_data_buffer, file_field_ids(ifield), & | ||
diag_field%get_mask(), diag_field%get_weight(), & | ||
bounds, .False., Time=this%current_model_time) | ||
call this%allocate_diag_field_output_buffers(input_data_buffer, ifield) | ||
error_string = this%fms_diag_do_reduction(input_data_buffer, ifield, & | ||
diag_field%get_mask(), diag_field%get_weight(), & | ||
bounds, .False., Time=this%current_model_time) | ||
if (trim(error_string) .ne. "") call mpp_error(FATAL, "Field:"//trim(diag_field%get_varname()//& | ||
" -"//trim(error_string))) | ||
endif calling_math | ||
!> Clean up, clean up, everybody everywhere | ||
if (associated(diag_field)) nullify(diag_field) | ||
enddo field_loop | ||
!> Clean up, clean up, everybody do your share | ||
if (allocated(file_field_ids)) deallocate(file_field_ids) | ||
endif field_outer_if | ||
enddo file_loop | ||
|
||
call this%fms_diag_do_io() | ||
" -"//trim(error_string))) | ||
endif has_input_buff | ||
There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. Can we add an else here and error out with a message saying that |
||
|
||
enddo file_loop | ||
endif | ||
!> Clean up, clean up, everybody do your share | ||
if (allocated(file_ids)) deallocate(file_ids) | ||
if (associated(diag_field)) nullify(diag_field) | ||
enddo field_loop | ||
|
||
call this%fms_diag_do_io() | ||
#endif | ||
|
||
end subroutine fms_diag_send_complete | ||
|
Oops, something went wrong.
Add this suggestion to a batch that can be applied as a single commit.
This suggestion is invalid because no changes were made to the code.
Suggestions cannot be applied while the pull request is closed.
Suggestions cannot be applied while viewing a subset of changes.
Only one suggestion per line can be applied in a batch.
Add this suggestion to a batch that can be applied as a single commit.
Applying suggestions on deleted lines is not supported.
You must change the existing code in this line in order to create a valid suggestion.
Outdated suggestions cannot be applied.
This suggestion has been applied or marked resolved.
Suggestions cannot be applied from pending reviews.
Suggestions cannot be applied on multi-line comments.
Suggestions cannot be applied while the pull request is queued to merge.
Suggestion cannot be applied right now. Please check back later.
There was a problem hiding this comment.
Choose a reason for hiding this comment
The reason will be displayed to describe this comment to others. Learn more.
I think this if statement can be moved to outside the file loop.
allocate_diag_field_output_buffers
andfms_diag_do_reduction
are looping through the number of buffers for the field. The number of buffers for the field is the same as the number of files the field is in, so i think we will be doing the reductions twice here.