Skip to content
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 1x1 and 2x2 blocksize error #618

Merged
merged 1 commit into from
Jul 16, 2021
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
1 change: 1 addition & 0 deletions cicecore/cicedynB/analysis/ice_diagnostics.F90
Original file line number Diff line number Diff line change
Expand Up @@ -825,6 +825,7 @@ subroutine runtime_diags (dt)
enddo
endif
endif
psalt(n) = c0
if (vice(i,j,iblk) /= c0) psalt(n) = work2(i,j,iblk)/vice(i,j,iblk)
pTsfc(n) = trcr(i,j,nt_Tsfc,iblk) ! ice/snow sfc temperature
pevap(n) = evap(i,j,iblk)*dt/rhoi ! sublimation/condensation
Expand Down
136 changes: 13 additions & 123 deletions cicecore/cicedynB/infrastructure/comm/mpi/ice_gather_scatter.F90
Original file line number Diff line number Diff line change
Expand Up @@ -636,6 +636,7 @@ subroutine gather_global_ext_dbl(ARRAY_G, ARRAY, dst_task, src_dist, spc_val)
else
special_value = spval_dbl
endif
ARRAY_G = special_value

nx = nx_global + 2*nghost
ny = ny_global + 2*nghost
Expand Down Expand Up @@ -744,92 +745,7 @@ subroutine gather_global_ext_dbl(ARRAY_G, ARRAY, dst_task, src_dist, spc_val)
endif
endif

!*** fill land blocks with special values

else if (src_dist%blockLocation(n) == 0) then

this_block = get_block(n,n)

! interior
do j=this_block%jlo,this_block%jhi
do i=this_block%ilo,this_block%ihi
ARRAY_G(this_block%i_glob(i)+nghost, &
this_block%j_glob(j)+nghost) = special_value
end do
end do

#ifdef CICE_IN_NEMO
!echmod: this code is temporarily wrapped for nemo pending further testing elsewhere
! fill ghost cells
if (this_block%jblock == 1) then
! south block
do j=1, nghost
do i=this_block%ilo,this_block%ihi
ARRAY_G(this_block%i_glob(i)+nghost,j) = special_value
end do
end do
if (this_block%iblock == 1) then
! southwest corner
do j=1, nghost
do i=1, nghost
ARRAY_G(i,j) = special_value
end do
end do
endif
endif
if (this_block%jblock == nblocks_y) then
! north block
do j=1, nghost
do i=this_block%ilo,this_block%ihi
ARRAY_G(this_block%i_glob(i)+nghost, &
ny_global + nghost + j) = special_value
end do
end do
if (this_block%iblock == nblocks_x) then
! northeast corner
do j=1, nghost
do i=1, nghost
ARRAY_G(nx-i+1, ny-j+1) = special_value
end do
end do
endif
endif
if (this_block%iblock == 1) then
! west block
do j=this_block%jlo,this_block%jhi
do i=1, nghost
ARRAY_G(i,this_block%j_glob(j)+nghost) = special_value
end do
end do
if (this_block%jblock == nblocks_y) then
! northwest corner
do j=1, nghost
do i=1, nghost
ARRAY_G(i, ny-j+1) = special_value
end do
end do
endif
endif
if (this_block%iblock == nblocks_x) then
! east block
do j=this_block%jlo,this_block%jhi
do i=1, nghost
ARRAY_G(nx_global + nghost + i, &
this_block%j_glob(j)+nghost) = special_value
end do
end do
if (this_block%jblock == 1) then
! southeast corner
do j=1, nghost
do i=1, nghost
ARRAY_G( nx-i+1,j) = special_value
end do
end do
endif
endif
#endif

endif
endif ! src_dist%blockLocation

end do

Expand Down Expand Up @@ -939,7 +855,7 @@ subroutine gather_global_ext_dbl(ARRAY_G, ARRAY, dst_task, src_dist, spc_val)
!
!-----------------------------------------------------------------------

else
else ! master task

