Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Make AtomView rely on struct-of-array format #68

Merged
merged 5 commits into from
Mar 22, 2023
Merged
Show file tree
Hide file tree
Changes from 4 commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
1 change: 1 addition & 0 deletions docs/src/apireference.md
Original file line number Diff line number Diff line change
Expand Up @@ -38,6 +38,7 @@ element

```@docs
Atom
AtomView
FlexibleSystem
AbstractSystem
atomic_system
Expand Down
44 changes: 37 additions & 7 deletions src/atomview.jl
Original file line number Diff line number Diff line change
Expand Up @@ -2,15 +2,45 @@
# A simple view datastructure for atoms of struct of array systems
#
export AtomView
struct AtomView{FS <: AbstractSystem}
system::FS

"""
AtomView{S<:AbstractSystem}

Species type for atoms of systems implemented as struct-of-arrays.
Can be queried with the same API than for other species, like [`Atom`](@ref).

See [FastSystem](@ref Struct-of-Arrays-/-FastSystem) for an example of system
using `AtomView` as its species type.

## Example
```jldoctest; setup=:(using AtomsBase, Unitful; atoms = Atom[:C => [0.25, 0.25, 0.25]u"Å", :C => [0.75, 0.75, 0.75]u"Å"]; box = [[1, 0, 0], [0, 1, 0], [0, 0, 1]]u"Å"; boundary_conditions = [Periodic(), Periodic(), DirichletZero()])
julia> system = FastSystem(atoms, box, boundary_conditions);

julia> atom = system[2]
AtomView(C, atomic_number = 6, atomic_mass = 12.011 u):
position : [0.75,0.75,0.75]u"Å"

julia> atom isa AtomView{typeof(system)}
true

