Skip to content
This repository was archived by the owner on Aug 28, 2024. It is now read-only.

Commit

Permalink
Fixed a bunch of issues, but the NetCDF file is still not quite there
Browse files Browse the repository at this point in the history
  • Loading branch information
daminton committed Dec 3, 2022
1 parent 17e7dcb commit 52a595a
Show file tree
Hide file tree
Showing 2 changed files with 26 additions and 21 deletions.
45 changes: 25 additions & 20 deletions src/encounter/encounter_io.f90
Original file line number Diff line number Diff line change
Expand Up @@ -55,6 +55,7 @@ module subroutine encounter_io_initialize_output(self, param)
real(SP) :: sfill
logical :: fileExists
character(len=STRMAX) :: errmsg
integer(I4B), dimension(2), parameter :: collider_dimension = [1,2]

dfill = ieee_value(dfill, IEEE_QUIET_NAN)
sfill = ieee_value(sfill, IEEE_QUIET_NAN)
Expand Down Expand Up @@ -82,7 +83,7 @@ module subroutine encounter_io_initialize_output(self, param)

call check( nf90_def_var(self%ncid, TIME_DIMNAME, self%out_type, self%time_dimid, self%time_varid), "encounter_io_initialize_output nf90_def_var time_varid" )
call check( nf90_def_var(self%ncid, ENCID_DIMNAME, NF90_INT, self%encid_dimid, self%encid_varid), "encounter_io_initialize_output nf90_def_var encid_varid" )
call check( nf90_def_var(self%ncid, COLLIDER_DIMNAME, NF90_INT, self%collider_dimid, self%encid_varid), "encounter_io_initialize_output nf90_def_var collider_varid" )
call check( nf90_def_var(self%ncid, COLLIDER_DIMNAME, NF90_INT, self%collider_dimid, self%collider_varid), "encounter_io_initialize_output nf90_def_var collider_varid" )
call check( nf90_def_var(self%ncid, NENC_VARNAME, NF90_INT, self%time_dimid, self%nenc_varid), "encounter_io_initialize_output nf90_def_var nenc_varid" )
call check( nf90_def_var(self%ncid, NAME_VARNAME, NF90_CHAR, [self%str_dimid, self%collider_dimid, self%encid_dimid], self%name_varid), "encounter_io_initialize_output nf90_def_var name_varid" )
call check( nf90_def_var(self%ncid, ID_DIMNAME, NF90_INT, [self%collider_dimid, self%encid_dimid, self%time_dimid], self%id_varid), "encounter_io_initialize_output nf90_def_var id_varid" )
Expand All @@ -99,6 +100,7 @@ module subroutine encounter_io_initialize_output(self, param)

! Take the file out of define mode
call check( nf90_enddef(self%ncid), "encounter_io_initialize_output nf90_enddef" )
call check( nf90_put_var(self%ncid, self%collider_varid, collider_dimension, start=[1], count=[2]), "encounter_io_initialize_output nf90_put_var collider_varid" )

return

