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

Memory leak fix: pointers in interpolator and horiz_interp #1100

Merged
merged 5 commits into from
Jan 20, 2023
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
14 changes: 7 additions & 7 deletions horiz_interp/horiz_interp_bicubic.F90
Original file line number Diff line number Diff line change
Expand Up @@ -736,13 +736,13 @@ subroutine horiz_interp_bicubic_del( Interp )

type (horiz_interp_type), intent(inout) :: Interp

if(associated(Interp%rat_x)) deallocate ( Interp%rat_x )
if(associated(Interp%rat_y)) deallocate ( Interp%rat_y )
if(associated(Interp%lon_in)) deallocate ( Interp%lon_in )
if(associated(Interp%lat_in)) deallocate ( Interp%lat_in )
if(associated(Interp%i_lon)) deallocate ( Interp%i_lon )
if(associated(Interp%j_lat)) deallocate ( Interp%j_lat )
if(associated(Interp%wti)) deallocate ( Interp%wti )
if(allocated(Interp%rat_x)) deallocate ( Interp%rat_x )
if(allocated(Interp%rat_y)) deallocate ( Interp%rat_y )
if(allocated(Interp%lon_in)) deallocate ( Interp%lon_in )
if(allocated(Interp%lat_in)) deallocate ( Interp%lat_in )
if(allocated(Interp%i_lon)) deallocate ( Interp%i_lon )
if(allocated(Interp%j_lat)) deallocate ( Interp%j_lat )
if(allocated(Interp%wti)) deallocate ( Interp%wti )

end subroutine horiz_interp_bicubic_del

Expand Down
8 changes: 4 additions & 4 deletions horiz_interp/horiz_interp_bilinear.F90
Original file line number Diff line number Diff line change
Expand Up @@ -1240,10 +1240,10 @@ subroutine horiz_interp_bilinear_del( Interp )
!! have allocated arrays. The returned variable will contain
!! deallocated arrays

if(associated(Interp%wti)) deallocate(Interp%wti)
if(associated(Interp%wtj)) deallocate(Interp%wtj)
if(associated(Interp%i_lon)) deallocate(Interp%i_lon)
if(associated(Interp%j_lat)) deallocate(Interp%j_lat)
if(allocated(Interp%wti)) deallocate(Interp%wti)
if(allocated(Interp%wtj)) deallocate(Interp%wtj)
if(allocated(Interp%i_lon)) deallocate(Interp%i_lon)
if(allocated(Interp%j_lat)) deallocate(Interp%j_lat)

end subroutine horiz_interp_bilinear_del

