-
Notifications
You must be signed in to change notification settings - Fork 373
[IO] Adding functions to estimate max byte size of all fields in a pool #1217
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
base: master
Are you sure you want to change the base?
Changes from 17 commits
4abf490
362c44e
68271e8
7fba511
f94b526
d438c29
6ac4a3b
7b4aa9b
2107e91
d5972cd
1a784a5
ea8d1c8
2d36c22
54e6776
a12be08
cd764aa
04e2291
cbcb236
cfd894e
File filter
Filter by extension
Conversations
Jump to
Diff view
Diff view
There are no files selected for viewing
| Original file line number | Diff line number | Diff line change |
|---|---|---|
|
|
@@ -3084,7 +3084,7 @@ subroutine write_stream(manager, stream, blockID, timeLevel, mgLevel, forceWrite | |
| character (len=StrKIND) :: now_string, time_string | ||
| character (len=StrKIND) :: temp_filename, actualWhen | ||
| character (len=StrKIND) :: err_string | ||
| logical :: ringing_alarm, recordSeek, swapRecords | ||
| logical :: ringing_alarm, recordSeek, swapRecords, use_cdf2 | ||
| logical :: clobberRecords, clobberFiles, truncateFiles | ||
| integer :: maxRecords, tempRecord | ||
| integer :: local_ierr, threadNum | ||
|
|
@@ -3182,6 +3182,7 @@ subroutine write_stream(manager, stream, blockID, timeLevel, mgLevel, forceWrite | |
| ! | ||
| ! Build stream from pools of fields and attributes | ||
| ! | ||
| call check_max_var_size(stream % field_pool, manager % allFields, timeLevel, mgLevel, use_cdf2) | ||
| allocate(stream % stream) | ||
| call MPAS_createStream(stream % stream, manager % ioContext, stream % filename, stream % io_type, MPAS_IO_WRITE, & | ||
| precision=stream % precision, clobberRecords=clobberRecords, & | ||
|
|
@@ -4325,6 +4326,7 @@ end subroutine gen_random | |
| timeLevel = 1 | ||
| end if | ||
|
|
||
|
|
||
| select case (info % fieldType) | ||
| case (MPAS_POOL_REAL) | ||
| select case (info % nDims) | ||
|
|
@@ -4495,6 +4497,123 @@ subroutine update_stream(stream, allFields, timeLevelIn, mgLevelIn, ierr) !{{{ | |
| end subroutine update_stream !}}} | ||
|
|
||
|
|
||
| subroutine check_max_var_size(field_pool, allFields, timeLevelIn, mgLevelIn, use_cdf2) !{{{ | ||
| use iso_c_binding, only : c_sizeof | ||
|
|
||
| implicit none | ||
|
|
||
| type (mpas_pool_type), intent(inout) :: field_pool | ||
| type (MPAS_Pool_type), intent(in) :: allFields | ||
| integer, intent(in) :: timeLevelIn | ||
| integer, intent(in) :: mgLevelIn | ||
| logical, intent(out) :: use_cdf2 | ||
|
|
||
| type (MPAS_Pool_iterator_type) :: itr | ||
| type (mpas_pool_field_info_type) :: info | ||
| integer :: timeLevel | ||
|
|
||
| type (field5DReal), pointer :: real5d | ||
| type (field4DReal), pointer :: real4d | ||
| type (field3DReal), pointer :: real3d | ||
| type (field2DReal), pointer :: real2d | ||
| type (field1DReal), pointer :: real1d | ||
| type (field0DReal), pointer :: real0d | ||
|
|
||
| type (field3DInteger), pointer :: int3d | ||
| type (field2DInteger), pointer :: int2d | ||
| type (field1DInteger), pointer :: int1d | ||
| type (field0DInteger), pointer :: int0d | ||
|
|
||
| type (field1DChar), pointer :: char1d | ||
| type (field0DChar), pointer :: char0d | ||
|
|
||
| integer :: max_field_size, field_size | ||
| integer :: tmp_int | ||
| real(kind=RKIND) :: tmp_real | ||
| character :: tmp_char | ||
| integer(kind=I8KIND) :: max_field_bytes, field_bytes | ||
| integer(kind=I8KIND), parameter :: int_size = c_sizeof(tmp_int) | ||
| integer(kind=I8KIND), parameter :: real_size = c_sizeof(tmp_real) | ||
|
||
| integer(kind=I8KIND), parameter :: char_size = c_sizeof(tmp_char) | ||
| integer(kind=I8KIND), parameter :: cdf2_max_bytes = 2**32 - 4 | ||
| character(len=StrKIND):: message | ||
|
|
||
| call mpas_pool_begin_iteration(field_pool) | ||
| max_field_bytes = 0 | ||
| field_bytes = 0 | ||
| do while ( mpas_pool_get_next_member(field_pool, itr) ) | ||
|
|
||
| if (itr % memberType == MPAS_POOL_CONFIG) then | ||
|
|
||
| ! To avoid accidentally matching in case statements below... | ||
| info % fieldType = -1 | ||
| call mpas_pool_get_field_info(allFields, itr % memberName, info) | ||
| ! Set time level to read | ||
mgduda marked this conversation as resolved.
Outdated
Show resolved
Hide resolved
|
||
| if (info % nTimeLevels >= timeLevelIn) then | ||
| timeLevel = timeLevelIn | ||
| else | ||
| timeLevel = 1 | ||
| end if | ||
| call mpas_log_write('In check_max_var_size, field '//trim(itr % memberName)//' ndims: $i',intArgs=(/info % nDims/)) | ||
| select case (info % fieldType) | ||
| case (MPAS_POOL_REAL) | ||
| select case (info % nDims) | ||
| case (0) | ||
| field_size = 1 | ||
| case (1) | ||
| call mpas_pool_get_field(allFields, itr % memberName, real1d, timeLevel) | ||
| field_size = global_dim_size(real1d % block, real1d % dimNames, real1d % dimSizes, info % nDims) | ||
| case (2) | ||
| call mpas_pool_get_field(allFields, itr % memberName, real2d, timeLevel) | ||
| field_size = global_dim_size(real2d % block, real2d % dimNames, real2d % dimSizes, info % nDims) | ||
| case (3) | ||
| call mpas_pool_get_field(allFields, itr % memberName, real3d, timeLevel) | ||
| field_size = global_dim_size(real3d % block, real3d % dimNames, real3d % dimSizes, info % nDims) | ||
| case (4) | ||
| call mpas_pool_get_field(allFields, itr % memberName, real4d, timeLevel) | ||
| field_size = global_dim_size(real4d % block, real4d % dimNames, real4d % dimSizes, info % nDims) | ||
| case (5) | ||
| call mpas_pool_get_field(allFields, itr % memberName, real5d, timeLevel) | ||
| field_size = global_dim_size(real5d % block, real5d % dimNames, real5d % dimSizes, info % nDims) | ||
| end select | ||
| field_bytes = int(field_size, kind=I8KIND) * real_size | ||
|
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. Perhaps 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. Isn't there a type conversion involved whether it's in this function or inside 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. I've now changed this. The only type conversion happens in |
||
| case (MPAS_POOL_INTEGER) | ||
| select case (info % nDims) | ||
| case (0) | ||
| field_size = 1 | ||
| case (1) | ||
| call mpas_pool_get_field(allFields, itr % memberName, int1d, timeLevel) | ||
| field_size = global_dim_size(int1d % block, int1d % dimNames, int1d % dimSizes, info % nDims) | ||
| case (2) | ||
| call mpas_pool_get_field(allFields, itr % memberName, int2d, timeLevel) | ||
| field_size = global_dim_size(int2d % block, int2d % dimNames, int2d % dimSizes, info % nDims) | ||
| case (3) | ||
| call mpas_pool_get_field(allFields, itr % memberName, int3d, timeLevel) | ||
| field_size = global_dim_size(int3d % block, int3d % dimNames, int3d % dimSizes, info % nDims) | ||
| end select | ||
| field_bytes = int(field_size, kind=I8KIND) * int_size | ||
| case (MPAS_POOL_CHARACTER) | ||
| select case (info % nDims) | ||
| case (0) | ||
| field_size = 1 | ||
| case (1) | ||
| ! call mpas_pool_get_field(allFields, itr % memberName, char1d, timeLevel) | ||
| call mpas_log_write('In check_max_var_size, unsupported type field1DChar.', MPAS_LOG_ERR) | ||
| end select | ||
| field_bytes = int(field_size, kind=I8KIND) * char_size | ||
| end select | ||
| max_field_bytes = merge(field_bytes, max_field_bytes, field_bytes > max_field_bytes) | ||
| write(message,fmt='(A,i14,A,i18,A,i18)') 'check_max_var_size.. field_bytes =',field_bytes,' max_field_bytes =',max_field_bytes, & | ||
| & ' cdf2_max_bytes =',cdf2_max_bytes | ||
| call mpas_log_write(message) | ||
| end if | ||
| end do | ||
| use_cdf2 = cdf2_max_bytes > max_field_bytes | ||
| call mpas_log_write('In check_max_var_size use_cdf2: $l',logicArgs=(/use_cdf2/)) | ||
|
|
||
| end subroutine check_max_var_size !}}} | ||
|
|
||
|
|
||
| !----------------------------------------------------------------------- | ||
| ! routine stream_active_pkg_check | ||
| ! | ||
|
|
@@ -4845,6 +4964,51 @@ logical function is_decomposed_dim(dimName) !{{{ | |
|
|
||
| end function is_decomposed_dim !}}} | ||
|
|
||
|
|
||
|
|
||
| integer function global_dim_size(block, dimNames, dimSizes, nDim) !{{{ | ||
| use iso_c_binding, only : c_sizeof | ||
mgduda marked this conversation as resolved.
Outdated
Show resolved
Hide resolved
|
||
|
|
||
| implicit none | ||
|
|
||
| character(len=*), intent(in), dimension(:) :: dimNames | ||
| integer, intent(in), dimension(:) :: dimSizes | ||
| type(block_type), intent(in) :: block | ||
| integer, intent(in) :: nDim | ||
| integer, pointer :: block_dim_size | ||
| integer :: iDim, sum_block_dim_size | ||
| character(len=StrKIND):: message | ||
|
|
||
| global_dim_size = 1 | ||
| do iDim = 1, nDim | ||
mgduda marked this conversation as resolved.
Outdated
Show resolved
Hide resolved
|
||
| if ( is_decomposed_dim(dimNames(iDim))) then | ||
| if (trim(dimNames(iDim)) == 'nCells') then | ||
| call mpas_pool_get_dimension(block % dimensions, 'nCellsSolve', block_dim_size) | ||
| else if (trim(dimNames(iDim)) == 'nEdges') then | ||
| call mpas_pool_get_dimension(block % dimensions, 'nEdgesSolve', block_dim_size) | ||
| else if (trim(dimNames(iDim)) == 'nVertices') then | ||
| call mpas_pool_get_dimension(block % dimensions, 'nVerticesSolve', block_dim_size) | ||
| else | ||
| global_dim_size = -1 | ||
| end if | ||
|
|
||
| call mpas_dmpar_sum_int(block % domain % dminfo, block_dim_size, sum_block_dim_size) | ||
| call mpas_log_write('----- from global_dim_size ----- Dimname is decomposed '//trim(dimNames(iDim))//' local size $i global size $i ',intArgs=(/block_dim_size, sum_block_dim_size/)) | ||
|
|
||
| global_dim_size = global_dim_size * sum_block_dim_size | ||
| else | ||
| global_dim_size = global_dim_size * dimSizes(iDim) | ||
mgduda marked this conversation as resolved.
Outdated
Show resolved
Hide resolved
|
||
| call mpas_log_write('----- from global_dim_size ----- Dimname is not decomposed '//trim(dimNames(iDim))//' Size $i',intArgs=(/dimSizes(iDim)/)) | ||
| end if | ||
| end do | ||
|
|
||
| call mpas_log_write('----- from global_dim_size: cumulative global_dim_size $i',intArgs=(/global_dim_size/)) | ||
|
|
||
|
|
||
| end function global_dim_size !}}} | ||
|
|
||
|
|
||
|
|
||
|
|
||
| !----------------------------------------------------------------------- | ||
| ! routine prewrite_reindex | ||
|
|
||
Uh oh!
There was an error while loading. Please reload this page.