allocate(snd_request(nblocks_tot), &
snd_status (MPI_STATUS_SIZE, nblocks_tot))
Expand All @@ -960,7 +876,7 @@ subroutine gather_global_ext_dbl(ARRAY_G, ARRAY, dst_task, src_dist, spc_val)
call MPI_WAITALL(nsends, snd_request, snd_status, ierr)
deallocate(snd_request, snd_status)

endif
endif ! master task

if (add_mpi_barriers) then
call ice_barrier()
Expand Down Expand Up @@ -1028,8 +944,9 @@ subroutine gather_global_ext_int(ARRAY_G, ARRAY, dst_task, src_dist, spc_val)
if (present(spc_val)) then
special_value = spc_val
else
special_value = 0 !MHRI NOTE: 0,1,-999,???
special_value = -9999
endif
ARRAY_G = special_value

nx = nx_global + 2*nghost
ny = ny_global + 2*nghost
Expand Down Expand Up @@ -1138,21 +1055,7 @@ subroutine gather_global_ext_int(ARRAY_G, ARRAY, dst_task, src_dist, spc_val)
endif
endif

!*** fill land blocks with special values

else if (src_dist%blockLocation(n) == 0) then

this_block = get_block(n,n)

! interior
do j=this_block%jlo,this_block%jhi
do i=this_block%ilo,this_block%ihi
ARRAY_G(this_block%i_glob(i)+nghost, &
this_block%j_glob(j)+nghost) = special_value
end do
end do

endif
endif ! src_dist%blockLocation

end do

Expand Down Expand Up @@ -1262,7 +1165,7 @@ subroutine gather_global_ext_int(ARRAY_G, ARRAY, dst_task, src_dist, spc_val)
!
!-----------------------------------------------------------------------

else
else ! master task

allocate(snd_request(nblocks_tot), &
snd_status (MPI_STATUS_SIZE, nblocks_tot))
Expand All @@ -1283,7 +1186,7 @@ subroutine gather_global_ext_int(ARRAY_G, ARRAY, dst_task, src_dist, spc_val)
call MPI_WAITALL(nsends, snd_request, snd_status, ierr)
deallocate(snd_request, snd_status)

endif
endif ! master task

if (add_mpi_barriers) then
call ice_barrier()
Expand Down Expand Up @@ -1353,6 +1256,7 @@ subroutine gather_global_ext_log(ARRAY_G, ARRAY, dst_task, src_dist, spc_val)
else
special_value = .false. !MHRI NOTE: .true./.false. ???
endif
ARRAY_G = special_value

nx = nx_global + 2*nghost
ny = ny_global + 2*nghost
Expand Down Expand Up @@ -1461,21 +1365,7 @@ subroutine gather_global_ext_log(ARRAY_G, ARRAY, dst_task, src_dist, spc_val)
endif
endif

!*** fill land blocks with special values

else if (src_dist%blockLocation(n) == 0) then

this_block = get_block(n,n)

! interior
do j=this_block%jlo,this_block%jhi
do i=this_block%ilo,this_block%ihi
ARRAY_G(this_block%i_glob(i)+nghost, &
this_block%j_glob(j)+nghost) = special_value
end do
end do

endif
endif ! src_dist%blockLocation

end do

Expand Down Expand Up @@ -1585,7 +1475,7 @@ subroutine gather_global_ext_log(ARRAY_G, ARRAY, dst_task, src_dist, spc_val)
!
!-----------------------------------------------------------------------

else
else ! master task

allocate(snd_request(nblocks_tot), &
snd_status (MPI_STATUS_SIZE, nblocks_tot))
Expand All @@ -1606,7 +1496,7 @@ subroutine gather_global_ext_log(ARRAY_G, ARRAY, dst_task, src_dist, spc_val)
call MPI_WAITALL(nsends, snd_request, snd_status, ierr)
deallocate(snd_request, snd_status)

endif
endif ! master task

if (add_mpi_barriers) then
call ice_barrier()
Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -373,6 +373,7 @@ subroutine gather_global_ext_dbl(ARRAY_G, ARRAY, dst_task, src_dist, spc_val)
else
special_value = spval_dbl
endif
ARRAY_G = special_value

