CliMA / Oceananigans.jl

🌊 Julia software for fast, friendly, flexible, ocean-flavored fluid dynamics on CPUs and GPUs
https://clima.github.io/OceananigansDocumentation/stable
MIT License
984 stars 193 forks source link

Illegal memory access on CPU and GPU with `LagrangianParticles` + `HydrostaticFreeSurfaceModel` + `LatitudeLongitudeGrid` #3852

Open ali-ramadhan opened 2 hours ago

ali-ramadhan commented 2 hours ago

I'm trying to add some particles to a hydrostatic model on a lat-lon grid, but ran into some CUDA memory issues. After reducing down to a MWE I noticed that it also segfaults on the CPU.

The MWE seems to be sensitive to the exact grid. Some lat-lon ranges lead to illegal memory accesses and others do not. I could not find a pattern though.

On the CPU the segfault seems to occur after ~2 iterations. On the GPU after ~29 iterations.

The particles are initialized within the domain and without any dynamics the particles should stay perfectly still. So I'm not sure where the illegal memory access is happening, but should be easy to debug on the CPU?

MWE:

using Oceananigans
using Oceananigans.Architectures: on_architecture

arch = GPU()

H = 100

Nλ = 100
Nφ = 200
Nz = 60

grid = LatitudeLongitudeGrid(
    arch;
    topology = (Bounded, Bounded, Bounded),
    size = (Nλ, Nφ, Nz),
    longitude = (0.79, 1.23),
    latitude = (-1.96, -1.12),
    z = (-H, 0),
    halo = (4, 4, 4)
)

Np = 100 # Number of particles

particles = LagrangianParticles(
    x = on_architecture(arch, 1 * ones(Np)),
    y = on_architecture(arch, -1.5 * ones(Np)),
    z = on_architecture(arch, -H/10 .* ones(Np))
)

model = HydrostaticFreeSurfaceModel(;
    grid,
    particles
)

for n in 1:100
    @info "Iteration $n..."
    time_step!(model, 0.1)
end

CPU segfault:

