Skip to content

Commit

Permalink
+add h to drifters interface (#408)
Browse files Browse the repository at this point in the history
This commit brings the drifters interface up-to-date with the current version of the drifters package, which requires h (layer thickness) to calculate the vertical movement of particles. The interfaces in the code and in config_src/external are updated to pass this information to the drifters package.
  • Loading branch information
cspencerjones authored Jul 25, 2023
1 parent 147ddf1 commit d5ba107
Show file tree
Hide file tree
Showing 2 changed files with 44 additions and 17 deletions.
44 changes: 31 additions & 13 deletions config_src/external/drifters/MOM_particles.F90
Original file line number Diff line number Diff line change
Expand Up @@ -11,28 +11,29 @@ module MOM_particles_mod
implicit none ; private

public particles, particles_run, particles_init, particles_save_restart, particles_end
public particles_to_k_space, particles_to_z_space

contains

!> Initializes particles container "parts"
subroutine particles_init(parts, Grid, Time, dt, u, v)
subroutine particles_init(parts, Grid, Time, dt, u, v, h)
! Arguments
type(particles), pointer, intent(out) :: parts !< Container for all types and memory
type(ocean_grid_type), target, intent(in) :: Grid !< Grid type from parent model
type(time_type), intent(in) :: Time !< Time type from parent model
real, intent(in) :: dt !< particle timestep [s]
real, dimension(:,:,:), intent(in) :: u !< Zonal velocity field [m s-1]
real, dimension(:,:,:), intent(in) :: v !< Meridional velocity field [m s-1]

real, intent(in) :: dt !< particle timestep in seconds [T ~> s]
real, dimension(:,:,:),intent(in) :: u !< Zonal velocity field [L T-1 ~> m s-1]
real, dimension(:,:,:),intent(in) :: v !< Meridional velocity field [L T-1 ~> m s-1]
real, dimension(:,:,:),intent(in) :: h !< Thickness of each layer [H ~> m or kg m-2]
end subroutine particles_init

!> The main driver the steps updates particles
subroutine particles_run(parts, time, uo, vo, ho, tv, stagger)
! Arguments
type(particles), pointer :: parts !< Container for all types and memory
type(time_type), intent(in) :: time !< Model time
real, dimension(:,:,:), intent(in) :: uo !< Ocean zonal velocity [m s-1]
real, dimension(:,:,:), intent(in) :: vo !< Ocean meridional velocity [m s-1]
real, dimension(:,:,:), intent(in) :: uo !< Ocean zonal velocity [L T-1 ~>m s-1]
real, dimension(:,:,:), intent(in) :: vo !< Ocean meridional velocity [L T-1~> m s-1]
real, dimension(:,:,:), intent(in) :: ho !< Ocean layer thickness [H ~> m or kg m-2]
type(thermo_var_ptrs), intent(in) :: tv !< structure containing pointers to available thermodynamic fields
integer, optional, intent(in) :: stagger !< Flag for whether velocities are staggered
Expand All @@ -41,21 +42,38 @@ end subroutine particles_run


!>Save particle locations (and sometimes other vars) to restart file
subroutine particles_save_restart(parts, temp, salt)
subroutine particles_save_restart(parts, h, temp, salt)
! Arguments
type(particles), pointer :: parts !< Container for all types and memory
real, dimension(:,:,:), optional, intent(in) :: temp !< Optional container for temperature
real, dimension(:,:,:), optional, intent(in) :: salt !< Optional container for salinity
real, dimension(:,:,:),intent(in) :: h !< Thickness of each layer [H ~> m or kg m-2]
real, dimension(:,:,:), optional, intent(in) :: temp !< Optional container for temperature [C ~> degC]
real, dimension(:,:,:), optional, intent(in) :: salt !< Optional container for salinity [S ~> ppt]

end subroutine particles_save_restart

!> Deallocate all memory and disassociated pointer
subroutine particles_end(parts, temp, salt)
subroutine particles_end(parts, h, temp, salt)
! Arguments
type(particles), pointer :: parts !< Container for all types and memory
real, dimension(:,:,:), optional, intent(in) :: temp !< Optional container for temperature
real, dimension(:,:,:), optional, intent(in) :: salt !< Optional container for salinity
real, dimension(:,:,:),intent(in) :: h !< Thickness of each layer [H ~> m or kg m-2]
real, dimension(:,:,:), optional, intent(in) :: temp !< Optional container for temperature [C ~> degC]
real, dimension(:,:,:), optional, intent(in) :: salt !< Optional container for salinity [S ~> ppt]

end subroutine particles_end

subroutine particles_to_k_space(parts, h)
! Arguments
type(particles), pointer :: parts !< Container for all types and memory
real, dimension(:,:,:),intent(in) :: h !< Thickness of layers [H ~> m or kg m-2]

end subroutine particles_to_k_space


subroutine particles_to_z_space(parts, h)
! Arguments
type(particles), pointer :: parts !< Container for all types and memory
real, dimension(:,:,:),intent(in) :: h !< Thickness of layers [H ~> m or kg m-2]

end subroutine particles_to_z_space

end module MOM_particles_mod
17 changes: 13 additions & 4 deletions src/core/MOM.F90
Original file line number Diff line number Diff line change
Expand Up @@ -162,7 +162,7 @@ module MOM
use MOM_offline_main, only : offline_advection_layer, offline_transport_end
use MOM_ice_shelf, only : ice_shelf_CS, ice_shelf_query, initialize_ice_shelf
use MOM_particles_mod, only : particles, particles_init, particles_run, particles_save_restart, particles_end

use MOM_particles_mod, only : particles_to_k_space, particles_to_z_space
implicit none ; private

#include <MOM_memory.h>
Expand Down Expand Up @@ -1541,6 +1541,10 @@ subroutine step_MOM_thermo(CS, G, GV, US, u, v, h, tv, fluxes, dtdia, &

call preAle_tracer_diagnostics(CS%tracer_Reg, G, GV)

if (CS%use_particles) then
call particles_to_z_space(CS%particles, h)
endif

if (CS%debug) then
call MOM_state_chksum("Pre-ALE ", u, v, h, CS%uh, CS%vh, G, GV, US, omit_corners=.true.)
call hchksum(tv%T,"Pre-ALE T", G%HI, haloshift=1, omit_corners=.true., scale=US%C_to_degC)
Expand Down Expand Up @@ -1588,6 +1592,11 @@ subroutine step_MOM_thermo(CS, G, GV, US, u, v, h, tv, fluxes, dtdia, &
call cpu_clock_end(id_clock_ALE)
endif ! endif for the block "if ( CS%use_ALE_algorithm )"


if (CS%use_particles) then
call particles_to_k_space(CS%particles, h)
endif

dynamics_stencil = min(3, G%Domain%nihalo, G%Domain%njhalo)
call create_group_pass(pass_uv_T_S_h, u, v, G%Domain, halo=dynamics_stencil)
if (associated(tv%T)) &
Expand Down Expand Up @@ -3232,7 +3241,7 @@ subroutine finish_MOM_initialization(Time, dirs, CS, restart_CSp)
G => CS%G ; GV => CS%GV ; US => CS%US

if (CS%use_particles) then
call particles_init(CS%particles, G, CS%Time, CS%dt_therm, CS%u, CS%v)
call particles_init(CS%particles, G, CS%Time, CS%dt_therm, CS%u, CS%v, CS%h)
endif

! Write initial conditions
Expand Down Expand Up @@ -3935,7 +3944,7 @@ subroutine save_MOM6_internal_state(CS, dirs, time, stamp_time)

! Could call save_restart(CS%restart_CSp) here

if (CS%use_particles) call particles_save_restart(CS%particles)
if (CS%use_particles) call particles_save_restart(CS%particles, CS%h)

end subroutine save_MOM6_internal_state

Expand Down Expand Up @@ -3978,7 +3987,7 @@ subroutine MOM_end(CS)
endif

if (CS%use_particles) then
call particles_end(CS%particles)
call particles_end(CS%particles, CS%h)
deallocate(CS%particles)
endif

Expand Down

0 comments on commit d5ba107

Please sign in to comment.