nx = nx_global + 2*nghost
ny = ny_global + 2*nghost
Expand Down Expand Up @@ -477,16 +478,7 @@ subroutine gather_global_ext_dbl(ARRAY_G, ARRAY, dst_task, src_dist, spc_val)
endif
endif

else !*** fill land blocks with special values

do j=this_block%jlo,this_block%jhi
do i=this_block%ilo,this_block%ihi
ARRAY_G(this_block%i_glob(i)+nghost, &
this_block%j_glob(j)+nghost) = special_value
end do
end do

endif
endif ! src_dist%blockLocation

end do

Expand Down Expand Up @@ -537,8 +529,9 @@ subroutine gather_global_ext_int(ARRAY_G, ARRAY, dst_task, src_dist, spc_val)
if (present(spc_val)) then
special_value = spc_val
else
special_value = 0 !MHRI: 0,1,999,-999 ??
special_value = -9999
endif
ARRAY_G = special_value

nx = nx_global + 2*nghost
ny = ny_global + 2*nghost
Expand Down Expand Up @@ -643,16 +636,7 @@ subroutine gather_global_ext_int(ARRAY_G, ARRAY, dst_task, src_dist, spc_val)
endif
endif

else !*** fill land blocks with special values

do j=this_block%jlo,this_block%jhi
do i=this_block%ilo,this_block%ihi
ARRAY_G(this_block%i_glob(i)+nghost, &
this_block%j_glob(j)+nghost) = special_value
end do
end do

endif
endif ! src_dist%blockLocation

end do

Expand Down Expand Up @@ -705,6 +689,7 @@ subroutine gather_global_ext_log(ARRAY_G, ARRAY, dst_task, src_dist, spc_val)
else
special_value = .false. !MHRI: true/false
endif
ARRAY_G = special_value

nx = nx_global + 2*nghost
ny = ny_global + 2*nghost
Expand Down Expand Up @@ -809,16 +794,7 @@ subroutine gather_global_ext_log(ARRAY_G, ARRAY, dst_task, src_dist, spc_val)
endif
endif

else !*** fill land blocks with special values

do j=this_block%jlo,this_block%jhi
do i=this_block%ilo,this_block%ihi
ARRAY_G(this_block%i_glob(i)+nghost, &
this_block%j_glob(j)+nghost) = special_value
end do
end do

endif
endif ! src_dist%blockLocation

end do

Expand Down
2 changes: 1 addition & 1 deletion cicecore/cicedynB/infrastructure/ice_domain.F90
Original file line number Diff line number Diff line change
Expand Up @@ -441,7 +441,7 @@ subroutine init_domain_distribution(KMTG,ULATG)
!----------------------------------------------------------------------

if (distribution_wght == 'latitude') then
flat = NINT(abs(ULATG*rad_to_deg), int_kind) ! linear function
flat = max(NINT(abs(ULATG*rad_to_deg), int_kind),1) ! linear function
else
flat = 1
endif
Expand Down
1 change: 1 addition & 0 deletions configuration/scripts/tests/reprosum_suite.ts
Original file line number Diff line number Diff line change
Expand Up @@ -8,4 +8,5 @@ logbfb gx3 1x20x5x29x80 dsectrobin,diag1,short,reprosum l
logbfb gx3 8x2x8x10x20 droundrobin,diag1,reprosum logbfb_gx3_4x2x25x29x4_diag1_dslenderX2_reprosum
logbfb gx3 6x2x50x58x1 droundrobin,diag1,reprosum logbfb_gx3_4x2x25x29x4_diag1_dslenderX2_reprosum
logbfb gx3 6x2x4x29x18 dspacecurve,diag1,maskhalo,reprosum logbfb_gx3_4x2x25x29x4_diag1_dslenderX2_reprosum
logbfb gx3 17x2x1x1x800 droundrobin,diag1,maskhalo,reprosum logbfb_gx3_4x2x25x29x4_diag1_dslenderX2_reprosum
#logbfb gx3 8x2x8x10x20 droundrobin,diag1 logbfb_gx3_4x2x25x29x4_diag1_dslenderX2