symba_drift.f90 Source File


This file depends on

sourcefile~~symba_drift.f90~~EfferentGraph sourcefile~symba_drift.f90 symba_drift.f90 sourcefile~symba_module.f90 symba_module.f90 sourcefile~symba_drift.f90->sourcefile~symba_module.f90 sourcefile~swiftest_module.f90 swiftest_module.f90 sourcefile~symba_drift.f90->sourcefile~swiftest_module.f90 sourcefile~symba_module.f90->sourcefile~swiftest_module.f90 sourcefile~helio_module.f90 helio_module.f90 sourcefile~symba_module.f90->sourcefile~helio_module.f90 sourcefile~operator_module.f90 operator_module.f90 sourcefile~swiftest_module.f90->sourcefile~operator_module.f90 sourcefile~base_module.f90 base_module.f90 sourcefile~swiftest_module.f90->sourcefile~base_module.f90 sourcefile~collision_module.f90 collision_module.f90 sourcefile~swiftest_module.f90->sourcefile~collision_module.f90 sourcefile~lambda_function_module.f90 lambda_function_module.f90 sourcefile~swiftest_module.f90->sourcefile~lambda_function_module.f90 sourcefile~encounter_module.f90 encounter_module.f90 sourcefile~swiftest_module.f90->sourcefile~encounter_module.f90 sourcefile~walltime_module.f90 walltime_module.f90 sourcefile~swiftest_module.f90->sourcefile~walltime_module.f90 sourcefile~io_progress_bar_module.f90 io_progress_bar_module.f90 sourcefile~swiftest_module.f90->sourcefile~io_progress_bar_module.f90 sourcefile~netcdf_io_module.f90 netcdf_io_module.f90 sourcefile~swiftest_module.f90->sourcefile~netcdf_io_module.f90 sourcefile~solver_module.f90 solver_module.f90 sourcefile~swiftest_module.f90->sourcefile~solver_module.f90 sourcefile~helio_module.f90->sourcefile~swiftest_module.f90 sourcefile~whm_module.f90 whm_module.f90 sourcefile~helio_module.f90->sourcefile~whm_module.f90 sourcefile~coarray_module.f90 coarray_module.f90 sourcefile~base_module.f90->sourcefile~coarray_module.f90 sourcefile~collision_module.f90->sourcefile~base_module.f90 sourcefile~collision_module.f90->sourcefile~encounter_module.f90 sourcefile~encounter_module.f90->sourcefile~base_module.f90 sourcefile~encounter_module.f90->sourcefile~netcdf_io_module.f90 sourcefile~walltime_module.f90->sourcefile~base_module.f90 sourcefile~io_progress_bar_module.f90->sourcefile~base_module.f90 sourcefile~netcdf_io_module.f90->sourcefile~base_module.f90 sourcefile~solver_module.f90->sourcefile~base_module.f90 sourcefile~solver_module.f90->sourcefile~lambda_function_module.f90 sourcefile~whm_module.f90->sourcefile~swiftest_module.f90

Contents

Source Code


Source Code

! Copyright 2024 - The Minton Group at Purdue University
! This file is part of Swiftest.
! Swiftest is free software: you can redistribute it and/or modify it under the terms of the GNU General Public License 
! as published by the Free Software Foundation, either version 3 of the License, or (at your option) any later version.
! Swiftest is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty 
! of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.
! You should have received a copy of the GNU General Public License along with Swiftest. 
! If not, see: https://www.gnu.org/licenses. 

  submodule (symba) s_symba_drift
   use swiftest
contains

   module subroutine symba_drift_pl(self, nbody_system, param, dt)
      !! author: David A. Minton
      !!
      !! Wrapper function used to call the body drift routine from a symba_pl structure
      implicit none
      ! Arguments
      class(symba_pl),              intent(inout) :: self   
         !! Helio massive body object
      class(swiftest_nbody_system), intent(inout) :: nbody_system 
         !! Swiftest nbody system object
      class(swiftest_parameters),   intent(in)    :: param  
         !! Current run configuration parameters 
      real(DP),                     intent(in)    :: dt     
         !! Stepsize

      if (self%nbody == 0) return
      associate(pl => self, npl => self%nbody)
         select type(nbody_system)
         class is (symba_nbody_system)
            pl%lmask(1:npl) = pl%status(1:npl) /= INACTIVE .and. pl%levelg(1:npl) == nbody_system%irec
            call helio_drift_body(pl, nbody_system, param, dt)
            pl%lmask(1:npl) = pl%status(1:npl) /= INACTIVE 
         end select
      end associate

      return
   end subroutine symba_drift_pl


   module subroutine symba_drift_tp(self, nbody_system, param, dt)
      !! author: David A. Minton
      !!
      !! Wrapper function used to call the body drift routine from a symba_pl structure
      implicit none
      ! Arguments
      class(symba_tp),              intent(inout) :: self   !! Helio massive body object
      class(swiftest_nbody_system), intent(inout) :: nbody_system !! Swiftest nbody system object
      class(swiftest_parameters),   intent(in)    :: param  !! Current run configuration parameters 
      real(DP),                     intent(in)    :: dt     !! Stepsize

      if (self%nbody == 0) return
      associate (tp => self, ntp => self%nbody)
         select type(nbody_system)
         class is (symba_nbody_system)
            tp%lmask(1:ntp) = tp%status(1:ntp) /= INACTIVE .and. tp%levelg(1:ntp) == nbody_system%irec
            call helio_drift_body(tp, nbody_system, param, dt)
            tp%lmask(1:ntp) = tp%status(1:ntp) /= INACTIVE 
         end select
      end associate

      return
   end subroutine symba_drift_tp

end submodule s_symba_drift