-
Notifications
You must be signed in to change notification settings - Fork 2
A. Basic usage: supercells, stripes and flakes, all out of the box!
Here we use only the honeytools
and honeyplots
high-level modules:
use honeytools ! We import only the top-level API here
use honeyplots ! And the bindings to matplotlib and gnuplot
use stdlib_strings, only: str => to_string ! Just for I/O
implicit none
type(xy_lattice) :: hlattice
type(unit_cell) :: mylayout
integer :: i
logical,allocatable :: NN(:,:),NNN(:,:)
Where in hlattice
we would store the generated x,y coordinates, their "A" or "B" sublattice-character and integer indices to build consistently the nearest and next-nearest neighbors matrices. Anyhow you probably don't care, just declare a type(xy_lattice
variable.
You would probably care instead about unit-cell initialization, which is done by passing a predefined orientation field (armchair
or zigzag
, they are constant expressions exported by the library, in most cases you don't want to customize your own orientation object, so we provide these two 🙃). You can instead change the lattice spacing and the origin, by passing the size
and origin
optional parameters (default values are 1
and [0,0]
).
Finally NN
and NNN
are two logical masks two identify which lattice points are to be connected by a nearest-neighbors or next-nearest-neighbors hoppings (or something else[1]). These are accepted by the plotting routines too, to allow easier visual exploration.
So, let's start with the examples!
mylayout = unit_cell(orientation=armchair)
do i = 1,5
hlattice = get_flake(radius=i,layout=mylayout)
if(i==1)then
call xy_print(hlattice)
print*
print*, "NN matrix:"
do j = 1,hlattice%size
write(*,*) (NN(j,k), k=1,hlattice%size)
enddo
call plot(hlattice,backend='gnuplot',set_terminal='dumb')
endif
call xy_nearest_neighbors(lattice=hlattice,nn_mask=NN)
call plot(hlattice,NN,figure_name='flake'//str(i)//'.svg')
enddo
Would print to standard output:
site # 1 [x,y]: -0.99999999999999989 0.0000000000000000
site # 2 [x,y]: -1.9999999999999998 2.2204460492503131E-016
site # 3 [x,y]: -2.5000000000000000 -0.86602540378443849
site # 4 [x,y]: -2.0000000000000004 -1.7320508075688770
site # 5 [x,y]: -1.0000000000000007 -1.7320508075688776
site # 6 [x,y]: -0.50000000000000000 -0.86602540378443882
site # 7 [x,y]: -0.99999999999999989 1.7320508075688772
site # 8 [x,y]: -1.9999999999999998 1.7320508075688774
site # 9 [x,y]: -2.5000000000000000 0.86602540378443871
site # 10 [x,y]: -0.50000000000000000 0.86602540378443837
site # 11 [x,y]: 0.50000000000000011 -0.86602540378443860
site # 12 [x,y]: -0.50000000000000044 -2.5980762113533156
site # 13 [x,y]: 0.49999999999999928 -2.5980762113533160
site # 14 [x,y]: 1.0000000000000000 -1.7320508075688774
site # 15 [x,y]: 0.50000000000000011 0.86602540378443860
site # 16 [x,y]: 1.0000000000000000 -2.4492935982947064E-016
site # 17 [x,y]: 0.50000000000000011 2.5980762113533160
site # 18 [x,y]: -0.49999999999999978 2.5980762113533160
site # 19 [x,y]: 1.0000000000000000 1.7320508075688770
site # 20 [x,y]: 2.0000000000000000 0.0000000000000000
site # 21 [x,y]: 1.9999999999999993 -1.7320508075688776
site # 22 [x,y]: 2.5000000000000000 -0.86602540378443882
site # 23 [x,y]: 2.0000000000000000 1.7320508075688772
site # 24 [x,y]: 2.5000000000000000 0.86602540378443837
NN matrix:
F T F F F T F F F T F F F F F F F F F F F F F F
T F T F F F F F T F F F F F F F F F F F F F F F
F T F T F F F F F F F F F F F F F F F F F F F F
F F T F T F F F F F F F F F F F F F F F F F F F
F F F T F T F F F F F T F F F F F F F F F F F F
T F F F T F F F F F T F F F F F F F F F F F F F
F F F F F F F T F T F F F F F F F T F F F F F F
F F F F F F T F T F F F F F F F F F F F F F F F
F T F F F F F T F F F F F F F F F F F F F F F F
T F F F F F T F F F F F F F T F F F F F F F F F
F F F F F T F F F F F F F T F T F F F F F F F F
F F F F T F F F F F F F T F F F F F F F F F F F
F F F F F F F F F F F T F T F F F F F F F F F F
F F F F F F F F F F T F T F F F F F F F T F F F
F F F F F F F F F T F F F F F T F F T F F F F F
F F F F F F F F F F T F F F T F F F F T F F F F
F F F F F F F F F F F F F F F F F T T F F F F F
F F F F F F T F F F F F F F F F T F F F F F F F
F F F F F F F F F F F F F F T F T F F F F F T F
F F F F F F F F F F F F F F F T F F F F F T F T
F F F F F F F F F F F F F T F F F F F F F T F F
F F F F F F F F F F F F F F F F F F F T T F F F
F F F F F F F F F F F F F F F F F F T F F F F T
F F F F F F F F F F F F F F F F F F F T F F T F
> Gnuplot GUI popping up...
3 +-------------------------------------+
| G G |
| |
2 |-+ |
| G G G G |
| |
1 |-+ |
| G G G G |
| |
0 |-+ G G G |
y | G |
| |
| G G G G |
-1 |-+ |
| |
| G G G G |
-2 |-+ |
| |
| + + G + G + + |
-3 +-------------------------------------+
-3 -2 -1 0 1 2 3
x
And save to disk:
radius=1 |
radius=2 |
radius=3 |
radius=4 |
radius=5 |
radius=6 |
mylayout = unit_cell(armchair,size=1,origin=[-5,-18])
hlattice = get_supercell(5,3,mylayout)
call xy_next_nearest_neighbors(hlattice,NNN,NN)
call plot(hlattice,NN,NNN,figure_name='supercell_armchair.svg')
mylayout = unit_cell(zigzag,size=3,origin=[-9,-9])
hlattice = get_supercell(3,5,mylayout)
call xy_next_nearest_neighbors(hlattice,NNN,NN)
call plot(hlattice,NN,NNN,figure_name='supercell_zigzag.svg')
Would draw also next-neighbor links!
supercell_armchair.svg |
supercell_zigzag.svg |
---|---|
mylayout = unit_cell(armchair)
hlattice = get_stripe(5,9,mylayout)
call xy_next_nearest_neighbors(hlattice,NNN,NN)
call plot(hlattice,NN,NNN,figure_name='stripe_armchair.svg')
mylayout = unit_cell(zigzag)
hlattice = get_stripe(7,5,mylayout)
call xy_next_nearest_neighbors(hlattice,NNN,NN)
call plot(hlattice,NN,NNN,figure_name='stripe_zigzag.svg')
stripe_armchair.svg |
stripe_zigzag.svg |
---|---|
do i = 2,5 ! triangles of size = 1 cannot exist
hlattice = get_triangle(i,unit_cell(zigzag))
call xy_next_nearest_neighbors(hlattice,NNN,NN)
call plot(hlattice,NN,NNN,figure_name='triangle'//str(i)//'.svg')
enddo
size=2 |
size=3 |
size=4 |
size=5 |
---|---|---|---|
[1]: We also provide facilities to build nth-order shells, just call the xy_shells
subroutine and build yourself the mask as shell_table == distance_set(n)
(these are the two intent(out)
variables in the call). Maybe I'll add an explicit example in the future.