Expand Down Expand Up @@ -168,25 +170,28 @@ module subroutine encounter_io_write_frame(self, iu, param)
i = iu%ienc_frame
n = int(self%nenc, kind=I4B)
call check( nf90_set_fill(iu%ncid, nf90_nofill, old_mode), "encounter_io_write_frame_base nf90_set_fill" )
call check( nf90_put_var(iu%ncid, iu%time_varid, self%t, start=[i]), "netcdf_write_hdr_system nf90_put_var time_varid" )
call check( nf90_put_var(iu%ncid, iu%xhx_varid, self%x1(1, :), start=[1, 1, i], count=[1,n,1]), "netcdf_write_frame_base nf90_put_var xhx_varid 1" )
call check( nf90_put_var(iu%ncid, iu%xhy_varid, self%x1(2, :), start=[1, 1, i], count=[1,n,1]), "netcdf_write_frame_base nf90_put_var xhy_varid 1" )
call check( nf90_put_var(iu%ncid, iu%xhz_varid, self%x1(3, :), start=[1, 1, i], count=[1,n,1]), "netcdf_write_frame_base nf90_put_var xhz_varid 1" )
call check( nf90_put_var(iu%ncid, iu%xhx_varid, self%x2(1, :), start=[2, 1, i], count=[1,n,1]), "netcdf_write_frame_base nf90_put_var xhx_varid 2" )
call check( nf90_put_var(iu%ncid, iu%xhy_varid, self%x2(2, :), start=[2, 1, i], count=[1,n,1]), "netcdf_write_frame_base nf90_put_var xhy_varid 2" )
call check( nf90_put_var(iu%ncid, iu%xhz_varid, self%x2(3, :), start=[2, 1, i], count=[1,n,1]), "netcdf_write_frame_base nf90_put_var xhz_varid 2" )
call check( nf90_put_var(iu%ncid, iu%vhx_varid, self%v1(1, :), start=[1, 1, i], count=[1,n,1]), "netcdf_write_frame_base nf90_put_var vhx_varid 1" )
call check( nf90_put_var(iu%ncid, iu%vhy_varid, self%v1(2, :), start=[1, 1, i], count=[1,n,1]), "netcdf_write_frame_base nf90_put_var vhy_varid 1" )
call check( nf90_put_var(iu%ncid, iu%vhz_varid, self%v1(3, :), start=[1, 1, i], count=[1,n,1]), "netcdf_write_frame_base nf90_put_var vhz_varid 1" )
call check( nf90_put_var(iu%ncid, iu%vhx_varid, self%v2(1, :), start=[2, 1, i], count=[1,n,1]), "netcdf_write_frame_base nf90_put_var vhx_varid 2" )
call check( nf90_put_var(iu%ncid, iu%vhy_varid, self%v2(2, :), start=[2, 1, i], count=[1,n,1]), "netcdf_write_frame_base nf90_put_var vhy_varid 2" )
call check( nf90_put_var(iu%ncid, iu%vhz_varid, self%v2(3, :), start=[2, 1, i], count=[1,n,1]), "netcdf_write_frame_base nf90_put_var vhz_varid 2" )
call check( nf90_put_var(iu%ncid, iu%name_varid, self%name1(:), start=[1, 1, i], count=[NAMELEN,1,1]), "netcdf_write_frame_base nf90_put_var name 1" )
call check( nf90_put_var(iu%ncid, iu%name_varid, self%name2(:), start=[1, 2, i], count=[NAMELEN,1,1]), "netcdf_write_frame_base nf90_put_var name 2" )
call check( nf90_put_var(iu%ncid, iu%Gmass_varid, self%Gmass1(:), start=[1, 1, i], count=[1,n,1]), "netcdf_write_frame_base nf90_put_var Gmass 1" )
call check( nf90_put_var(iu%ncid, iu%Gmass_varid, self%Gmass2(:), start=[2, 1, i], count=[1,n,1]), "netcdf_write_frame_base nf90_put_var Gmass 2" )
call check( nf90_put_var(iu%ncid, iu%radius_varid, self%radius1(:), start=[1, 1, i], count=[1,n,1]), "netcdf_write_frame_base nf90_put_var radius 1" )
call check( nf90_put_var(iu%ncid, iu%radius_varid, self%radius2(:), start=[2, 1, i], count=[1,n,1]), "netcdf_write_frame_base nf90_put_var radius 2" )
call check( nf90_put_var(iu%ncid, iu%time_varid, self%t, start=[i]), "encounter_io_write_frame nf90_put_var time_varid" )
call check( nf90_put_var(iu%ncid, iu%encid_varid, [(i, i=1,n)], start=[i], count=[n]), "encounter_io_write_frame nf90_put_var encid_varid" )