[ Info: Iteration 1...
[ Info: Iteration 2...

[503062] signal (11.1): Segmentation fault
in expression starting at /home/alir/atdepth/Oceananigans.jl/particles_error.jl:35
advect_particle at /home/alir/atdepth/Oceananigans.jl/src/Models/LagrangianParticleTracking/lagrangian_particle_advection.jl:0 [inlined]
macro expansion at /home/alir/atdepth/Oceananigans.jl/src/Models/LagrangianParticleTracking/lagrangian_particle_advection.jl:177 [inlined]
cpu__advect_particles! at /home/alir/.julia/packages/KernelAbstractions/491pi/src/macros.jl:291 [inlined]
cpu__advect_particles! at ./none:0
__thread_run at /home/alir/.julia/packages/KernelAbstractions/491pi/src/cpu.jl:144
unknown function (ip: 0x7c0090512182)
_jl_invoke at /cache/build/builder-amdci4-4/julialang/julia-release-1-dot-10/src/gf.c:2895 [inlined]
ijl_apply_generic at /cache/build/builder-amdci4-4/julialang/julia-release-1-dot-10/src/gf.c:3077
__run at /home/alir/.julia/packages/KernelAbstractions/491pi/src/cpu.jl:111
unknown function (ip: 0x7c009050feb3)
_jl_invoke at /cache/build/builder-amdci4-4/julialang/julia-release-1-dot-10/src/gf.c:2895 [inlined]
ijl_apply_generic at /cache/build/builder-amdci4-4/julialang/julia-release-1-dot-10/src/gf.c:3077
#_#16 at /home/alir/.julia/packages/KernelAbstractions/491pi/src/cpu.jl:46
_jl_invoke at /cache/build/builder-amdci4-4/julialang/julia-release-1-dot-10/src/gf.c:2895 [inlined]
ijl_apply_generic at /cache/build/builder-amdci4-4/julialang/julia-release-1-dot-10/src/gf.c:3077
jl_apply at /cache/build/builder-amdci4-4/julialang/julia-release-1-dot-10/src/julia.h:1982 [inlined]
do_apply at /cache/build/builder-amdci4-4/julialang/julia-release-1-dot-10/src/builtins.c:768
Kernel at /home/alir/.julia/packages/KernelAbstractions/491pi/src/cpu.jl:39
_jl_invoke at /cache/build/builder-amdci4-4/julialang/julia-release-1-dot-10/src/gf.c:2895 [inlined]
ijl_apply_generic at /cache/build/builder-amdci4-4/julialang/julia-release-1-dot-10/src/gf.c:3077
advect_lagrangian_particles! at /home/alir/atdepth/Oceananigans.jl/src/Models/LagrangianParticleTracking/lagrangian_particle_advection.jl:193
step_lagrangian_particles! at /home/alir/atdepth/Oceananigans.jl/src/Models/LagrangianParticleTracking/LagrangianParticleTracking.jl:143 [inlined]
step_lagrangian_particles! at /home/alir/atdepth/Oceananigans.jl/src/Models/HydrostaticFreeSurfaceModels/HydrostaticFreeSurfaceModels.jl:107 [inlined]
#time_step!#8 at /home/alir/atdepth/Oceananigans.jl/src/TimeSteppers/quasi_adams_bashforth_2.jl:124
time_step! at /home/alir/atdepth/Oceananigans.jl/src/TimeSteppers/quasi_adams_bashforth_2.jl:76
unknown function (ip: 0x7c00a0f12fbd)
_jl_invoke at /cache/build/builder-amdci4-4/julialang/julia-release-1-dot-10/src/gf.c:2895 [inlined]
ijl_apply_generic at /cache/build/builder-amdci4-4/julialang/julia-release-1-dot-10/src/gf.c:3077
top-level scope at /home/alir/atdepth/Oceananigans.jl/particles_error.jl:37
jl_toplevel_eval_flex at /cache/build/builder-amdci4-4/julialang/julia-release-1-dot-10/src/toplevel.c:925
jl_toplevel_eval_flex at /cache/build/builder-amdci4-4/julialang/julia-release-1-dot-10/src/toplevel.c:877
ijl_toplevel_eval_in at /cache/build/builder-amdci4-4/julialang/julia-release-1-dot-10/src/toplevel.c:985
eval at ./boot.jl:385 [inlined]
include_string at ./loading.jl:2076
_jl_invoke at /cache/build/builder-amdci4-4/julialang/julia-release-1-dot-10/src/gf.c:2895 [inlined]
ijl_apply_generic at /cache/build/builder-amdci4-4/julialang/julia-release-1-dot-10/src/gf.c:3077
_include at ./loading.jl:2136
include at ./client.jl:489
unknown function (ip: 0x7c00f54ff855)
_jl_invoke at /cache/build/builder-amdci4-4/julialang/julia-release-1-dot-10/src/gf.c:2895 [inlined]
ijl_apply_generic at /cache/build/builder-amdci4-4/julialang/julia-release-1-dot-10/src/gf.c:3077
jl_apply at /cache/build/builder-amdci4-4/julialang/julia-release-1-dot-10/src/julia.h:1982 [inlined]
do_call at /cache/build/builder-amdci4-4/julialang/julia-release-1-dot-10/src/interpreter.c:126
eval_value at /cache/build/builder-amdci4-4/julialang/julia-release-1-dot-10/src/interpreter.c:223
eval_stmt_value at /cache/build/builder-amdci4-4/julialang/julia-release-1-dot-10/src/interpreter.c:174 [inlined]
eval_body at /cache/build/builder-amdci4-4/julialang/julia-release-1-dot-10/src/interpreter.c:617
jl_interpret_toplevel_thunk at /cache/build/builder-amdci4-4/julialang/julia-release-1-dot-10/src/interpreter.c:775
jl_toplevel_eval_flex at /cache/build/builder-amdci4-4/julialang/julia-release-1-dot-10/src/toplevel.c:934
jl_toplevel_eval_flex at /cache/build/builder-amdci4-4/julialang/julia-release-1-dot-10/src/toplevel.c:877
jl_toplevel_eval_flex at /cache/build/builder-amdci4-4/julialang/julia-release-1-dot-10/src/toplevel.c:877
jl_toplevel_eval_flex at /cache/build/builder-amdci4-4/julialang/julia-release-1-dot-10/src/toplevel.c:877
ijl_toplevel_eval_in at /cache/build/builder-amdci4-4/julialang/julia-release-1-dot-10/src/toplevel.c:985
eval at ./boot.jl:385 [inlined]
eval_user_input at /cache/build/builder-amdci4-4/julialang/julia-release-1-dot-10/usr/share/julia/stdlib/v1.10/REPL/src/REPL.jl:150
repl_backend_loop at /cache/build/builder-amdci4-4/julialang/julia-release-1-dot-10/usr/share/julia/stdlib/v1.10/REPL/src/REPL.jl:246
#start_repl_backend#46 at /cache/build/builder-amdci4-4/julialang/julia-release-1-dot-10/usr/share/julia/stdlib/v1.10/REPL/src/REPL.jl:231
start_repl_backend at /cache/build/builder-amdci4-4/julialang/julia-release-1-dot-10/usr/share/julia/stdlib/v1.10/REPL/src/REPL.jl:228
_jl_invoke at /cache/build/builder-amdci4-4/julialang/julia-release-1-dot-10/src/gf.c:2895 [inlined]
ijl_apply_generic at /cache/build/builder-amdci4-4/julialang/julia-release-1-dot-10/src/gf.c:3077
#run_repl#59 at /cache/build/builder-amdci4-4/julialang/julia-release-1-dot-10/usr/share/julia/stdlib/v1.10/REPL/src/REPL.jl:389
run_repl at /cache/build/builder-amdci4-4/julialang/julia-release-1-dot-10/usr/share/julia/stdlib/v1.10/REPL/src/REPL.jl:375
jfptr_run_repl_91805.1 at /home/alir/.julia/juliaup/julia-1.10.5+0.x64.linux.gnu/lib/julia/sys.so (unknown line)
_jl_invoke at /cache/build/builder-amdci4-4/julialang/julia-release-1-dot-10/src/gf.c:2895 [inlined]
ijl_apply_generic at /cache/build/builder-amdci4-4/julialang/julia-release-1-dot-10/src/gf.c:3077
#1013 at ./client.jl:432
jfptr_YY.1013_82772.1 at /home/alir/.julia/juliaup/julia-1.10.5+0.x64.linux.gnu/lib/julia/sys.so (unknown line)
_jl_invoke at /cache/build/builder-amdci4-4/julialang/julia-release-1-dot-10/src/gf.c:2895 [inlined]
ijl_apply_generic at /cache/build/builder-amdci4-4/julialang/julia-release-1-dot-10/src/gf.c:3077
jl_apply at /cache/build/builder-amdci4-4/julialang/julia-release-1-dot-10/src/julia.h:1982 [inlined]
jl_f__call_latest at /cache/build/builder-amdci4-4/julialang/julia-release-1-dot-10/src/builtins.c:812
#invokelatest#2 at ./essentials.jl:892 [inlined]
invokelatest at ./essentials.jl:889 [inlined]
run_main_repl at ./client.jl:416
exec_options at ./client.jl:333
_start at ./client.jl:552
jfptr__start_82798.1 at /home/alir/.julia/juliaup/julia-1.10.5+0.x64.linux.gnu/lib/julia/sys.so (unknown line)
_jl_invoke at /cache/build/builder-amdci4-4/julialang/julia-release-1-dot-10/src/gf.c:2895 [inlined]
ijl_apply_generic at /cache/build/builder-amdci4-4/julialang/julia-release-1-dot-10/src/gf.c:3077
jl_apply at /cache/build/builder-amdci4-4/julialang/julia-release-1-dot-10/src/julia.h:1982 [inlined]
true_main at /cache/build/builder-amdci4-4/julialang/julia-release-1-dot-10/src/jlapi.c:582
jl_repl_entrypoint at /cache/build/builder-amdci4-4/julialang/julia-release-1-dot-10/src/jlapi.c:731
main at /cache/build/builder-amdci4-4/julialang/julia-release-1-dot-10/cli/loader_exe.c:58
unknown function (ip: 0x7c00f758ce07)
__libc_start_main at /usr/lib/libc.so.6 (unknown line)
unknown function (ip: 0x4010b8)
Allocations: 67298744 (Pool: 67235612; Big: 63132); GC: 66
fish: Job 1, 'julia --project' terminated by signal SIGSEGV (Address boundary error)

GPU illegal memory access:

[ Info: Skipping precompilation since __precompile__(false). Importing Oceananigans [9e8cae18-63c1-5223-a75c-80ca9d6e9a09].
[ Info: Iteration 1...
[ Info: Iteration 2...
[ Info: Iteration 3...
[ Info: Iteration 4...
[ Info: Iteration 5...
[ Info: Iteration 6...
[ Info: Iteration 7...
[ Info: Iteration 8...
[ Info: Iteration 9...
[ Info: Iteration 10...
[ Info: Iteration 11...
[ Info: Iteration 12...
[ Info: Iteration 13...
[ Info: Iteration 14...
[ Info: Iteration 15...
[ Info: Iteration 16...
[ Info: Iteration 17...
[ Info: Iteration 18...
[ Info: Iteration 19...
[ Info: Iteration 20...
[ Info: Iteration 21...
[ Info: Iteration 22...
[ Info: Iteration 23...
[ Info: Iteration 24...
[ Info: Iteration 25...
[ Info: Iteration 26...
[ Info: Iteration 27...
[ Info: Iteration 28...
[ Info: Iteration 29...
ERROR: LoadError: CUDA error: an illegal memory access was encountered (code 700, ERROR_ILLEGAL_ADDRESS)
Stacktrace:
  [1] throw_api_error(res::CUDA.cudaError_enum)
    @ CUDA ~/.julia/packages/CUDA/2kjXI/lib/cudadrv/libcuda.jl:30
  [2] check
    @ ~/.julia/packages/CUDA/2kjXI/lib/cudadrv/libcuda.jl:37 [inlined]
  [3] cuStreamGetCaptureInfo
    @ ~/.julia/packages/CUDA/2kjXI/lib/utils/call.jl:34 [inlined]
  [4] capture_status(stream::CUDA.CuStream)
    @ CUDA ~/.julia/packages/CUDA/2kjXI/lib/cudadrv/graph.jl:174
  [5] is_capturing
    @ ~/.julia/packages/CUDA/2kjXI/lib/cudadrv/graph.jl:179 [inlined]
  [6] convert(::Type{CUDA.CuPtr{Float64}}, managed::CUDA.Managed{CUDA.DeviceMemory})
    @ CUDA ~/.julia/packages/CUDA/2kjXI/src/memory.jl:539
  [7] unsafe_convert
    @ ~/.julia/packages/CUDA/2kjXI/src/array.jl:434 [inlined]
  [8] #pointer#1123
    @ ~/.julia/packages/CUDA/2kjXI/src/array.jl:392 [inlined]
  [9] pointer (repeats 2 times)
    @ ~/.julia/packages/CUDA/2kjXI/src/array.jl:384 [inlined]
 [10] unsafe_convert(::Type{CUDA.CuDeviceArray{Float64, 3, 1}}, a::CUDA.CuArray{Float64, 3, CUDA.DeviceMemory})
    @ CUDA ~/.julia/packages/CUDA/2kjXI/src/array.jl:454
 [11] adapt_storage
    @ ~/.julia/packages/CUDA/2kjXI/src/compiler/execution.jl:162 [inlined]
 [12] adapt_structure
    @ ~/.julia/packages/Adapt/7T9au/src/Adapt.jl:57 [inlined]
 [13] adapt
    @ ~/.julia/packages/Adapt/7T9au/src/Adapt.jl:40 [inlined]
 [14] cudaconvert
    @ ~/.julia/packages/CUDA/2kjXI/src/compiler/execution.jl:198 [inlined]
 [15] map
    @ ./tuple.jl:293 [inlined]
 [16] macro expansion
    @ ~/.julia/packages/CUDA/2kjXI/src/compiler/execution.jl:110 [inlined]
 [17] #launch_heuristic#1200
    @ ~/.julia/packages/CUDA/2kjXI/src/gpuarrays.jl:17 [inlined]
 [18] launch_heuristic
    @ ~/.julia/packages/CUDA/2kjXI/src/gpuarrays.jl:15 [inlined]
 [19] gpu_call(::GPUArrays.var"#6#7", ::CUDA.CuArray{Float64, 3, CUDA.DeviceMemory}, ::Float64; target::CUDA.CuArray{Float64, 3, CUDA.DeviceMemory}, elements::Nothing, threads::Nothing, blocks::Nothing, name::Nothing)
    @ GPUArrays ~/.julia/packages/GPUArrays/qt4ax/src/device/execution.jl:61
 [20] gpu_call
    @ ~/.julia/packages/GPUArrays/qt4ax/src/device/execution.jl:34 [inlined]
 [21] fill!
    @ ~/.julia/packages/GPUArrays/qt4ax/src/host/construction.jl:14 [inlined]
 [22] fill!
    @ ~/atdepth/Oceananigans.jl/src/Fields/field.jl:407 [inlined]
 [23] initialize_free_surface_state!(state::Oceananigans.Models.HydrostaticFreeSurfaceModels.SplitExplicitState{…}, η::Field{…}, timestepper::Oceananigans.Models.HydrostaticFreeSurfaceModels.ForwardBackwardScheme)
    @ Oceananigans.Models.HydrostaticFreeSurfaceModels ~/atdepth/Oceananigans.jl/src/Models/HydrostaticFreeSurfaceModels/split_explicit_free_surface_kernels.jl:243
 [24] #apply_regionally!#56
    @ ~/atdepth/Oceananigans.jl/src/Utils/multi_region_transformation.jl:121 [inlined]
 [25] apply_regionally!
    @ ~/atdepth/Oceananigans.jl/src/Utils/multi_region_transformation.jl:118 [inlined]
 [26] macro expansion
    @ ~/atdepth/Oceananigans.jl/src/Utils/multi_region_transformation.jl:233 [inlined]
 [27] macro expansion
    @ ~/atdepth/Oceananigans.jl/src/Models/HydrostaticFreeSurfaceModels/split_explicit_free_surface_kernels.jl:325 [inlined]
 [28] macro expansion
    @ ~/atdepth/Oceananigans.jl/src/Utils/multi_region_transformation.jl:250 [inlined]
 [29] split_explicit_free_surface_step!(free_surface::SplitExplicitFreeSurface{…}, model::HydrostaticFreeSurfaceModel{…}, Δt::Float64, χ::Float64)
    @ Oceananigans.Models.HydrostaticFreeSurfaceModels ~/atdepth/Oceananigans.jl/src/Models/HydrostaticFreeSurfaceModels/split_explicit_free_surface_kernels.jl:324
 [30] ab2_step_free_surface!
    @ ~/atdepth/Oceananigans.jl/src/Models/HydrostaticFreeSurfaceModels/split_explicit_free_surface_kernels.jl:295 [inlined]
 [31] ab2_step!
    @ ~/atdepth/Oceananigans.jl/src/Models/HydrostaticFreeSurfaceModels/hydrostatic_free_surface_ab2_step.jl:23 [inlined]
 [32] time_step!(model::HydrostaticFreeSurfaceModel{…}, Δt::Float64; callbacks::Vector{…}, euler::Bool)
    @ Oceananigans.TimeSteppers ~/atdepth/Oceananigans.jl/src/TimeSteppers/quasi_adams_bashforth_2.jl:114
 [33] time_step!(model::HydrostaticFreeSurfaceModel{…}, Δt::Float64)
    @ Oceananigans.TimeSteppers ~/atdepth/Oceananigans.jl/src/TimeSteppers/quasi_adams_bashforth_2.jl:76
 [34] top-level scope
    @ ~/atdepth/Oceananigans.jl/particles_error.jl:37
 [35] include(fname::String)
    @ Base.MainInclude ./client.jl:489
 [36] top-level scope
    @ REPL[1]:1
in expression starting at /home/alir/atdepth/Oceananigans.jl/particles_error.jl:35
Some type information was truncated. Use `show(err)` to see complete types.
ali-ramadhan commented 2 hours ago

Running with --check-bounds=yes on the CPU provides a strong hint:

at index [39914881, -59303136, 54]

Yeah that'll do it lol.


[ Info: Iteration 1...
[ Info: Iteration 2...
ERROR: LoadError: BoundsError: attempt to access 109×208×68 OffsetArray(::Array{Float64, 3}, -3:105, -3:204, -3:64) with eltype Float64 with indices -3:105×-3:204×-3:64 at index [39914881, -59303136, 54]
Stacktrace:
  [1] throw_boundserror(A::OffsetArrays.OffsetArray{Float64, 3, Array{Float64, 3}}, I::Tuple{Int64, Int64, Int64})
    @ Base ./abstractarray.jl:737
  [2] checkbounds
    @ ./abstractarray.jl:702 [inlined]
  [3] getindex
    @ ~/.julia/packages/OffsetArrays/hwmnB/src/OffsetArrays.jl:422 [inlined]
  [4] getindex
    @ ~/atdepth/Oceananigans.jl/src/Fields/field.jl:401 [inlined]
  [5] _interpolate
    @ ~/atdepth/Oceananigans.jl/src/Fields/interpolate.jl:295 [inlined]
  [6] interpolate
    @ ~/atdepth/Oceananigans.jl/src/Fields/interpolate.jl:245 [inlined]
  [7] advect_particle
    @ ~/atdepth/Oceananigans.jl/src/Models/LagrangianParticleTracking/lagrangian_particle_advection.jl:113 [inlined]
  [8] macro expansion
    @ ~/atdepth/Oceananigans.jl/src/Models/LagrangianParticleTracking/lagrangian_particle_advection.jl:177 [inlined]
  [9] cpu__advect_particles!
    @ ~/.julia/packages/KernelAbstractions/491pi/src/macros.jl:291 [inlined]
 [10] cpu__advect_particles!(__ctx__::KernelAbstractions.CompilerMetadata{…}, particles::StructArrays.StructVector{…}, restitution::Float64, grid::LatitudeLongitudeGrid{…}, Δt::Float64, velocities::@NamedTuple{…})
    @ Oceananigans.Models.LagrangianParticleTracking ./none:0
 [11] __thread_run(tid::Int64, len::Int64, rem::Int64, obj::KernelAbstractions.Kernel{…}, ndrange::Nothing, iterspace::KernelAbstractions.NDIteration.NDRange{…}, args::Tuple{…}, dynamic::KernelAbstractions.NDIteration.NoDynamicCheck)
    @ KernelAbstractions ~/.julia/packages/KernelAbstractions/491pi/src/cpu.jl:144
 [12] __run(obj::KernelAbstractions.Kernel{…}, ndrange::Nothing, iterspace::KernelAbstractions.NDIteration.NDRange{…}, args::Tuple{…}, dynamic::KernelAbstractions.NDIteration.NoDynamicCheck, static_threads::Bool)
    @ KernelAbstractions ~/.julia/packages/KernelAbstractions/491pi/src/cpu.jl:111
 [13] (::KernelAbstractions.Kernel{…})(::StructArrays.StructVector{…}, ::Vararg{…}; ndrange::Nothing, workgroupsize::Nothing)
    @ KernelAbstractions ~/.julia/packages/KernelAbstractions/491pi/src/cpu.jl:46
 [14] (::KernelAbstractions.Kernel{…})(::StructArrays.StructVector{…}, ::Vararg{…})
    @ KernelAbstractions ~/.julia/packages/KernelAbstractions/491pi/src/cpu.jl:39
 [15] advect_lagrangian_particles!(particles::LagrangianParticles{…}, model::HydrostaticFreeSurfaceModel{…}, Δt::Float64)
    @ Oceananigans.Models.LagrangianParticleTracking ~/atdepth/Oceananigans.jl/src/Models/LagrangianParticleTracking/lagrangian_particle_advection.jl:193
 [16] step_lagrangian_particles!
    @ ~/atdepth/Oceananigans.jl/src/Models/LagrangianParticleTracking/LagrangianParticleTracking.jl:143 [inlined]
 [17] step_lagrangian_particles!
    @ ~/atdepth/Oceananigans.jl/src/Models/HydrostaticFreeSurfaceModels/HydrostaticFreeSurfaceModels.jl:107 [inlined]
 [18] time_step!(model::HydrostaticFreeSurfaceModel{…}, Δt::Float64; callbacks::Vector{…}, euler::Bool)
    @ Oceananigans.TimeSteppers ~/atdepth/Oceananigans.jl/src/TimeSteppers/quasi_adams_bashforth_2.jl:124
 [19] time_step!(model::HydrostaticFreeSurfaceModel{…}, Δt::Float64)
    @ Oceananigans.TimeSteppers ~/atdepth/Oceananigans.jl/src/TimeSteppers/quasi_adams_bashforth_2.jl:76
 [20] top-level scope
    @ ~/atdepth/Oceananigans.jl/particles_error.jl:37
 [21] include(fname::String)
    @ Base.MainInclude ./client.jl:489
 [22] top-level scope
    @ REPL[1]:1
in expression starting at /home/alir/atdepth/Oceananigans.jl/particles_error.jl:35
Some type information was truncated. Use `show(err)` to see complete types.
jagoosw commented 52 minutes ago

~I have had errors like this and found it came from the particle leaving the domain and then the boundary is enforcer moving it back in but it had moved so far out when it went back "in" it was out of the domain the other way~

I just realised there is no velocity in the MWE so this can't be what's going on.

ali-ramadhan commented 4 minutes ago

I think these lines should be using ξnode, ηnode, and rnode:

https://github.com/CliMA/Oceananigans.jl/blob/7cbf013cb6bed2bd7cef0f4d8e5f04c078e50ee0/src/Models/LagrangianParticleTracking/lagrangian_particle_advection.jl#L136-L142

I'll open a PR with a fix tomorrow. Should probably also add a test for particle advection on a lat-lon grid.


Some debug printing inside advect_particle with 1 particle:

[ Info: Iteration 1...
[ Info: X=(1.0, -1.5, -10.0), I=(47, 109, 53)
[ Info: (before) X⁺=(1.0, -1.5, -10.0)
(iᴿ, jᴿ, kᴿ) = (101, 201, 61)
(xᴸ, yᴸ, zᴸ) = (87813.63270401207, -217942.05622333512, -100.0)
(xᴿ, yᴿ, zᴿ) = (136722.49142523398, -124538.3178419058, 0.0)
(x⁺, y⁺, z⁺) = (175626.26540802413, -249075.1356838116, -10.0)
[ Info: (after) X⁺=(175626.26540802413, -249075.1356838116, -10.0)
[ Info: Iteration 2...
[ Info: X=(175626.26540802413, -249075.1356838116, -10.0), I=(39914880, -59303137, 53)
ERROR: LoadError: BoundsError: attempt to access 109×208×68 OffsetArray(::Array{Float64, 3}, -3:105, -3:204, -3:64) with eltype Float64 with indices -3:105×-3:204×-3:64 at index [39914881, -59303136, 54]