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

Commit

Permalink
Merge branch 'debug' into IntelAdvisor
Browse files Browse the repository at this point in the history
  • Loading branch information
daminton committed Oct 26, 2021
2 parents 81d12c8 + 17d6c57 commit a5e38ca
Show file tree
Hide file tree
Showing 4 changed files with 7 additions and 9 deletions.
7 changes: 1 addition & 6 deletions src/encounter/encounter_check.f90
Original file line number Diff line number Diff line change
Expand Up @@ -189,7 +189,7 @@ module subroutine encounter_check_all_plplm(param, nplm, nplt, xplm, vplm, xplt,
call move_alloc(itmp, index1)
allocate(itmp(nenc+plmplt_nenc))
itmp(1:nenc) = index2(1:nenc)
itmp(nenc+1:nenc+plmplt_nenc) = plmplt_index2(1:plmplt_nenc)
itmp(nenc+1:nenc+plmplt_nenc) = plmplt_index2(1:plmplt_nenc) + nplm ! Be sure to shift these indices back to their natural range
call move_alloc(itmp, index2)
allocate(ltmp(nenc+plmplt_nenc))
ltmp(1:nenc) = lvdotr(1:nenc)
Expand Down Expand Up @@ -543,11 +543,6 @@ subroutine encounter_check_all_sort_and_sweep_plplm(nplm, nplt, xplm, vplm, xplt
call timer%stop()
write(*,*) "plplm check : ",timer%count_stop_step - timer%count_start_step

! Shift the tiny body indices back to their natural range
index2(:) = index2(:) + nplm

call timer%reset()
call timer%start()
call encounter_check_reduce_broadphase(ntot, nenc, index1, index2, lencounter, lvdotr)
call timer%stop()
write(*,*) "plplm reduce: ",timer%count_stop_step - timer%count_start_step
Expand Down
3 changes: 1 addition & 2 deletions src/netcdf/netcdf.f90
Original file line number Diff line number Diff line change
Expand Up @@ -68,7 +68,7 @@ module function netcdf_get_old_t_final_system(self, param) result(old_t_final)
real(DP), dimension(1) :: val
real(DP) :: KE_orb_orig, KE_spin_orig, PE_orig, Ltmp

call param%nciu%open(param, readonly=.true.)
call param%nciu%open(param)
call check( nf90_inquire_dimension(param%nciu%ncid, param%nciu%time_dimid, len=itmax) )
call check( nf90_inquire_dimension(param%nciu%ncid, param%nciu%id_dimid, len=idmax) )
allocate(vals(idmax))
Expand Down Expand Up @@ -110,7 +110,6 @@ module function netcdf_get_old_t_final_system(self, param) result(old_t_final)
end if

deallocate(vals)
call param%nciu%close()

return
end function netcdf_get_old_t_final_system
Expand Down
5 changes: 5 additions & 0 deletions src/symba/symba_encounter_check.f90
Original file line number Diff line number Diff line change
Expand Up @@ -35,6 +35,8 @@ module function symba_encounter_check_pl(self, param, system, dt, irec) result(l
nplm = pl%nplm
nplt = npl - nplm

call pl%set_renc(irec)

if (nplt == 0) then
call encounter_check_all_plpl(param, npl, pl%xh, pl%vh, pl%renc, dt, lvdotr, index1, index2, nenc)
else
Expand Down Expand Up @@ -116,6 +118,8 @@ module function symba_encounter_check(self, param, system, dt, irec) result(lany
nenc_enc = count(lencmask(:))
if (nenc_enc == 0) return

call pl%set_renc(irec)

allocate(encidx(nenc_enc))
allocate(lencounter(nenc_enc))
encidx(:) = pack([(k, k = 1, self%nenc)], lencmask(:))
Expand Down Expand Up @@ -203,6 +207,7 @@ module function symba_encounter_check_tp(self, param, system, dt, irec) result(l
if (self%nbody == 0) return

associate(tp => self, ntp => self%nbody, pl => system%pl, npl => system%pl%nbody)
call pl%set_renc(irec)
call encounter_check_all_pltp(param, npl, ntp, pl%xh, pl%vh, tp%xh, tp%vh, pl%renc, dt, lvdotr, index1, index2, nenc)

lany_encounter = nenc > 0
Expand Down
1 change: 0 additions & 1 deletion src/symba/symba_step.f90
Original file line number Diff line number Diff line change
Expand Up @@ -173,7 +173,6 @@ module recursive subroutine symba_step_recur_system(self, param, t, ireci)
nloops = NTENC
end if
do j = 1, nloops
call system%pl%set_renc(irecp)
lencounter = plplenc_list%encounter_check(param, system, dtl, irecp) .or. pltpenc_list%encounter_check(param, system, dtl, irecp)

call plplenc_list%kick(system, dth, irecp, 1)
Expand Down

0 comments on commit a5e38ca

Please sign in to comment.