1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241
|
program spglib_example
use, intrinsic :: iso_c_binding
implicit none
call run_examples()
contains
subroutine run_examples()
print '("Rutile:")'
call example_spg_get_dataset_rutile()
print '("")'
print '("FCC:")'
call example_spg_get_dataset_FCC()
print '("")'
print '("Wurtzite:")'
call example_spg_get_dataset_wurtzite()
end subroutine run_examples
subroutine example_spg_get_dataset_rutile()
integer(c_int) :: num_atom
real(c_double) :: symprec = 1e-6
real(c_double) :: lattice(3, 3)
real(c_double) :: positions(3, 6)
integer(c_int) :: atom_types(6)
integer(c_int) :: mesh(3) = [4, 4, 6]
integer(c_int) :: shifted(3)
integer(c_int) :: grid_point(3, 96)
integer(c_int) :: map(96)
integer(c_int) :: is_time_reversal = 1
real :: a
num_atom = 6
lattice(1, :) = [4.0, 0.0, 0.0]
lattice(2, :) = [0.0, 4.0, 0.0]
lattice(3, :) = [0.0, 0.0, 2.6]
positions(:, 1) = [0.0, 0.0, 0.0]
positions(:, 2) = [0.5, 0.5, 0.5]
positions(:, 3) = [0.3, 0.3, 0.0]
positions(:, 4) = [0.7, 0.7, 0.0]
positions(:, 5) = [0.2, 0.8, 0.5]
positions(:, 6) = [0.8, 0.2, 0.5]
atom_types(:) = [1, 1, 2, 2, 2, 2]
shifted(:) = [1, 1, 1]
call show_syminfo(num_atom, lattice, symprec, atom_types, positions)
call show_grid_info(num_atom, lattice, symprec, atom_types, positions, &
mesh, shifted, is_time_reversal)
end subroutine example_spg_get_dataset_rutile
subroutine example_spg_get_dataset_FCC()
integer(c_int) :: num_atom
real(c_double) :: symprec = 1e-6
real(c_double) :: lattice(3, 3)
real(c_double) :: positions(3, 1)
integer(c_int) :: atom_types(1)
integer(c_int) :: mesh(3) = [4, 4, 4]
integer(c_int) :: shifted(3)
integer(c_int) :: grid_point(3, 64)
integer(c_int) :: map(64)
integer(c_int) :: is_time_reversal = 1
num_atom = 1
lattice(1, :) = [0.0, 2.0, 2.0]
lattice(2, :) = [2.0, 0.0, 2.0]
lattice(3, :) = [2.0, 2.0, 0.0]
positions(:, 1) = [0.0, 0.0, 0.0]
atom_types(1) = 1
shifted = [0, 0, 0]
call show_syminfo(num_atom, lattice, symprec, atom_types, positions)
call show_grid_info(num_atom, lattice, symprec, atom_types, positions, &
mesh, shifted, is_time_reversal)
end subroutine example_spg_get_dataset_FCC
subroutine example_spg_get_dataset_wurtzite
integer(c_int) :: num_atom
real(c_double) :: symprec = 1e-6
real(c_double) :: lattice(3, 3)
real(c_double) :: positions(3, 4)
integer(c_int) :: atom_types(4)
integer(c_int) :: mesh(3) = [6, 6, 4]
integer(c_int) :: shifted(3)
integer(c_int) :: grid_point(3, 144)
integer(c_int) :: map(144)
integer(c_int) :: is_time_reversal = 1
num_atom = 4
lattice(:, 1) = [3.111, -1.5555, 0.0]
lattice(:, 2) = [0.0, 2.6942050311733885, 0.0]
lattice(:, 3) = [0.0, 0.0, 4.988]
positions(:, 1) = [1.0/3, 2.0/3, 0.0]
positions(:, 2) = [2.0/3, 1.0/3, 0.5]
positions(:, 3) = [1.0/3, 2.0/3, 0.6181]
positions(:, 4) = [2.0/3, 1.0/3, 0.1181]
atom_types(:) = [1, 1, 2, 2]
shifted(:) = [0, 0, 1]
call show_syminfo(num_atom, lattice, symprec, atom_types, positions)
call show_grid_info(num_atom, lattice, symprec, atom_types, positions, &
mesh, shifted, is_time_reversal)
end subroutine example_spg_get_dataset_wurtzite
subroutine show_syminfo(num_atom, lattice, symprec, atom_types, positions)
use spglib_f08, only: &
SpglibDataset, &
SpglibSpacegroupType, &
spg_get_dataset, &
spg_get_spacegroup_type, &
spg_get_ir_reciprocal_mesh
! Arguments ------------------------------------
! scalars
integer(c_int), intent(in) :: num_atom
real(c_double), intent(in) :: symprec
! arrays
integer(c_int), intent(in) :: atom_types(num_atom)
real(c_double), intent(in) :: lattice(3, 3)
real(c_double), intent(in) :: positions(3, num_atom)
! Local variables-------------------------------
! scalars
integer :: i, j, space_group, num_sym, indent
character(len=21) :: international
character(len=7) :: schoenflies
character(len=30) :: space
character(len=27) :: wyckoffs
type(SpglibDataset) :: dset
type(SpglibSpacegroupType) :: spgtype
!**************************************************************************
! The allocatable components of dset get deallocated on scope exit
dset = spg_get_dataset(lattice, positions, atom_types, num_atom, symprec)
num_sym = dset%n_operations
wyckoffs = "abcdefghijklmnopqrstuvwxyz"
space = " "
indent = 1
if (dset%spacegroup_number /= 0) then
spgtype = spg_get_spacegroup_type(dset%hall_number)
print('(a, "space_group: ", i3)'), space(1:indent*2), dset%spacegroup_number
print('(a, "international: ", a, a)'), space(1:indent*2), trim(dset%international_symbol)
print('(a, "schoenflies: ", a)'), space(1:indent*2), trim(spgtype%schoenflies)
else
print '("Space group could not be found,")'
end if
print('(a, "atom-type:")'), space(1:indent*2)
do i = 1, num_atom
print('(a, "- { type: ", i3, "}")'), space(1:indent*2), atom_types(i)
end do
print('(a, "real-basis:")'), space(1:indent*2)
do i = 1, 3
print('(a, "- [", f19.14, ", ", f19.14, ", ", f19.14, "]")'), space(1:indent*2), lattice(i, :)
end do
print('(a, "position:")'), space(1:indent*2)
do i = 1, num_atom
print('(a, "- [", f17.14, ", ", f17.14, ", ", f17.14, "]")'), space(1:indent*2), positions(:, i)
end do
print('(a, "operation:")'), space(1:indent*2)
do i = 1, num_sym
print('(a, "- rotation: #", i4)'), space(1:indent*2), i
do j = 1, 3
print('(a, " - [", i3,",", i3,",", i3,"]")'), space(1:indent*2), dset%rotations(:, j, i)
end do
print('(a, " translation: [ ", f10.7,", ", f10.7,", ", f10.7,"]")'), space(1:indent*2), dset%translations(:, i)
end do
print "(a,'wyckoffs:')", space(1:indent*2)
do i = 1, dset%n_atoms
print('(a, "- ", i0, ": ", a)'), space(1:indent*2), i, wyckoffs(dset%wyckoffs(i) + 1:dset%wyckoffs(i) + 1)
end do
print "(a,'equivalent_atoms:')", space(1:indent*2)
do i = 1, dset%n_atoms
print('(a, "- ", i0, ": ", i0)'), space(1:indent*2), i, dset%equivalent_atoms(i) + 1
end do
end subroutine show_syminfo
subroutine show_grid_info(num_atom, lattice, symprec, atom_types, positions, &
mesh, shifted, is_time_reversal)
use spglib_f08, only: spg_get_spacegroup_type, spg_get_ir_reciprocal_mesh
! Arguments ------------------------------------
! scalars
integer(c_int), intent(in) :: num_atom, is_time_reversal
real(c_double), intent(in) :: symprec
! arrays
integer(c_int), intent(in) :: atom_types(num_atom)
integer(c_int), intent(in) :: mesh(3), shifted(3)
real(c_double), intent(in) :: lattice(3, 3)
real(c_double), intent(in) :: positions(3, num_atom)
! Local variables-------------------------------
! scalars
integer :: i, counter, weight, indent
integer :: num_ir_grid
character(len=30) :: space
! arrays
integer(c_int) :: grid_point(3, mesh(1)*mesh(2)*mesh(3))
integer(c_int) :: map(mesh(1)*mesh(2)*mesh(3))
!**************************************************************************
num_ir_grid = spg_get_ir_reciprocal_mesh(grid_point, map, &
& mesh, shifted, is_time_reversal, lattice, positions, &
& atom_types, num_atom, symprec)
space = " "
indent = 1
print('(a, "reciprocal-mesh: [", i3, ", ", i3, ", ", i3, "]")'), space(1:indent*2), mesh(:)
print('(a, "- shifted: [", i3, ", ", i3, ", ", i3, "]")'), space(1:indent*2), shifted(:)
print('(a, "- is_time_reversal: ", i3)'), space(1:indent*2), is_time_reversal
print('(a, "- num-ir-grid-point:", i4)'), space(1:indent*2), num_ir_grid
print('(a, "- ir-grid-point:")'), space(1:indent*2)
counter = 0
do i = 1, product(mesh)
if (i - 1 == map(i)) then
! Ad-hoc and intuitive implementation of weight. Not optimal for very large size
weight = count(map == i - 1)
counter = counter + 1
print '(i5, 3(x,f9.6), 2x,i4,2x, f9.6, 3(x,i4))', counter, &
& real(shifted + 2*grid_point(:, i))/real(2*mesh), map(i) + 1, &
& real(weight)/real(product(mesh)), grid_point(:, i)
end if
end do
end subroutine show_grid_info
end program spglib_example
|