Skip to content

Commit

Permalink
Add plain-text and html-style graphical representation of structures …
Browse files Browse the repository at this point in the history
…to AtomsBase (#70)
  • Loading branch information
mfherbst authored Apr 8, 2023
1 parent e1cabff commit 47aca33
Show file tree
Hide file tree
Showing 6 changed files with 148 additions and 4 deletions.
6 changes: 4 additions & 2 deletions Project.toml
Original file line number Diff line number Diff line change
@@ -1,11 +1,13 @@
name = "AtomsBase"
uuid = "a963bdd2-2df7-4f54-a1ee-49d51e6be12a"
authors = ["JuliaMolSim community"]
version = "0.3.2"
version = "0.3.3"

[deps]
LinearAlgebra = "37e2e46d-f89d-539d-b4ee-838fcccc9c8e"
PeriodicTable = "7b2266bf-644c-5ea3-82d8-af4bbd25a884"
Printf = "de0858da-6303-5e67-8744-51eddeeeb8d7"
Requires = "ae029012-a4dd-5104-9daa-d747884805df"
StaticArrays = "90137ffa-7385-5640-81b9-e52037218182"
Unitful = "1986cc42-f94f-5a68-af5c-568840ba703d"
UnitfulAtomic = "a7773ee8-282e-5fa2-be4e-bd808c38a91a"
Expand All @@ -18,8 +20,8 @@ UnitfulAtomic = "1"
julia = "1.6"

[extras]
Test = "8dfed614-e22c-5e08-85e1-65c5234f0b40"
Pkg = "44cfe95a-1eb2-52ea-b672-e2afdf69b78f"
Test = "8dfed614-e22c-5e08-85e1-65c5234f0b40"
Unitful = "1986cc42-f94f-5a68-af5c-568840ba703d"
UnitfulAtomic = "a7773ee8-282e-5fa2-be4e-bd808c38a91a"

Expand Down
1 change: 1 addition & 0 deletions docs/src/apireference.md
Original file line number Diff line number Diff line change
Expand Up @@ -21,6 +21,7 @@ periodicity
species_type
atomkeys
hasatomkey
visualize_ascii
```

## Species / atom properties
Expand Down
10 changes: 10 additions & 0 deletions src/AtomsBase.jl
Original file line number Diff line number Diff line change
Expand Up @@ -2,13 +2,23 @@ module AtomsBase
using Unitful
using UnitfulAtomic
using StaticArrays
using Requires

include("interface.jl")
include("properties.jl")
include("visualize_ascii.jl")
include("show.jl")
include("flexible_system.jl")
include("atomview.jl")
include("atom.jl")
include("fast_system.jl")

function __init__()
@require AtomsView="ee286e10-dd2d-4ff2-afcb-0a3cd50c8041" begin
function Base.show(io::IO, mime::MIME"text/html", system::FlexibleSystem)
write(io, AtomsView.visualize_structure(system, mime))
end
end
end

end
9 changes: 7 additions & 2 deletions src/show.jl
Original file line number Diff line number Diff line change
Expand Up @@ -56,13 +56,18 @@ function show_system(io::IO, ::MIME"text/plain", system::AbstractSystem{D}) wher
extra_line = true
@printf io " %-17s : %s\n" string(k) string(v)
end

# TODO Better would be some ascii-graphical representation of the structure
if length(system) < 10
extra_line && println(io)
for atom in system
println(io, " ", atom)
end
extra_line = true
end

ascii = visualize_ascii(system)
if !isempty(ascii)
extra_line && println(io)
println(io, " ", replace(ascii, "\n" => "\n "))
end
end

Expand Down
112 changes: 112 additions & 0 deletions src/visualize_ascii.jl
Original file line number Diff line number Diff line change
@@ -0,0 +1,112 @@
using UnitfulAtomic
using LinearAlgebra

export visualize_ascii

"""
Build an ASCII representation of the passed atomistic structure. The string may
be empty if the passed structure could not be represented
(structure not supported or invalid).
"""
function visualize_ascii(system::AbstractSystem{D}) where {D}
# Heavily inspired by the ascii art plot algorithm of GPAW
# See output.py in the GPAW sources

# Unit cell matrix (vectors column-by-column) and plotting box (xyz)
cell = austrip.(reduce(hcat, bounding_box(system)))
box = Vector(diag(cell))
shift = zero(box)

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
# by lumping each cartesian component on the diagonal
box = sum.(eachrow(cell))

# Shift centre of the original unit cell to the centre of the orthorhomic cell
centre_atoms = austrip.(sum(position(system)) / length(system))
shift = box / 2 - centre_atoms

plot_box = false
end

# Normalise positions
normpos = [@. box * mod((shift + austrip(p)) / box, 1.0)
for p in position(system)]

scaling = 1.3
sx = nothing
sy = nothing
sz = nothing
canvas = nothing
while scaling > 0.1
if D == 2
scaled = scaling .* [box[1], 0.0, box[2]]
else
scaled = scaling .* box .* (1.0, 0.25, 0.5)
end

sx, sy, sz = round.(Int, scaled)
canvas = fill(' ', sx + sy + 4, sy + sz + 1)
all(size(canvas) .≤ 100) && break
scaling *= 0.9
end

if D == 2
projector = Diagonal([sx sz] ./ box)
elseif D == 3
projector = [sx sy 0; 0 sy sz] * Diagonal(1 ./ box)
end
pos2d = [1 .+ round.(Int, projector * p .+ eps(Float64)) for p in normpos]

# Draw box onto canvas
if plot_box
# 7 Corners:
canvas[2 + sy, 1 + sy ] = '.'
canvas[2 + sx + sy, 1 + sy ] = '.'
canvas[2 + sx + sy, 1 + sy + sz] = '.'
canvas[2 + sy, 1 + sy + sz] = '.'
canvas[2, 1 ] = '*'
canvas[2 + sx, 1 ] = '*'
canvas[2, 1 + sz] = '*'
if D < 3 # Better use a *
canvas[2 + sx + sy, 1 + sy + sz] = '*'
end

for y in 1:sy-1 # Bars along y
canvas[2 + y, 1 + y ] = '/'
canvas[2 + y + sx, 1 + y ] = '/'
canvas[2 + y, 1 + y + sz] = '/'
end

for z in 1:sz-1
canvas[2, 1 + z] = '|'
canvas[2 + sy, 1 + sy + z] = '|'
canvas[2 + sx + sy, 1 + sy + z] = '|'
end

for x in 1:sx-1 # Bars along x
canvas[2 + x, 1 ] = '-'
canvas[2 + x + sy, 1 + sy ] = '-'
canvas[2 + x + sy, 1 + sy + sz] = '-'
end
end

depth2d = Inf * ones(size(canvas)) # Keep track of things covering each other
for (iatom, symbol) in enumerate(atomic_symbol(system))
x, y = pos2d[iatom]
for (i, c) in enumerate(string(symbol))
if normpos[iatom][2] < depth2d[x + i, y]
canvas[x + i, y] = c
depth2d[x + i, y] = normpos[iatom][2]
end
end
end

join(reverse([join(col) for col in eachcol(canvas)]), "\n")
end
visualize_ascii(::AbstractSystem{1}) = ""
14 changes: 14 additions & 0 deletions test/printing.jl
Original file line number Diff line number Diff line change
Expand Up @@ -22,3 +22,17 @@ using Test
show(stdout, MIME("text/plain"), fast_system)
show(stdout, MIME("text/plain"), fast_system[1])
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))

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

2 comments on commit 47aca33

@mfherbst
Copy link
Member Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

@JuliaRegistrator
Copy link

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Registration pull request created: JuliaRegistries/General/81263

After the above pull request is merged, it is recommended that a tag is created on this repository for the registered package version.

This will be done automatically if the Julia TagBot GitHub Action is installed, or can be done manually through the github interface, or via:

git tag -a v0.3.3 -m "<description of version>" 47aca331c24f1b75f39cf499e2d0070af588708f
git push origin v0.3.3

Please sign in to comment.