Expand Down
24 changes: 12 additions & 12 deletions horiz_interp/horiz_interp_conserve.F90
Original file line number Diff line number Diff line change
Expand Up @@ -821,7 +821,7 @@ subroutine horiz_interp_conserve_version1 ( Interp, data_in, data_out, verbose,
if (present(mask_in)) then
call data_sum ( data_in(is:ie,js:je), Interp%area_src(is:ie,js:je), &
fis, fie, fjs,fje, dwtsum, wtsum, arsum, mask_in(is:ie,js:je) )
else if( ASSOCIATED(Interp%mask_in) ) then
else if( allocated(Interp%mask_in) ) then
call data_sum ( data_in(is:ie,js:je), Interp%area_src(is:ie,js:je), &
fis, fie, fjs,fje, dwtsum, wtsum, arsum, Interp%mask_in(is:ie,js:je) )
else
Expand Down Expand Up @@ -924,18 +924,18 @@ subroutine horiz_interp_conserve_del ( Interp )

select case(Interp%version)
case (1)
if(associated(Interp%area_src)) deallocate(Interp%area_src)
if(associated(Interp%area_dst)) deallocate(Interp%area_dst)
if(associated(Interp%facj)) deallocate(Interp%facj)
if(associated(Interp%jlat)) deallocate(Interp%jlat)
if(associated(Interp%faci)) deallocate(Interp%faci)
if(associated(Interp%ilon)) deallocate(Interp%ilon)
if(allocated(Interp%area_src)) deallocate(Interp%area_src)
if(allocated(Interp%area_dst)) deallocate(Interp%area_dst)
if(allocated(Interp%facj)) deallocate(Interp%facj)
if(allocated(Interp%jlat)) deallocate(Interp%jlat)
if(allocated(Interp%faci)) deallocate(Interp%faci)
if(allocated(Interp%ilon)) deallocate(Interp%ilon)
case (2)
if(associated(Interp%i_src)) deallocate(Interp%i_src)
if(associated(Interp%j_src)) deallocate(Interp%j_src)
if(associated(Interp%i_dst)) deallocate(Interp%i_dst)
if(associated(Interp%j_dst)) deallocate(Interp%j_dst)
if(associated(Interp%area_frac_dst)) deallocate(Interp%area_frac_dst)
if(allocated(Interp%i_src)) deallocate(Interp%i_src)
if(allocated(Interp%j_src)) deallocate(Interp%j_src)
if(allocated(Interp%i_dst)) deallocate(Interp%i_dst)
if(allocated(Interp%j_dst)) deallocate(Interp%j_dst)
if(allocated(Interp%area_frac_dst)) deallocate(Interp%area_frac_dst)
end select

end subroutine horiz_interp_conserve_del
Expand Down
10 changes: 5 additions & 5 deletions horiz_interp/horiz_interp_spherical.F90
Original file line number Diff line number Diff line change
Expand Up @@ -208,7 +208,7 @@ subroutine horiz_interp_spherical_new(Interp, lon_in,lat_in,lon_out,lat_out, &
endif

! allocate memory to data type
if(ASSOCIATED(Interp%i_lon)) then
if(allocated(Interp%i_lon)) then
if(size(Interp%i_lon,1) .NE. map_dst_xsize .OR. &
size(Interp%i_lon,2) .NE. map_dst_ysize ) call mpp_error(FATAL, &
'horiz_interp_spherical_mod: size(Interp%i_lon(:),1) .NE. map_dst_xsize .OR. '// &
Expand Down Expand Up @@ -503,10 +503,10 @@ subroutine horiz_interp_spherical_del( Interp )
!! must have allocated arrays. The returned variable will
!! contain deallocated arrays.

if(associated(Interp%src_dist)) deallocate(Interp%src_dist)
if(associated(Interp%num_found)) deallocate(Interp%num_found)
if(associated(Interp%i_lon)) deallocate(Interp%i_lon)
if(associated(Interp%j_lat)) deallocate(Interp%j_lat)
if(allocated(Interp%src_dist)) deallocate(Interp%src_dist)
if(allocated(Interp%num_found)) deallocate(Interp%num_found)
if(allocated(Interp%i_lon)) deallocate(Interp%i_lon)
if(allocated(Interp%j_lat)) deallocate(Interp%j_lat)

end subroutine horiz_interp_spherical_del

Expand Down
95 changes: 48 additions & 47 deletions horiz_interp/horiz_interp_type.F90
Original file line number Diff line number Diff line change
Expand Up @@ -53,26 +53,26 @@ module horiz_interp_type_mod
!<PUBLICTYPE >
!> @ingroup horiz_interp_type_mod
type horiz_interp_type
real, dimension(:,:), pointer :: faci =>NULL() !< weights for conservative scheme
real, dimension(:,:), pointer :: facj =>NULL() !< weights for conservative scheme
integer, dimension(:,:), pointer :: ilon =>NULL() !< indices for conservative scheme
integer, dimension(:,:), pointer :: jlat =>NULL() !< indices for conservative scheme
real, dimension(:,:), pointer :: area_src =>NULL() !< area of the source grid
real, dimension(:,:), pointer :: area_dst =>NULL() !< area of the destination grid
real, dimension(:,:,:), pointer :: wti =>NULL() !< weights for bilinear interpolation
real, dimension(:,:), allocatable :: faci !< weights for conservative scheme
real, dimension(:,:), allocatable :: facj !< weights for conservative scheme
integer, dimension(:,:), allocatable :: ilon !< indices for conservative scheme
integer, dimension(:,:), allocatable :: jlat !< indices for conservative scheme
real, dimension(:,:), allocatable :: area_src !< area of the source grid
real, dimension(:,:), allocatable :: area_dst !< area of the destination grid
real, dimension(:,:,:), allocatable :: wti !< weights for bilinear interpolation
!! wti ist used for derivative "weights" in bicubic
real, dimension(:,:,:), pointer :: wtj =>NULL() !< weights for bilinear interpolation
real, dimension(:,:,:), allocatable :: wtj !< weights for bilinear interpolation
!! wti ist used for derivative "weights" in bicubic
integer, dimension(:,:,:), pointer :: i_lon =>NULL() !< indices for bilinear interpolation
integer, dimension(:,:,:), allocatable :: i_lon !< indices for bilinear interpolation
!! and spherical regrid
integer, dimension(:,:,:), pointer :: j_lat =>NULL() !< indices for bilinear interpolation
integer, dimension(:,:,:), allocatable :: j_lat !< indices for bilinear interpolation
!! and spherical regrid
real, dimension(:,:,:), pointer :: src_dist =>NULL() !< distance between destination grid and
real, dimension(:,:,:), allocatable :: src_dist !< distance between destination grid and
!! neighbor source grid.
logical, dimension(:,:), pointer :: found_neighbors =>NULL() !< indicate whether destination grid
logical, dimension(:,:), allocatable :: found_neighbors !< indicate whether destination grid
!! has some source grid around it.
real :: max_src_dist
integer, dimension(:,:), pointer :: num_found => NULL()
integer, dimension(:,:), allocatable :: num_found
integer :: nlon_src !< size of source grid
integer :: nlat_src !< size of source grid
integer :: nlon_dst !< size of destination grid
Expand All @@ -82,26 +82,26 @@ module horiz_interp_type_mod
!! =2, bilinear interpolation
!! =3, spherical regrid
!! =4, bicubic regrid
real, dimension(:,:), pointer :: rat_x =>NULL() !< the ratio of coordinates of the dest grid
real, dimension(:,:), allocatable :: rat_x !< the ratio of coordinates of the dest grid
!! (x_dest -x_src_r)/(x_src_l -x_src_r)
!! and (y_dest -y_src_r)/(y_src_l -y_src_r)
real, dimension(:,:), pointer :: rat_y =>NULL() !< the ratio of coordinates of the dest grid
real, dimension(:,:), allocatable :: rat_y !< the ratio of coordinates of the dest grid
!! (x_dest -x_src_r)/(x_src_l -x_src_r)
!! and (y_dest -y_src_r)/(y_src_l -y_src_r)
real, dimension(:), pointer :: lon_in =>NULL() !< the coordinates of the source grid
real, dimension(:), pointer :: lat_in =>NULL() !< the coordinates of the source grid
real, dimension(:), allocatable :: lon_in !< the coordinates of the source grid
real, dimension(:), allocatable :: lat_in !< the coordinates of the source grid
logical :: I_am_initialized=.false.
integer :: version !< indicate conservative
!! interpolation version with value 1 or 2
!--- The following are for conservative interpolation scheme version 2 ( through xgrid)
integer :: nxgrid !< number of exchange grid
!! between src and dst grid.
integer, dimension(:), pointer :: i_src=>NULL() !< indices in source grid.
integer, dimension(:), pointer :: j_src=>NULL() !< indices in source grid.
integer, dimension(:), pointer :: i_dst=>NULL() !< indices in destination grid.
integer, dimension(:), pointer :: j_dst=>NULL() !< indices in destination grid.
real, dimension(:), pointer :: area_frac_dst=>NULL() !< area fraction in destination grid.
real, dimension(:,:), pointer :: mask_in=>NULL()
integer, dimension(:), allocatable :: i_src !< indices in source grid.
integer, dimension(:), allocatable :: j_src !< indices in source grid.
integer, dimension(:), allocatable :: i_dst !< indices in destination grid.
integer, dimension(:), allocatable :: j_dst !< indices in destination grid.
real, dimension(:), allocatable :: area_frac_dst !< area fraction in destination grid.
real, dimension(:,:), allocatable :: mask_in
end type
!</PUBLICTYPE>

Expand Down Expand Up @@ -181,43 +181,44 @@ subroutine stats ( dat, low, high, avg, miss, missing_value, mask )
end subroutine stats

!######################################################################################################################
!> @brief horiz_interp_type_eq creates a copy of the horiz_interp_type object
subroutine horiz_interp_type_eq(horiz_interp_out, horiz_interp_in)
type(horiz_interp_type), intent(inout) :: horiz_interp_out
type(horiz_interp_type), intent(in) :: horiz_interp_in
type(horiz_interp_type), intent(inout) :: horiz_interp_out !< Output object being set
type(horiz_interp_type), intent(in) :: horiz_interp_in !< Input object being copied

if(.not.horiz_interp_in%I_am_initialized) then
call mpp_error(FATAL,'horiz_interp_type_eq: horiz_interp_type variable on right hand side is unassigned')
endif

horiz_interp_out%faci => horiz_interp_in%faci
horiz_interp_out%facj => horiz_interp_in%facj
horiz_interp_out%ilon => horiz_interp_in%ilon
horiz_interp_out%jlat => horiz_interp_in%jlat
horiz_interp_out%area_src => horiz_interp_in%area_src
horiz_interp_out%area_dst => horiz_interp_in%area_dst
horiz_interp_out%wti => horiz_interp_in%wti
horiz_interp_out%wtj => horiz_interp_in%wtj
horiz_interp_out%i_lon => horiz_interp_in%i_lon
horiz_interp_out%j_lat => horiz_interp_in%j_lat
horiz_interp_out%src_dist => horiz_interp_in%src_dist
horiz_interp_out%found_neighbors => horiz_interp_in%found_neighbors
horiz_interp_out%faci = horiz_interp_in%faci
horiz_interp_out%facj = horiz_interp_in%facj
horiz_interp_out%ilon = horiz_interp_in%ilon
horiz_interp_out%jlat = horiz_interp_in%jlat
horiz_interp_out%area_src = horiz_interp_in%area_src
horiz_interp_out%area_dst = horiz_interp_in%area_dst
horiz_interp_out%wti = horiz_interp_in%wti
horiz_interp_out%wtj = horiz_interp_in%wtj
horiz_interp_out%i_lon = horiz_interp_in%i_lon
horiz_interp_out%j_lat = horiz_interp_in%j_lat
horiz_interp_out%src_dist = horiz_interp_in%src_dist
if (allocated(horiz_interp_in%found_neighbors)) horiz_interp_out%found_neighbors = horiz_interp_in%found_neighbors
horiz_interp_out%max_src_dist = horiz_interp_in%max_src_dist
horiz_interp_out%num_found => horiz_interp_in%num_found
horiz_interp_out%num_found = horiz_interp_in%num_found
horiz_interp_out%nlon_src = horiz_interp_in%nlon_src
horiz_interp_out%nlat_src = horiz_interp_in%nlat_src
horiz_interp_out%nlon_dst = horiz_interp_in%nlon_dst
horiz_interp_out%nlat_dst = horiz_interp_in%nlat_dst
horiz_interp_out%interp_method = horiz_interp_in%interp_method
horiz_interp_out%rat_x => horiz_interp_in%rat_x
horiz_interp_out%rat_y => horiz_interp_in%rat_y
horiz_interp_out%lon_in => horiz_interp_in%lon_in
horiz_interp_out%lat_in => horiz_interp_in%lat_in
horiz_interp_out%rat_x = horiz_interp_in%rat_x
horiz_interp_out%rat_y = horiz_interp_in%rat_y
horiz_interp_out%lon_in = horiz_interp_in%lon_in
horiz_interp_out%lat_in = horiz_interp_in%lat_in
horiz_interp_out%I_am_initialized = .true.
horiz_interp_out%i_src => horiz_interp_in%i_src
horiz_interp_out%j_src => horiz_interp_in%j_src
horiz_interp_out%i_dst => horiz_interp_in%i_dst
horiz_interp_out%j_dst => horiz_interp_in%j_dst
horiz_interp_out%area_frac_dst => horiz_interp_in%area_frac_dst
horiz_interp_out%i_src = horiz_interp_in%i_src
horiz_interp_out%j_src = horiz_interp_in%j_src
horiz_interp_out%i_dst = horiz_interp_in%i_dst
horiz_interp_out%j_dst = horiz_interp_in%j_dst
horiz_interp_out%area_frac_dst = horiz_interp_in%area_frac_dst
if(horiz_interp_in%interp_method == CONSERVE) then
horiz_interp_out%version = horiz_interp_in%version
if(horiz_interp_in%version==2) horiz_interp_out%nxgrid = horiz_interp_in%nxgrid
Expand Down
Loading