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

Fixed bug where model seg faults for dyn_npes<npes #1393

Merged
merged 1 commit into from
Apr 14, 2017
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
6 changes: 5 additions & 1 deletion components/cam/src/dynamics/se/dyn_comp.F90
Original file line number Diff line number Diff line change
Expand Up @@ -50,6 +50,8 @@ Module dyn_comp
! !REVISION HISTORY:
!
! JPE 06.05.31: created
! Aaron Donahue 17.04.11: Fixed bug in write_grid_mapping which caused
! a segmentation fault when dyn_npes<npes
!
!----------------------------------------------------------------------

Expand Down Expand Up @@ -444,7 +446,9 @@ subroutine write_grid_mapping(par, elem)
ierr = pio_def_var(nc, 'element_corners', PIO_INT, (/dim1,dim2/),vid)

ierr = pio_enddef(nc)
call createmetadata(par, elem, subelement_corners)
if (iam<par%nprocs) then
call createmetadata(par, elem, subelement_corners)
end if

jj=0
do cc=0,3
Expand Down
9 changes: 8 additions & 1 deletion components/cam/src/dynamics/se/gravity_waves_sources.F90
Original file line number Diff line number Diff line change
Expand Up @@ -100,13 +100,18 @@ subroutine compute_frontogenesis(frontgf,frontga,tl,elem,ederiv,hybrid,nets,nete
! Change by Santos, 10 Aug 2011:
! Integrated into gravity_waves_sources module, several arguments made global
! to prevent repeated allocation/initialization
! Change by Aaron Donahue, April 2017:
! Fixed bug where boundary information was called for processors not associated
! with dynamics when dyn_npes<npes
!
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
use physical_constants, only : kappa
use derivative_mod, only : gradient_sphere, ugradv_sphere
use edge_mod, only : edgevpack, edgevunpack
use bndry_mod, only : bndry_exchangev
use dyn_comp, only : hvcoord
use spmd_utils, only : iam
use parallel_mod, only : par
implicit none
type (hybrid_t) , intent(in) :: hybrid
type (element_t) , intent(inout), target :: elem(:)
Expand Down Expand Up @@ -150,7 +155,9 @@ subroutine compute_frontogenesis(frontgf,frontga,tl,elem,ederiv,hybrid,nets,nete
call edgeVpack(edge3, frontgf(:,:,:,ie),nlev,0,ie)
call edgeVpack(edge3, gradth(:,:,:,:,ie),2*nlev,nlev,ie)
enddo
call bndry_exchangeV(hybrid,edge3)
if (iam<par%nprocs) then
call bndry_exchangeV(hybrid,edge3)
end if
do ie=nets,nete
call edgeVunpack(edge3, frontgf(:,:,:,ie),nlev,0,ie)
call edgeVunpack(edge3, gradth(:,:,:,:,ie),2*nlev,nlev,ie)
Expand Down