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

Disable printing for cells with negative vectors #82

Merged
merged 4 commits into from
Sep 17, 2023
Merged
Show file tree
Hide file tree
Changes from all 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
2 changes: 1 addition & 1 deletion Project.toml
Original file line number Diff line number Diff line change
@@ -1,7 +1,7 @@
name = "AtomsBase"
uuid = "a963bdd2-2df7-4f54-a1ee-49d51e6be12a"
authors = ["JuliaMolSim community"]
version = "0.3.4"
version = "0.3.5"

[deps]
LinearAlgebra = "37e2e46d-f89d-539d-b4ee-838fcccc9c8e"
Expand Down
2 changes: 0 additions & 2 deletions docs/make.jl
Original file line number Diff line number Diff line change
Expand Up @@ -6,7 +6,6 @@ DocMeta.setdocmeta!(AtomsBase, :DocTestSetup, :(using AtomsBase); recursive=true
makedocs(;
modules=[AtomsBase],
authors="JuliaMolSim community",
repo="https://github.com/JuliaMolSim/AtomsBase.jl/blob/{commit}{path}#{line}",
sitename="AtomsBase.jl",
format=Documenter.HTML(;
prettyurls=get(ENV, "CI", "false") == "true",
Expand All @@ -22,7 +21,6 @@ makedocs(;
"apireference.md"
],
checkdocs=:exports,
strict=true,
)

deploydocs(;
Expand Down
5 changes: 4 additions & 1 deletion src/visualize_ascii.jl
Original file line number Diff line number Diff line change
Expand Up @@ -16,11 +16,11 @@ function visualize_ascii(system::AbstractSystem{D}) where {D}
cell = austrip.(reduce(hcat, bounding_box(system)))
box = Vector(diag(cell))
shift = zero(box)
plot_box = D > 1

is_right_handed = det(cell) > 0
is_right_handed || return ""

plot_box = D > 1
is_orthorhombic = isdiag(cell)
if !is_orthorhombic
# Build an orthorhombic cell inscribing the actual unit cell
Expand All @@ -34,6 +34,9 @@ function visualize_ascii(system::AbstractSystem{D}) where {D}
plot_box = false
end

# If one of the box coordinates is negative
any(box .≤ 0) && return ""

# Normalise positions
normpos = [@. box * mod((shift + austrip(p)) / box, 1.0)
for p in position(system)]
Expand Down
32 changes: 22 additions & 10 deletions test/printing.jl
Original file line number Diff line number Diff line change
Expand Up @@ -24,15 +24,27 @@ using Test
end

@testset "Test ASCII representation of structures" begin
atoms = [:Si => [0.0, -0.125, 0.0],
:C => [0.125, 0.0, 0.0]]
box = [[10, 0.0, 0.0], [0.0, 5, 0.0], [0.0, 0.0, 7]]u"Å"
system = periodic_system(atoms, box; fractional=true)
println(visualize_ascii(system))
@testset "3D standard system" begin
atoms = [:Si => [0.0, -0.125, 0.0],
:C => [0.125, 0.0, 0.0]]
box = [[10, 0.0, 0.0], [0.0, 5, 0.0], [0.0, 0.0, 7]]u"Å"
system = periodic_system(atoms, box; fractional=true)
println(visualize_ascii(system))
end

@testset "2D standard system" begin
atoms = [:Si => [0.0, -0.125],
:C => [0.125, 0.0]]
box = [[10, 0.0], [0.0, 5]]u"Å"
system = periodic_system(atoms, box; fractional=true)
println(visualize_ascii(system))
end

atoms = [:Si => [0.0, -0.125],
:C => [0.125, 0.0]]
box = [[10, 0.0], [0.0, 5]]u"Å"
system = periodic_system(atoms, box; fractional=true)
println(visualize_ascii(system))
@testset "3D with negative unit cell" begin
atoms = [:Si => [0.75, 0.75, 0.75],
:Si => [0.0, 0.0, 0.0]]
box = [[-2.73, -2.73, 0.0], [-2.73, 0.0, -2.73], [0.0, -2.73, -2.73]]u"Å"
system = periodic_system(atoms, box; fractional=true)
println(visualize_ascii(system))
end
end