call check( nf90_put_var(iu%ncid, iu%nenc_varid, self%nenc, start=[i]), "encounter_io_frame nf90_put_var nenc_varid" )
call check( nf90_put_var(iu%ncid, iu%name_varid, self%name1(1:n), start=[1, 1, i], count=[NAMELEN,1,1]), "netcdf_write_frame_base nf90_put_var name 1" )
call check( nf90_put_var(iu%ncid, iu%name_varid, self%name2(1:n), start=[1, 2, i], count=[NAMELEN,1,1]), "netcdf_write_frame_base nf90_put_var name 2" )
call check( nf90_put_var(iu%ncid, iu%xhx_varid, self%x1(1, 1:n), start=[1, 1, i], count=[1,n,1]), "netcdf_write_frame_base nf90_put_var xhx_varid 1" )
call check( nf90_put_var(iu%ncid, iu%xhy_varid, self%x1(2, 1:n), start=[1, 1, i], count=[1,n,1]), "netcdf_write_frame_base nf90_put_var xhy_varid 1" )
call check( nf90_put_var(iu%ncid, iu%xhz_varid, self%x1(3, 1:n), start=[1, 1, i], count=[1,n,1]), "netcdf_write_frame_base nf90_put_var xhz_varid 1" )
call check( nf90_put_var(iu%ncid, iu%xhx_varid, self%x2(1, 1:n), start=[2, 1, i], count=[1,n,1]), "netcdf_write_frame_base nf90_put_var xhx_varid 2" )
call check( nf90_put_var(iu%ncid, iu%xhy_varid, self%x2(2, 1:n), start=[2, 1, i], count=[1,n,1]), "netcdf_write_frame_base nf90_put_var xhy_varid 2" )
call check( nf90_put_var(iu%ncid, iu%xhz_varid, self%x2(3, 1:n), start=[2, 1, i], count=[1,n,1]), "netcdf_write_frame_base nf90_put_var xhz_varid 2" )
call check( nf90_put_var(iu%ncid, iu%vhx_varid, self%v1(1, 1:n), start=[1, 1, i], count=[1,n,1]), "netcdf_write_frame_base nf90_put_var vhx_varid 1" )
call check( nf90_put_var(iu%ncid, iu%vhy_varid, self%v1(2, 1:n), start=[1, 1, i], count=[1,n,1]), "netcdf_write_frame_base nf90_put_var vhy_varid 1" )
call check( nf90_put_var(iu%ncid, iu%vhz_varid, self%v1(3, 1:n), start=[1, 1, i], count=[1,n,1]), "netcdf_write_frame_base nf90_put_var vhz_varid 1" )
call check( nf90_put_var(iu%ncid, iu%vhx_varid, self%v2(1, 1:n), start=[2, 1, i], count=[1,n,1]), "netcdf_write_frame_base nf90_put_var vhx_varid 2" )
call check( nf90_put_var(iu%ncid, iu%vhy_varid, self%v2(2, 1:n), start=[2, 1, i], count=[1,n,1]), "netcdf_write_frame_base nf90_put_var vhy_varid 2" )
call check( nf90_put_var(iu%ncid, iu%vhz_varid, self%v2(3, 1:n), start=[2, 1, i], count=[1,n,1]), "netcdf_write_frame_base nf90_put_var vhz_varid 2" )
call check( nf90_put_var(iu%ncid, iu%Gmass_varid, self%Gmass1(1:n), start=[1, 1, i], count=[1,n,1]), "netcdf_write_frame_base nf90_put_var Gmass 1" )
call check( nf90_put_var(iu%ncid, iu%Gmass_varid, self%Gmass2(1:n), start=[2, 1, i], count=[1,n,1]), "netcdf_write_frame_base nf90_put_var Gmass 2" )
call check( nf90_put_var(iu%ncid, iu%radius_varid, self%radius1(1:n), start=[1, 1, i], count=[1,n,1]), "netcdf_write_frame_base nf90_put_var radius 1" )
call check( nf90_put_var(iu%ncid, iu%radius_varid, self%radius2(1:n), start=[2, 1, i], count=[1,n,1]), "netcdf_write_frame_base nf90_put_var radius 2" )


return
Expand Down
2 changes: 1 addition & 1 deletion src/symba/symba_util.f90
Original file line number Diff line number Diff line change
Expand Up @@ -925,7 +925,7 @@ module subroutine symba_util_resize_storage(self, nnew)
end if

if (nnew > nold) then
allocate(encounter_storage(2 * nnew) :: tmp)
allocate(encounter_storage(8 * nnew) :: tmp)
if (lmalloc) then
do i = 1, nold
if (allocated(self%encounter_history%frame(i)%item)) tmp%frame(i) = self%encounter_history%frame(i)%item
Expand Down

0 comments on commit 52a595a

Please sign in to comment.