Octopus
symmetries.F90
Go to the documentation of this file.
1!! Copyright (C) 2009 X. Andrade
2!!
3!! This program is free software; you can redistribute it and/or modify
4!! it under the terms of the GNU General Public License as published by
5!! the Free Software Foundation; either version 2, or (at your option)
6!! any later version.
7!!
8!! This program is distributed in the hope that it will be useful,
9!! but WITHOUT ANY WARRANTY; without even the implied warranty of
10!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
11!! GNU General Public License for more details.
12!!
13!! You should have received a copy of the GNU General Public License
14!! along with this program; if not, write to the Free Software
15!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
16!! 02110-1301, USA.
17!!
18
19#include "global.h"
20
22 use debug_oct_m
23 use global_oct_m
24 use iso_c_binding
25 use, intrinsic :: iso_fortran_env
27 use math_oct_m
29 use mpi_oct_m
31 use parser_oct_m
33 use string_oct_m
34 use space_oct_m
35 use spglib_f08
37
38 implicit none
39
40 private
41
42 public :: &
52
53 type symmetries_t
54 private
55 type(symm_op_t), allocatable, public :: ops(:)
56 integer, public :: nops = 0
57 real(real64) :: breakdir(1:3)
58 integer, public :: periodic_dim
59 integer :: space_group
60 logical :: any_non_spherical
61 logical :: symmetries_compute = .false.
62 character(len=6) :: group_name = ""
63 character(len=30) :: group_elements = ""
64 character(len=11) :: symbol = ""
65 character(len=7) :: schoenflies = ""
66
67 type(symm_op_t), allocatable, public :: non_symmorphic_ops(:)
68 integer, public :: nops_nonsymmorphic = 0
69 contains
70 procedure :: copy => symmetries_copy
71 generic :: assignment(=) => copy
73 end type symmetries_t
74
75 real(real64), parameter, public :: default_symprec = 1.e-5_real64
76 real(real64), public :: symprec
77
79 interface
80 subroutine symmetries_finite_init(natoms, types, positions, verbosity, point_group)
81 use, intrinsic :: iso_fortran_env
82 integer, intent(in) :: natoms
83 integer, intent(in) :: types
84 real(real64), intent(in) :: positions
85 integer, intent(in) :: verbosity
86 integer, intent(out) :: point_group
87 end subroutine symmetries_finite_init
88
89 subroutine symmetries_finite_get_group_name(point_group, name)
90 integer, intent(in) :: point_group
91 character(len=*), intent(out) :: name
93
94 subroutine symmetries_finite_get_group_elements(point_group, elements)
95 integer, intent(in) :: point_group
96 character(len=*), intent(out) :: elements
98
99 subroutine symmetries_finite_end()
100 end subroutine symmetries_finite_end
101 end interface
102
103 interface symmetries_t
104 module procedure symmetries_constructor
105 end interface symmetries_t
106
107contains
108
109 function symmetries_constructor(namespace, space, latt, n_sites, site_pos, site_type, spherical_site) result(this)
110 type(namespace_t), intent(in) :: namespace
111 class(space_t), intent(in) :: space
112 type(lattice_vectors_t), intent(in) :: latt
113 integer, intent(in) :: n_sites
114 real(real64), intent(in) :: site_pos(:,:)
115 integer, intent(in) :: site_type(1:n_sites)
116 logical, intent(in) :: spherical_site(:)
117 type(symmetries_t) :: this
118
119 integer :: dim4syms
120 integer :: idir, isite, isite_symm, iop, verbosity, point_group
121 real(real64), allocatable :: position(:, :)
122 type(block_t) :: blk
123 type(symm_op_t) :: tmpop
124 integer :: identity(3,3)
125 logical :: found_identity, is_supercell
126 logical :: def_sym_comp
127 real(real64) :: rsite(space%dim)
128
129 type(SpglibDataset) :: spg_dataset
130 type(SpglibSpacegroupType) :: spg_spacegroup
131
132 push_sub(symmetries_constructor)
133
134 ! if someone cares, they could try to analyze the symmetry point group of the individual species too
135 this%any_non_spherical = any(.not. spherical_site)
136
137 this%periodic_dim = space%periodic_dim
138
139 dim4syms = min(3, space%dim)
140
141 def_sym_comp = (n_sites < 100) .or. space%is_periodic()
142 def_sym_comp = def_sym_comp .and. space%dim == 3
143
144 !%Variable SymmetriesCompute
145 !%Type logical
146 !%Section Execution::Symmetries
147 !%Description
148 !% If disabled, <tt>Octopus</tt> will not compute
149 !% nor print the symmetries.
150 !%
151 !% By default, symmetries are computed when running in 3
152 !% dimensions for systems with less than 100 atoms.
153 !% For periodic systems, the default is always true, irrespective of the number of atoms.
154 !%End
155 call parse_variable(namespace, 'SymmetriesCompute', def_sym_comp, this%symmetries_compute)
157 !%Variable SymmetriesTolerance
158 !%Type float
159 !%Section Execution::Symmetries
160 !%Description
161 !% For periodic systems, this variable controls the tolerance used by the symmetry finder
162 !% (spglib) to find the spacegroup and symmetries of the crystal.
163 !%End
164 call parse_variable(namespace, 'SymmetriesTolerance', default_symprec, symprec)
166
167 if (this%symmetries_compute .and. space%dim /= 3) then
168 call messages_experimental('symmetries for non 3D systems', namespace=namespace)
169 end if
170
171 if (this%any_non_spherical .or. .not. this%symmetries_compute) then
172 call init_identity()
174 return
175 end if
176
177 ! In all cases, we must check that the grid respects the symmetries. --DAS
178
179 if (space%periodic_dim == 0) then
180
181 call init_identity()
183 ! for the moment symmetries are only used for information, so we compute them only on one node.
184 if (mpi_grp_is_root(mpi_world)) then
185 verbosity = -1
186
187 safe_allocate(position(1:3, 1:n_sites))
188
189 do isite = 1, n_sites
190 position(1:3, isite) = m_zero
191 position(1:dim4syms, isite) = site_pos(1:dim4syms, isite)
192 end do
193
194 if (this%symmetries_compute) then
195 call symmetries_finite_init(n_sites, site_type(1), position(1, 1), verbosity, point_group)
196 if (point_group > -1) then
197 call symmetries_finite_get_group_name(point_group, this%group_name)
198 call symmetries_finite_get_group_elements(point_group, this%group_elements)
199 else
200 this%group_name = ""
201 this%group_elements = ""
202 this%symbol = ""
203 this%schoenflies = ""
204 this%symmetries_compute = .false.
205 end if
207 end if
208 safe_deallocate_a(position)
209 end if
210
211 else
212
213 spg_dataset = symmetries_get_spg_dataset(namespace, space, latt, n_sites, site_pos, site_type)
214 if (spg_dataset%spglib_error /= 0) then
215 call init_identity()
217 return
218 end if
219
220 this%space_group = spg_dataset%spacegroup_number
221 this%symbol = spg_dataset%international_symbol
222 spg_spacegroup = spg_get_spacegroup_type(spg_dataset%hall_number)
223 this%schoenflies = spg_spacegroup%schoenflies
224
225 do iop = 1, spg_dataset%n_operations
226 ! transpose due to array order difference between C and fortran
227 spg_dataset%rotations(:,:,iop) = transpose(spg_dataset%rotations(:,:,iop))
228 ! sometimes spglib may return lattice vectors as 'fractional' translations
229 spg_dataset%translations(:, iop) = spg_dataset%translations(:, iop) - &
230 anint(spg_dataset%translations(:, iop) + m_half * tol_translation)
231 end do
232
233 ! we need to check that it is not a supercell, as in the QE routine (sgam_at)
234 ! they disable fractional translations if the identity has one, because the sym ops might not form a group.
235 ! spglib may return duplicate operations in this case!
236 is_supercell = (spg_dataset%n_operations > 48)
237 found_identity = .false.
238 identity = diagonal_matrix(3, 1)
239 do iop = 1, spg_dataset%n_operations
240 if (all(spg_dataset%rotations(1:3, 1:3, iop) == identity(1:3, 1:3))) then
241 found_identity = .true.
242 if (any(abs(spg_dataset%translations(1:3, iop)) > tol_translation )) then
243 is_supercell = .true.
244 write(message(1),'(a,3f12.6)') 'Identity has a fractional translation ', spg_dataset%translations(1:3, iop)
245 call messages_info(1, namespace=namespace)
246 end if
247 end if
248 end do
249 if (.not. found_identity) then
250 message(1) = "Symmetries internal error: Identity is missing from symmetry operations."
251 call messages_fatal(1, namespace=namespace)
252 end if
253
254 if (is_supercell) then
255 message(1) = "Disabling fractional translations. System appears to be a supercell."
256 call messages_info(1, namespace=namespace)
257 end if
258 ! actually, we do not use fractional translations regardless currently
259
260 ! this is a hack to get things working, this variable should be
261 ! eliminated and the direction calculated automatically from the
262 ! perturbations.
263
264 !%Variable SymmetryBreakDir
265 !%Type block
266 !%Section Mesh::Simulation Box
267 !%Description
268 !% This variable specifies a direction in which the symmetry of
269 !% the system will be broken. This is useful for generating <i>k</i>-point
270 !% grids when an external perturbation is applied.
271 !%End
272
273 this%breakdir(1:3) = m_zero
274
275 if (parse_block(namespace, 'SymmetryBreakDir', blk) == 0) then
276
277 do idir = 1, dim4syms
278 call parse_block_float(blk, 0, idir - 1, this%breakdir(idir))
279 end do
280
281 call parse_block_end(blk)
282
283 end if
284
285 safe_allocate(this%ops(1:spg_dataset%n_operations))
286 safe_allocate(this%non_symmorphic_ops(1:spg_dataset%n_operations))
287
288 ! check all operations and leave those that kept the symmetry-breaking
289 ! direction invariant and (for the moment) that do not have a translation
290 this%nops = 0
291 ! We also keep track of all non-symmorphic operations
292 this%nops_nonsymmorphic = 0
293 do iop = 1, spg_dataset%n_operations
294 ! If spglib returns tiny translations, we set them to zero to avoid numerical artifacts
295 if (all(abs(spg_dataset%translations(1:3, iop)) < tol_translation)) then
296 spg_dataset%translations(1:3, iop) = m_zero
297 end if
298
299 call symm_op_init(tmpop, spg_dataset%rotations(1:3, 1:3, iop), latt, dim4syms, spg_dataset%translations(1:3, iop))
300
301 ! If we do not specify SymmetryBreakDir, we want to use all symmetries
302 if( all(abs(this%breakdir) < m_epsilon) ) then
303 if (.not. symm_op_has_translation(tmpop, tol_translation)) then
304 this%nops = this%nops + 1
305 call symm_op_copy(tmpop, this%ops(this%nops))
306 else
307 this%nops_nonsymmorphic = this%nops_nonsymmorphic + 1
308 call symm_op_copy(tmpop, this%non_symmorphic_ops(this%nops_nonsymmorphic))
309 end if
310 else ! We only test symmorphic symmetries
311 if (symm_op_invariant_cart(tmpop, this%breakdir, symprec) .and. &
312 .not. symm_op_has_translation(tmpop, symprec)) then
313 this%nops = this%nops + 1
314 call symm_op_copy(tmpop, this%ops(this%nops))
315 end if
316 end if
317 end do
318
319
320 end if
321
322 call symmetries_write_info(this, space, namespace=namespace)
323
324 ! Checks if that the atomic coordinates are compatible with the symmetries
325 !
326 ! We want to use for instance that
327 !
328 ! \int dr f(Rr) V_isite(r) \nabla f(R(v)) = R\int dr f(r) V_isite(R*r) f(r)
329 !
330 ! and that the operator R should map the position of atom
331 ! isite to the position of some other atom isite_symm, so that
332 !
333 ! V_isite(R*r) = V_isite_symm(r)
334 !
335 do iop = 1, symmetries_number(this)
336 if (iop == symmetries_identity_index(this)) cycle
337
338 do isite = 1, n_sites
339 rsite = symm_op_apply_cart(this%ops(iop), site_pos(:, isite))
340
341 rsite = latt%fold_into_cell(rsite)
342
343 ! find isite_symm
344 do isite_symm = 1, n_sites
345 if (all(abs(rsite - site_pos(:, isite_symm)) < symprec)) exit
346 end do
347
348 if (isite_symm > n_sites) then
349 write(message(1),'(a,i6)') 'Internal error: could not find symmetric partner for atom number', isite
350 write(message(2),'(a,i3,a)') 'with symmetry operation number ', iop, '.'
351 call messages_fatal(2, namespace=namespace)
352 end if
353
354 end do
355 end do
356
357 ! Check non-symmorphic operations
358 do iop = 1, this%nops_nonsymmorphic
359 do isite = 1, n_sites
360 rsite = symm_op_apply_cart(this%non_symmorphic_ops(iop), site_pos(:, isite))
361
362 rsite = latt%fold_into_cell(rsite)
363
364 ! find isite_symm
365 do isite_symm = 1, n_sites
366 if (all(abs(rsite - site_pos(:, isite_symm)) < symprec)) exit
367 end do
368
369 if (isite_symm > n_sites) then
370 write(message(1),'(a,i6)') 'Internal error: could not find symmetric partner for atom number', isite
371 write(message(2),'(a,i3,a)') 'with non-symmorphic symmetry operation number ', iop, '.'
372 call messages_fatal(2, namespace=namespace)
373 end if
374
375 end do
376 end do
377
379
380 contains
381
382 subroutine init_identity()
383
384 push_sub(symmetries_init.init_identity)
385
386 safe_allocate(this%ops(1))
387 this%nops = 1
388 call symm_op_init(this%ops(1), diagonal_matrix(3, 1), latt, dim4syms)
389 this%nops_nonsymmorphic = 0
390 this%breakdir = m_zero
391 this%space_group = 1
392
393 pop_sub(symmetries_init.init_identity)
394
395 end subroutine init_identity
396
397 end function symmetries_constructor
398
399 ! -------------------------------------------------------------------------------
400 subroutine symmetries_copy(lhs, rhs)
401 class(symmetries_t), intent(out) :: lhs
402 class(symmetries_t), intent(in) :: rhs
403
404 integer :: iop
405
406 push_sub(symmetries_copy)
407
408 safe_allocate(lhs%ops(1:rhs%nops))
409 do iop = 1, rhs%nops
410 call symm_op_copy(rhs%ops(iop), lhs%ops(iop))
411 end do
412 lhs%nops = rhs%nops
413 safe_allocate(lhs%non_symmorphic_ops(1:rhs%nops_nonsymmorphic))
414 do iop = 1, rhs%nops_nonsymmorphic
415 call symm_op_copy(rhs%non_symmorphic_ops(iop), lhs%non_symmorphic_ops(iop))
416 end do
417 lhs%nops_nonsymmorphic = rhs%nops_nonsymmorphic
418 lhs%breakdir = rhs%breakdir
419 lhs%periodic_dim = rhs%periodic_dim
420 lhs%space_group = rhs%space_group
421 lhs%any_non_spherical = rhs%any_non_spherical
422 lhs%symmetries_compute = rhs%symmetries_compute
423 lhs%group_name = rhs%group_name
424 lhs%group_elements = rhs%group_elements
425 lhs%symbol = rhs%symbol
426 lhs%schoenflies = rhs%schoenflies
427
428 pop_sub(symmetries_copy)
429 end subroutine symmetries_copy
430
431 ! -------------------------------------------------------------------------------
432
433 subroutine symmetries_finalizer(this)
434 type(symmetries_t), intent(inout) :: this
435
436 push_sub(symmetries_finalizer)
437
438 safe_deallocate_a(this%ops)
439 safe_deallocate_a(this%non_symmorphic_ops)
440
441 pop_sub(symmetries_finalizer)
442 end subroutine symmetries_finalizer
443
444 ! -------------------------------------------------------------------------------
445
446 integer pure function symmetries_number(this) result(number)
447 type(symmetries_t), intent(in) :: this
448
449 number = this%nops
450 end function symmetries_number
451
452 ! -------------------------------------------------------------------------------
453
454 subroutine symmetries_apply_kpoint_red(this, iop, aa, bb)
455 type(symmetries_t), intent(in) :: this
456 integer, intent(in) :: iop
457 real(real64), intent(in) :: aa(1:3)
458 real(real64), intent(out) :: bb(1:3)
459
461
462 assert(iop /= 0 .and. abs(iop) <= this%nops)
463
464 !Time-reversal symmetry
465 if (iop < 0) then
466 bb(1:3) = symm_op_apply_transpose_red(this%ops(abs(iop)), -aa(1:3))
467 else
468 bb(1:3) = symm_op_apply_transpose_red(this%ops(abs(iop)), aa(1:3))
469 end if
470
472 end subroutine symmetries_apply_kpoint_red
473
474 ! -------------------------------------------------------------------------------
476 integer pure function symmetries_space_group_number(this) result(number)
477 type(symmetries_t), intent(in) :: this
478
479 number = this%space_group
481
482 ! -------------------------------------------------------------------------------
483
484 logical pure function symmetries_have_break_dir(this) result(have)
485 type(symmetries_t), intent(in) :: this
486
487 have = any(abs(this%breakdir(1:3)) > m_epsilon)
488 end function symmetries_have_break_dir
489
490 ! -------------------------------------------------------------------------------
491
492 integer pure function symmetries_identity_index(this) result(index)
493 type(symmetries_t), intent(in) :: this
494
495 integer :: iop
496
497 do iop = 1, this%nops
498 if (symm_op_is_identity(this%ops(iop))) then
499 index = iop
500 cycle
501 end if
502 end do
503
504 end function symmetries_identity_index
505
506 ! ---------------------------------------------------------
507 subroutine symmetries_write_info(this, space, iunit, namespace)
508 type(symmetries_t), intent(in) :: this
509 class(space_t), intent(in) :: space
510 integer, optional, intent(in) :: iunit
511 type(namespace_t), optional, intent(in) :: namespace
512
513 integer :: iop
514
515 push_sub(symmetries_write_info)
516
517 call messages_print_with_emphasis(msg='Symmetries', iunit=iunit, namespace=namespace)
518
519 if (this%any_non_spherical) then
520 message(1) = "Symmetries are disabled since non-spherically symmetric species may be present."
521 call messages_info(1,iunit = iunit, namespace=namespace)
522 call messages_print_with_emphasis(iunit=iunit, namespace=namespace)
523 pop_sub(symmetries_write_info)
524 return
525 end if
527 if (.not. this%symmetries_compute) then
528 message(1) = "Symmetries have been disabled by SymmetriesCompute = false."
529 call messages_info(1, iunit=iunit, namespace=namespace)
530 call messages_print_with_emphasis(iunit=iunit, namespace=namespace)
531 pop_sub(symmetries_write_info)
532 return
533 end if
534
535 if (space%periodic_dim == 0) then
536 ! At the moment only the root node has information about symetries of finite systems.
537 if (mpi_grp_is_root(mpi_world)) then
538 if (this%symmetries_compute) then
539 call messages_write('Symmetry elements : '//trim(this%group_elements), new_line = .true.)
540 call messages_write('Symmetry group : '//trim(this%group_name))
541 call messages_info(iunit=iunit, namespace=namespace)
542 end if
543 end if
544 else
545 write(message(1),'(a, i4)') 'Space group No. ', this%space_group
546 write(message(2),'(2a)') 'International: ', trim(this%symbol)
547 write(message(3),'(2a)') 'Schoenflies: ', trim(this%schoenflies)
548 call messages_info(3, iunit=iunit, namespace=namespace)
549
550 write(message(1), '(a,i5,a)') 'Info: The system has ', this%nops, ' symmorphic symmetries that can be used.'
551 call messages_info(1, iunit=iunit, namespace=namespace)
552
553 write(message(1),'(a7,a31,12x,a33)') 'Index', 'Rotation matrix', 'Fractional translations'
554 call messages_info(1, iunit=iunit, namespace=namespace)
555 do iop = 1, this%nops
556 ! list all operations and leave those that kept the symmetry-breaking
557 ! direction invariant and (for the moment) that do not have a translation
558 if (space%dim == 1) then
559 write(message(1),'(i5,1x,a,2x,1(1i4,2x),1f12.6)') iop, ':', symm_op_rotation_matrix_red(this%ops(iop)), &
560 symm_op_translation_vector_red(this%ops(iop))
561 end if
562 if (space%dim == 2) then
563 write(message(1),'(i5,1x,a,2x,2(2i4,2x),2f12.6)') iop, ':', symm_op_rotation_matrix_red(this%ops(iop)), &
564 symm_op_translation_vector_red(this%ops(iop))
565 end if
566 if (space%dim == 3) then
567 write(message(1),'(i5,1x,a,2x,3(3i4,2x),3f12.6)') iop, ':', symm_op_rotation_matrix_red(this%ops(iop)), &
568 symm_op_translation_vector_red(this%ops(iop))
569 end if
570 call messages_info(1, iunit=iunit, namespace=namespace)
571 end do
572
573 write(message(1), '(a,i5,a)') 'Info: The system also has ', this%nops_nonsymmorphic, &
574 ' nonsymmorphic symmetries (not used for electrons).'
575 call messages_info(1, iunit=iunit, namespace=namespace)
576 do iop = 1, this%nops_nonsymmorphic
577 ! list all operations and leave those that kept the symmetry-breaking
578 ! direction invariant and (for the moment) that do not have a translation
579 if (space%dim == 1) then
580 write(message(1),'(i5,1x,a,2x,1(1i4,2x),1f12.6)') iop, ':', symm_op_rotation_matrix_red(this%non_symmorphic_ops(iop)), &
581 symm_op_translation_vector_red(this%non_symmorphic_ops(iop))
582 end if
583 if (space%dim == 2) then
584 write(message(1),'(i5,1x,a,2x,2(2i4,2x),2f12.6)') iop, ':', symm_op_rotation_matrix_red(this%non_symmorphic_ops(iop)), &
585 symm_op_translation_vector_red(this%non_symmorphic_ops(iop))
586 end if
587 if (space%dim == 3) then
588 write(message(1),'(i5,1x,a,2x,3(3i4,2x),3f12.6)') iop, ':', symm_op_rotation_matrix_red(this%non_symmorphic_ops(iop)), &
589 symm_op_translation_vector_red(this%non_symmorphic_ops(iop))
590 end if
591 call messages_info(1, iunit=iunit, namespace=namespace)
592 end do
593 end if
594 call messages_print_with_emphasis(iunit=iunit, namespace=namespace)
595
596 pop_sub(symmetries_write_info)
597 end subroutine symmetries_write_info
598
599 ! ---------------------------------------------------------
601 subroutine symmetries_update_lattice_vectors(this, latt, dim)
602 type(symmetries_t), intent(inout) :: this
603 type(lattice_vectors_t), intent(in) :: latt
604 integer, intent(in) :: dim
605
606 integer :: iop
607
609
610 do iop = 1, this%nops
611 call symm_op_build_cartesian(this%ops(iop), latt, dim)
612 end do
613
615 end subroutine
616
618 type(spglibdataset) function symmetries_get_spg_dataset(namespace, space, latt, n_sites, site_pos, site_type) result(spg_dataset)
619 type(namespace_t), intent(in) :: namespace
620 class(space_t), intent(in) :: space
621 type(lattice_vectors_t), intent(in) :: latt
622 integer, intent(in) :: n_sites
623 real(real64), intent(in) :: site_pos(:,:)
624 integer, intent(in) :: site_type(1:n_sites)
625
626 integer :: isite, idir, dim4syms
627 real(real64), allocatable :: position(:, :)
628 real(real64) :: lattice(1:3, 1:3)
629
631
632 assert(space%is_periodic())
633
634 dim4syms = min(3, space%dim)
635
636 safe_allocate(position(1:3, 1:n_sites))
637
638 do isite = 1, n_sites
639 position(1:3,isite) = m_zero
640
641 ! Transform atomic positions to reduced coordinates
642 position(1:dim4syms, isite) = latt%cart_to_red(site_pos(1:dim4syms, isite))
643 position(1:dim4syms, isite) = position(1:dim4syms, isite) - m_half
644 do idir = 1, dim4syms
645 position(idir, isite) = position(idir, isite) - anint(position(idir, isite))
646 end do
647 position(1:dim4syms, isite) = position(1:dim4syms,isite) + m_half
648 end do
649
650 lattice = m_zero
651 !NOTE: Why "inverse matrix" ? (NTD)
652 ! get inverse matrix to extract reduced coordinates for spglib
653 lattice(1:space%periodic_dim, 1:space%periodic_dim) = latt%rlattice(1:space%periodic_dim, 1:space%periodic_dim)
654 ! transpose the lattice vectors for use in spglib as row-major matrix
655 lattice(:,:) = transpose(lattice(:,:))
656 ! spglib always assume 3D periodicity, so we need to set the lattice vectors along the non-periodic dimensions.
657 ! Here we choose a value that should guarantee that no incorrect symmetries are found along the non-periodic dimensions.
658 do idir = space%periodic_dim + 1, 3
659 if (idir <= space%dim) then
660 lattice(idir, idir) = m_two*(maxval(site_pos(idir, :)) - minval(site_pos(idir, :))) + m_one
661 else
662 lattice(idir, idir) = m_one
663 end if
664 end do
665
666 spg_dataset = spg_get_dataset(lattice, position, site_type, n_sites, symprec)
667 if (spg_dataset%spglib_error /= 0) then
668 message(1) = "Symmetry analysis failed in spglib. Disabling symmetries."
669 call messages_warning(1, namespace=namespace)
670
671 do isite = 1, n_sites
672 write(message(1),'(a,i6,a,3f12.6,a,3f12.6)') 'type ', site_type(isite), &
673 ' reduced coords ', position(:, isite), ' cartesian coords ', site_pos(:, isite)
674 call messages_info(1, namespace=namespace)
675 end do
676 end if
677 safe_deallocate_a(position)
678
680 end function symmetries_get_spg_dataset
681
682end module symmetries_oct_m
683
684!! Local Variables:
685!! mode: f90
686!! coding: utf-8
687!! End:
NOTE: unfortunately, these routines use global variables shared among them.
Definition: symmetries.F90:173
real(real64), parameter, public m_zero
Definition: global.F90:187
real(real64), parameter, public m_epsilon
Definition: global.F90:203
real(real64), parameter, public m_half
Definition: global.F90:193
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:115
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:160
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:420
subroutine, public messages_experimental(name, namespace)
Definition: messages.F90:1097
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:624
logical function mpi_grp_is_root(grp)
Is the current MPI process of grpcomm, root.
Definition: mpi.F90:430
type(mpi_grp_t), public mpi_world
Definition: mpi.F90:266
integer function, public parse_block(namespace, name, blk, check_varinfo_)
Definition: parser.F90:618
type(spglibspacegrouptype) function, public spg_get_spacegroup_type(hall_number)
Definition: spglib_f08.F90:415
subroutine, public symm_op_copy(inp, outp)
Definition: symm_op.F90:276
logical pure function, public symm_op_has_translation(this, prec)
Definition: symm_op.F90:293
real(real64), parameter, public tol_translation
Definition: symm_op.F90:183
subroutine, public symm_op_init(this, rot, latt, dim, trans)
Definition: symm_op.F90:190
subroutine, public symmetries_update_lattice_vectors(this, latt, dim)
Updates the symmetry operations when lattice vectors are updated.
Definition: symmetries.F90:695
subroutine, public symmetries_apply_kpoint_red(this, iop, aa, bb)
Definition: symmetries.F90:548
integer pure function, public symmetries_space_group_number(this)
Definition: symmetries.F90:570
integer pure function, public symmetries_identity_index(this)
Definition: symmetries.F90:586
type(symmetries_t) function symmetries_constructor(namespace, space, latt, n_sites, site_pos, site_type, spherical_site)
Definition: symmetries.F90:203
subroutine, public symmetries_write_info(this, space, iunit, namespace)
Definition: symmetries.F90:601
real(real64), public symprec
Definition: symmetries.F90:169
integer pure function, public symmetries_number(this)
Definition: symmetries.F90:540
subroutine symmetries_copy(lhs, rhs)
Definition: symmetries.F90:494
logical pure function, public symmetries_have_break_dir(this)
Definition: symmetries.F90:578
type(spglibdataset) function, public symmetries_get_spg_dataset(namespace, space, latt, n_sites, site_pos, site_type)
Returns the SpglibDataset for the system.
Definition: symmetries.F90:712
subroutine symmetries_finalizer(this)
Definition: symmetries.F90:527
subroutine init_identity()
Definition: symmetries.F90:476
int true(void)