julia> atomic_symbol(atom)
:C
```
"""
struct AtomView{S<:AbstractSystem}
system::S
index::Int
end
velocity(v::AtomView) = velocity(v.system, v.index)
position(v::AtomView) = position(v.system, v.index)
atomic_mass(v::AtomView) = atomic_mass(v.system, v.index)
atomic_symbol(v::AtomView) = atomic_symbol(v.system, v.index)
atomic_number(v::AtomView) = atomic_number(v.system, v.index)

function velocity(v::AtomView)
vel = velocity(v.system)
ismissing(vel) && return missing
return vel[v.index]
end
position(v::AtomView) = position(v.system)[v.index]
atomic_mass(v::AtomView) = atomic_mass(v.system)[v.index]
atomic_symbol(v::AtomView) = atomic_symbol(v.system)[v.index]
atomic_number(v::AtomView) = atomic_number(v.system)[v.index]
n_dimensions(v::AtomView) = n_dimensions(v.system)
element(atom::AtomView) = element(atomic_number(atom))

Expand Down
14 changes: 5 additions & 9 deletions src/fast_system.jl
Original file line number Diff line number Diff line change
Expand Up @@ -55,18 +55,14 @@ atomic_number(s::FastSystem) = s.atomic_number
atomic_mass(s::FastSystem) = s.atomic_mass
velocity(::FastSystem) = missing

position(s::FastSystem, i) = s.position[i]
mfherbst marked this conversation as resolved.
Show resolved Hide resolved
atomic_symbol(s::FastSystem, i) = s.atomic_symbol[i]
atomic_number(s::FastSystem, i) = s.atomic_number[i]
atomic_mass(s::FastSystem, i) = s.atomic_mass[i]
velocity(::FastSystem, i) = missing

# System property access
function Base.getindex(system::FastSystem, x::Symbol)
if x in (:bounding_box, :boundary_conditions)
mfherbst marked this conversation as resolved.
Show resolved Hide resolved
getfield(system, x)
if x === :bounding_box
bounding_box(system)
elseif x === :boundary_conditions
boundary_conditions(system)
else
throw(KeyError("Key $x not found"))
throw(KeyError(x))
end
end
Base.haskey(::FastSystem, x::Symbol) = x in (:bounding_box, :boundary_conditions)
Expand Down
6 changes: 4 additions & 2 deletions src/flexible_system.jl
Original file line number Diff line number Diff line change
Expand Up @@ -13,8 +13,10 @@ end

# System property access
function Base.getindex(system::FlexibleSystem, x::Symbol)
if x in (:bounding_box, :boundary_conditions)
getfield(system, x)
if x === :bounding_box
bounding_box(system)
elseif x === :boundary_conditions
boundary_conditions(system)
else
getindex(system.data, x)
end
Expand Down
9 changes: 9 additions & 0 deletions test/fast_system.jl
Original file line number Diff line number Diff line change
Expand Up @@ -2,6 +2,7 @@ using AtomsBase
using Test
using Unitful
using PeriodicTable
using StaticArrays

@testset "Fast system" begin
box = [[1, 0, 0], [0, 1, 0], [0, 0, 1]]u"m"
Expand All @@ -20,6 +21,8 @@ using PeriodicTable
@test !isinfinite(system)
@test element(system[1]) == element(:C)
@test keys(system) == (:bounding_box, :boundary_conditions)
@test haskey(system, :boundary_conditions)
@test system[:boundary_conditions][1] == Periodic()
@test atomkeys(system) == (:position, :atomic_symbol, :atomic_number, :atomic_mass)
@test keys(system[1]) == (:position, :atomic_symbol, :atomic_number, :atomic_mass)
@test hasatomkey(system, :atomic_symbol)
Expand All @@ -39,6 +42,12 @@ using PeriodicTable
:atomic_mass => atomic_mass(atoms[1]),
]

# check type stability
get_b_vector(syst) = bounding_box(syst)[2]
@test @inferred(get_b_vector(system)) === SVector{3}([0.0, 1.0, 0.0]u"m")
@test @inferred(position(system, 1)) === SVector{3}([0.25, 0.25, 0.25]u"m")
mfherbst marked this conversation as resolved.
Show resolved Hide resolved
@test ismissing(@inferred(velocity(system, 2)))

# Test AtomView
for method in (position, atomic_mass, atomic_symbol, atomic_number)
@test method(system[1]) == method(system, 1)
Expand Down
4 changes: 4 additions & 0 deletions test/interface.jl
Original file line number Diff line number Diff line change
Expand Up @@ -56,5 +56,9 @@ using PeriodicTable
@test ismissing(velocity(fast))
@test all(position(fast) .== position(flexible))
@test all(atomic_symbol(fast) .== atomic_symbol(flexible))

# type stability
get_z_periodicity(syst) = syst[:boundary_conditions][3]
@test @inferred(BoundaryCondition, get_z_periodicity(flexible)) === DirichletZero()
end
end
4 changes: 2 additions & 2 deletions test/printing.jl
Original file line number Diff line number Diff line change
Expand Up @@ -13,12 +13,12 @@ using Test

flexible_system = periodic_system(atoms, box; fractional=true, data=-12)
@test repr(flexible_system) == """
FlexibleSystem(CSi, periodic = TTT, bounding_box = [[10.0, 0.0, 0.0], [0.0, 5.0, 0.0], [0.0, 0.0, 7.0]]u"")"""
mfherbst marked this conversation as resolved.
Show resolved Hide resolved
FlexibleSystem(CSi, periodic = TTT, bounding_box = [[10.0, 0.0, 0.0], [0.0, 5.0, 0.0], [0.0, 0.0, 7.0]]u"Å")"""
show(stdout, MIME("text/plain"), flexible_system)

fast_system = FastSystem(flexible_system)
@test repr(fast_system) == """
FastSystem(CSi, periodic = TTT, bounding_box = [[10.0, 0.0, 0.0], [0.0, 5.0, 0.0], [0.0, 0.0, 7.0]]u"")"""
FastSystem(CSi, periodic = TTT, bounding_box = [[10.0, 0.0, 0.0], [0.0, 5.0, 0.0], [0.0, 0.0, 7.0]]u"Å")"""
show(stdout, MIME("text/plain"), fast_system)
show(stdout, MIME("text/plain"), fast_system[1])
end