46#include "../base/base_uses.f90"
52 CHARACTER(len=*),
PARAMETER,
PRIVATE :: moduleN =
'space_groups'
79 SUBROUTINE spgr_create(scoor, types, cell, gopt_env, eps_symmetry, pol, ranges, &
80 nparticle, n_atom, n_core, n_shell, iunit, print_atoms)
82 REAL(kind=
dp),
DIMENSION(:, :),
INTENT(IN) :: scoor
83 INTEGER,
DIMENSION(:),
INTENT(IN) :: types
84 TYPE(
cell_type),
INTENT(IN),
POINTER :: cell
86 REAL(kind=
dp),
INTENT(IN),
OPTIONAL :: eps_symmetry
87 REAL(kind=
dp),
DIMENSION(3),
INTENT(IN),
OPTIONAL :: pol
88 INTEGER,
DIMENSION(:, :),
INTENT(IN),
OPTIONAL :: ranges
89 INTEGER,
INTENT(IN),
OPTIONAL :: nparticle, n_atom, n_core, n_shell
90 INTEGER,
INTENT(IN) :: iunit
91 LOGICAL,
INTENT(IN) :: print_atoms
93 CHARACTER(LEN=*),
PARAMETER :: routinen =
'spgr_create'
95 CHARACTER(LEN=1000) :: buffer
96 INTEGER :: ierr, nchars, nop, tra_mat(3, 3)
98 INTEGER :: handle, i, j, n_sr_rep
99 INTEGER,
ALLOCATABLE,
DIMENSION(:) :: tmp_types
101 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:, :) :: tmp_coor
104 CALL timeset(routinen, handle)
106 spgr => gopt_env%spgr
107 cpassert(
ASSOCIATED(spgr))
112 IF (
PRESENT(nparticle))
THEN
113 cpassert(nparticle ==
SIZE(scoor, 2))
114 spgr%nparticle = nparticle
116 spgr%nparticle =
SIZE(scoor, 2)
119 IF (
PRESENT(n_atom))
THEN
121 ELSE IF (
PRESENT(n_core))
THEN
122 spgr%n_atom = spgr%nparticle - n_core
123 ELSE IF (
PRESENT(n_shell))
THEN
124 spgr%n_atom = spgr%nparticle - n_shell
126 spgr%n_atom = spgr%nparticle
129 IF (
PRESENT(n_core))
THEN
131 ELSE IF (
PRESENT(n_shell))
THEN
132 spgr%n_core = n_shell
135 IF (
PRESENT(n_shell))
THEN
136 spgr%n_shell = n_shell
137 ELSE IF (
PRESENT(n_core))
THEN
138 spgr%n_shell = n_core
141 IF (.NOT. (spgr%nparticle == (spgr%n_atom + spgr%n_shell)))
THEN
142 cpabort(
"spgr_create: nparticle not equal to natom + nshell.")
145 spgr%nparticle_sym = spgr%nparticle
146 spgr%n_atom_sym = spgr%n_atom
147 spgr%n_core_sym = spgr%n_core
148 spgr%n_shell_sym = spgr%n_shell
151 spgr%print_atoms = print_atoms
154 IF (
PRESENT(eps_symmetry))
THEN
155 spgr%eps_symmetry = eps_symmetry
159 IF (
PRESENT(pol))
THEN
165 ALLOCATE (spgr%lat(spgr%nparticle))
168 IF (
PRESENT(ranges))
THEN
169 n_sr_rep =
SIZE(ranges, 2)
171 DO j = ranges(1, i), ranges(2, i)
172 spgr%lat(j) = .false.
173 spgr%nparticle_sym = spgr%nparticle_sym - 1
174 IF (j <= spgr%n_atom)
THEN
175 spgr%n_atom_sym = spgr%n_atom_sym - 1
176 ELSE IF (j > spgr%n_atom .AND. j <= spgr%nparticle)
THEN
177 spgr%n_core_sym = spgr%n_core_sym - 1
178 spgr%n_shell_sym = spgr%n_shell_sym - 1
180 cpabort(
"Symmetry exclusion range larger than actual number of particles.")
186 ALLOCATE (tmp_coor(3, spgr%n_atom_sym), tmp_types(spgr%n_atom_sym))
189 DO i = 1, spgr%n_atom
190 IF (spgr%lat(i))
THEN
192 tmp_coor(:, j) = scoor(:, i)
193 tmp_types(j) = types(i)
198 NULLIFY (spgr%cell_ref)
200 CALL cell_copy(cell, spgr%cell_ref, tag=
"CELL_OPT_REF")
201 SELECT CASE (gopt_env%type_id)
203 CALL init_cell(spgr%cell_ref, hmat=cell%hmat)
205 CALL init_cell(spgr%cell_ref, hmat=gopt_env%h_ref)
207 cpabort(
"SPACE_GROUP_SYMMETRY is not compatible with md.")
211 ALLOCATE (spgr%atype(spgr%nparticle))
212 spgr%atype(1:spgr%nparticle) = types(1:spgr%nparticle)
214 spgr%n_operations = 0
219 spgr%space_group_number =
spg_get_international(spgr%international_symbol, transpose(cell%hmat), tmp_coor, tmp_types, &
220 spgr%n_atom_sym, eps_symmetry)
222 nchars =
strlcpy_c2f(buffer, spgr%international_symbol)
223 spgr%international_symbol = buffer(1:nchars)
224 IF (spgr%space_group_number == 0)
THEN
225 cpabort(
"Symmetry Library SPGLIB failed, most likely due a problem with the coordinates.")
229 spgr%n_atom_sym, eps_symmetry)
230 ALLOCATE (spgr%rotations(3, 3, nop), spgr%translations(3, nop))
231 ALLOCATE (spgr%eqatom(nop, spgr%nparticle))
232 ALLOCATE (spgr%lop(nop))
233 spgr%n_operations = nop
235 ierr =
spg_get_symmetry(spgr%rotations, spgr%translations, nop, transpose(cell%hmat), &
236 tmp_coor, tmp_types, spgr%n_atom_sym, eps_symmetry)
239 spgr%n_atom_sym, eps_symmetry)
242 spgr%schoenflies = buffer(1:nchars)
247 spgr%rotations, spgr%n_operations)
249 nchars =
strlcpy_c2f(buffer, spgr%pointgroup_symbol)
250 spgr%pointgroup_symbol = buffer(1:nchars)
253 cpabort(
"Symmetry library SPGLIB not available")
258 DEALLOCATE (tmp_coor, tmp_types)
260 CALL timestop(handle)
281 INTEGER,
INTENT(IN) :: iunit
283 CHARACTER(LEN=*),
PARAMETER :: routinen =
'identify_space_group'
285 INTEGER :: handle, i, k, n_atom, n_core, n_shell, &
286 n_sr_rep, nparticle, shell_index
287 INTEGER,
ALLOCATABLE,
DIMENSION(:) :: atype
288 INTEGER,
ALLOCATABLE,
DIMENSION(:, :) :: ranges
289 INTEGER,
DIMENSION(:),
POINTER :: tmp
290 LOGICAL :: print_atoms
291 REAL(kind=
dp) :: eps_symmetry
292 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:, :) :: scoord
293 REAL(kind=
dp),
DIMENSION(:),
POINTER :: pol
299 CALL timeset(routinen, handle)
308 NULLIFY (core_particles)
309 NULLIFY (shell_particles)
313 cpassert(
ASSOCIATED(cell))
315 CALL cp_subsys_get(subsys, particles=particles, shell_particles=shell_particles, core_particles=core_particles)
317 cpassert(
ASSOCIATED(particles))
318 n_atom = particles%n_els
320 IF (
ASSOCIATED(shell_particles))
THEN
321 n_shell = shell_particles%n_els
322 cpassert(
ASSOCIATED(core_particles))
323 n_core = subsys%core_particles%n_els
325 cpassert(n_core == n_shell)
326 ELSE IF (
ASSOCIATED(core_particles))
THEN
328 cpabort(
"Core particles should not be defined without corresponding shell particles.")
334 nparticle = n_atom + n_shell
335 ALLOCATE (scoord(3, nparticle), atype(nparticle))
337 shell_index = particles%els(i)%shell_index
338 IF (shell_index == 0)
THEN
339 CALL real_to_scaled(scoord(1:3, i), particles%els(i)%r(1:3), cell)
340 CALL get_atomic_kind(atomic_kind=particles%els(i)%atomic_kind, kind_number=atype(i))
342 CALL real_to_scaled(scoord(1:3, i), core_particles%els(shell_index)%r(1:3), cell)
343 CALL get_atomic_kind(atomic_kind=core_particles%els(shell_index)%atomic_kind, kind_number=atype(i))
344 k = n_atom + shell_index
345 CALL real_to_scaled(scoord(1:3, k), shell_particles%els(shell_index)%r(1:3), cell)
346 CALL get_atomic_kind(atomic_kind=shell_particles%els(shell_index)%atomic_kind, kind_number=atype(k))
354 IF (n_sr_rep > 0)
THEN
355 ALLOCATE (ranges(2, n_sr_rep))
358 ranges(:, i) = tmp(:)
360 CALL spgr_create(scoord, atype, cell, gopt_env, eps_symmetry=eps_symmetry, pol=pol(1:3), &
361 ranges=ranges, nparticle=nparticle, n_atom=n_atom, &
362 n_core=n_core, n_shell=n_shell, iunit=iunit, print_atoms=print_atoms)
365 CALL spgr_create(scoord, atype, cell, gopt_env, eps_symmetry=eps_symmetry, pol=pol(1:3), &
366 nparticle=nparticle, n_atom=n_atom, &
367 n_core=n_core, n_shell=n_shell, iunit=iunit, print_atoms=print_atoms)
371 spgr => gopt_env%spgr
374 CALL spgr_reduce_symm(spgr)
375 CALL spgr_rotations_subset(spgr)
377 DEALLOCATE (scoord, atype)
379 CALL timestop(handle)
393 TYPE(
spgr_type),
INTENT(INOUT),
POINTER :: spgr
394 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:, :), &
397 CHARACTER(LEN=*),
PARAMETER :: routinen =
'spgr_find_equivalent_atoms'
399 INTEGER :: handle, i, ia, ib, ir, j, natom, nop, &
401 REAL(kind=
dp) :: diff
402 REAL(kind=
dp),
DIMENSION(3) :: rb, ri, ro, tr
403 REAL(kind=
dp),
DIMENSION(3, 3) :: rot
405 CALL timeset(routinen, handle)
407 nop = spgr%n_operations
409 nshell = spgr%n_shell
411 IF (.NOT. (spgr%nparticle == (natom + nshell)))
THEN
412 cpabort(
"spgr_find_equivalent_atoms: nparticle not equal to natom + nshell.")
415 DO ia = 1, spgr%nparticle
416 spgr%eqatom(:, ia) = ia
421 IF (.NOT. spgr%lat(ia)) cycle
422 ri(1:3) = scoord(1:3, ia)
424 rot(1:3, 1:3) = spgr%rotations(1:3, 1:3, ir)
425 tr(1:3) = spgr%translations(1:3, ir)
427 IF (.NOT. spgr%lat(ib)) cycle
428 rb(1:3) = scoord(1:3, ib)
429 ro(1) = real(rot(1, 1),
dp)*rb(1) + real(rot(2, 1),
dp)*rb(2) + real(rot(3, 1),
dp)*rb(3) + tr(1)
430 ro(2) = real(rot(1, 2),
dp)*rb(1) + real(rot(2, 2),
dp)*rb(2) + real(rot(3, 2),
dp)*rb(3) + tr(2)
431 ro(3) = real(rot(1, 3),
dp)*rb(1) + real(rot(2, 3),
dp)*rb(2) + real(rot(3, 3),
dp)*rb(3) + tr(3)
432 ro(1) = ro(1) - real(nint(ro(1) - ri(1)),
dp)
433 ro(2) = ro(2) - real(nint(ro(2) - ri(2)),
dp)
434 ro(3) = ro(3) - real(nint(ro(3) - ri(3)),
dp)
435 diff = norm2(ri(:) - ro(:))
436 IF ((diff < spgr%eps_symmetry) .AND. (spgr%atype(ia) == spgr%atype(ib)))
THEN
437 spgr%eqatom(ir, ia) = ib
448 IF (.NOT. spgr%lat(ia)) cycle
449 ri(1:3) = scoord(1:3, ia)
451 rot(1:3, 1:3) = spgr%rotations(1:3, 1:3, ir)
452 tr(1:3) = spgr%translations(1:3, ir)
455 IF (.NOT. spgr%lat(ib)) cycle
456 rb(1:3) = scoord(1:3, ib)
457 ro(1) = real(rot(1, 1),
dp)*rb(1) + real(rot(2, 1),
dp)*rb(2) + real(rot(3, 1),
dp)*rb(3) + tr(1)
458 ro(2) = real(rot(1, 2),
dp)*rb(1) + real(rot(2, 2),
dp)*rb(2) + real(rot(3, 2),
dp)*rb(3) + tr(2)
459 ro(3) = real(rot(1, 3),
dp)*rb(1) + real(rot(2, 3),
dp)*rb(2) + real(rot(3, 3),
dp)*rb(3) + tr(3)
460 ro(1) = ro(1) - real(nint(ro(1) - ri(1)),
dp)
461 ro(2) = ro(2) - real(nint(ro(2) - ri(2)),
dp)
462 ro(3) = ro(3) - real(nint(ro(3) - ri(3)),
dp)
463 diff = norm2(ri(:) - ro(:))
464 IF ((diff < spgr%eps_symmetry) .AND. (spgr%atype(ia) == spgr%atype(ib)))
THEN
465 spgr%eqatom(ir, ia) = ib
473 CALL timestop(handle)
484 SUBROUTINE spgr_reduce_symm(spgr)
486 TYPE(
spgr_type),
INTENT(INOUT),
POINTER :: spgr
488 CHARACTER(LEN=*),
PARAMETER :: routinen =
'spgr_reduce_symm'
490 INTEGER :: handle, ia, ib, ir, ja, jb, nop, nops, &
492 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:) :: x, xold
493 REAL(kind=
dp),
DIMENSION(3) :: ri, ro
494 REAL(kind=
dp),
DIMENSION(3, 3) :: rot
496 CALL timeset(routinen, handle)
498 nop = spgr%n_operations
499 nparticle = spgr%nparticle
500 ALLOCATE (x(3*nparticle), xold(3*nparticle))
504 x(ja + 1) = x(ja + 1) + spgr%pol(1)
505 x(ja + 2) = x(ja + 2) + spgr%pol(2)
506 x(ja + 3) = x(ja + 3) + spgr%pol(3)
513 spgr%lop(ir) = .true.
514 rot(1:3, 1:3) = spgr%rotations(1:3, 1:3, ir)
516 IF (.NOT. spgr%lat(ia)) cycle
518 ri(1:3) = xold(ja + 1:ja + 3)
519 ro(1) = real(rot(1, 1),
dp)*ri(1) + real(rot(2, 1),
dp)*ri(2) + real(rot(3, 1),
dp)*ri(3)
520 ro(2) = real(rot(1, 2),
dp)*ri(1) + real(rot(2, 2),
dp)*ri(2) + real(rot(3, 2),
dp)*ri(3)
521 ro(3) = real(rot(1, 3),
dp)*ri(1) + real(rot(2, 3),
dp)*ri(2) + real(rot(3, 3),
dp)*ri(3)
522 x(ja + 1:ja + 3) = ro(1:3)
525 IF (.NOT. spgr%lat(ia)) cycle
526 ib = spgr%eqatom(ir, ia)
529 ro = x(jb + 1:jb + 3) - xold(ja + 1:ja + 3)
530 spgr%lop(ir) = (spgr%lop(ir) .AND. (abs(ro(1)) < spgr%eps_symmetry) &
531 .AND. (abs(ro(2)) < spgr%eps_symmetry) &
532 .AND. (abs(ro(3)) < spgr%eps_symmetry))
534 IF (spgr%lop(ir)) nops = nops + 1
537 spgr%n_reduced_operations = nops
540 CALL timestop(handle)
542 END SUBROUTINE spgr_reduce_symm
552 SUBROUTINE spgr_rotations_subset(spgr)
554 TYPE(
spgr_type),
INTENT(INOUT),
POINTER :: spgr
556 CHARACTER(LEN=*),
PARAMETER :: routinen =
'spgr_rotations_subset'
558 INTEGER :: handle, i, j
559 INTEGER,
DIMENSION(3, 3) :: d
560 LOGICAL,
ALLOCATABLE,
DIMENSION(:) :: mask
562 CALL timeset(routinen, handle)
564 ALLOCATE (mask(spgr%n_operations))
567 DO i = 1, spgr%n_operations
568 IF (.NOT. spgr%lop(i)) mask(i) = .false.
571 DO i = 1, spgr%n_operations - 1
572 IF (.NOT. mask(i)) cycle
573 DO j = i + 1, spgr%n_operations
574 IF (.NOT. mask(j)) cycle
575 d(:, :) = spgr%rotations(:, :, j) - spgr%rotations(:, :, i)
576 IF (sum(abs(d)) == 0) mask(j) = .false.
580 spgr%n_operations_subset = 0
581 DO i = 1, spgr%n_operations
582 IF (mask(i)) spgr%n_operations_subset = spgr%n_operations_subset + 1
585 ALLOCATE (spgr%rotations_subset(3, 3, spgr%n_operations_subset))
588 DO i = 1, spgr%n_operations
591 spgr%rotations_subset(:, :, j) = spgr%rotations(:, :, i)
596 CALL timestop(handle)
598 END SUBROUTINE spgr_rotations_subset
610 TYPE(
spgr_type),
INTENT(IN),
POINTER :: spgr
611 REAL(kind=
dp),
DIMENSION(:),
INTENT(INOUT) :: coord
613 CHARACTER(LEN=*),
PARAMETER :: routinen =
'spgr_apply_rotations_coord'
615 INTEGER :: handle, ia, ib, ir, ja, jb, nop, nops, &
617 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:) :: cold
618 REAL(kind=
dp),
DIMENSION(3) :: rf, ri, rn, ro, tr
619 REAL(kind=
dp),
DIMENSION(3, 3) :: rot
621 CALL timeset(routinen, handle)
623 ALLOCATE (cold(
SIZE(coord)))
626 nop = spgr%n_operations
627 nparticle = spgr%nparticle
628 nops = spgr%n_reduced_operations
632 IF (.NOT. spgr%lat(ia)) cycle
637 IF (.NOT. spgr%lop(ir)) cycle
638 ib = spgr%eqatom(ir, ia)
639 rot(1:3, 1:3) = spgr%rotations(1:3, 1:3, ir)
640 tr(1:3) = spgr%translations(1:3, ir)
643 ro(1) = real(rot(1, 1),
dp)*ri(1) + real(rot(2, 1),
dp)*ri(2) + real(rot(3, 1),
dp)*ri(3) + tr(1)
644 ro(2) = real(rot(1, 2),
dp)*ri(1) + real(rot(2, 2),
dp)*ri(2) + real(rot(3, 2),
dp)*ri(3) + tr(2)
645 ro(3) = real(rot(1, 3),
dp)*ri(1) + real(rot(2, 3),
dp)*ri(2) + real(rot(3, 3),
dp)*ri(3) + tr(3)
646 ro(1) = ro(1) - real(nint(ro(1) - rf(1)),
dp)
647 ro(2) = ro(2) - real(nint(ro(2) - rf(2)),
dp)
648 ro(3) = ro(3) - real(nint(ro(3) - rf(3)),
dp)
649 rn(1:3) = rn(1:3) + ro(1:3)
651 rn = rn/real(nops,
dp)
657 CALL timestop(handle)
671 TYPE(
spgr_type),
INTENT(IN),
POINTER :: spgr
672 REAL(kind=
dp),
DIMENSION(:),
INTENT(INOUT) :: force
674 CHARACTER(LEN=*),
PARAMETER :: routinen =
'spgr_apply_rotations_force'
676 INTEGER :: handle, ia, ib, ir, ja, jb, nop, nops, &
678 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:) :: fold
679 REAL(kind=
dp),
DIMENSION(3) :: ri, rn, ro
680 REAL(kind=
dp),
DIMENSION(3, 3) :: rot
682 CALL timeset(routinen, handle)
684 ALLOCATE (fold(
SIZE(force)))
687 nop = spgr%n_operations
688 nparticle = spgr%nparticle
689 nops = spgr%n_reduced_operations
693 IF (.NOT. spgr%lat(ia)) cycle
697 IF (.NOT. spgr%lop(ir)) cycle
698 ib = spgr%eqatom(ir, ia)
699 rot(1:3, 1:3) = spgr%rotations(1:3, 1:3, ir)
702 ro(1) = real(rot(1, 1),
dp)*ri(1) + real(rot(2, 1),
dp)*ri(2) + real(rot(3, 1),
dp)*ri(3)
703 ro(2) = real(rot(1, 2),
dp)*ri(1) + real(rot(2, 2),
dp)*ri(2) + real(rot(3, 2),
dp)*ri(3)
704 ro(3) = real(rot(1, 3),
dp)*ri(1) + real(rot(2, 3),
dp)*ri(2) + real(rot(3, 3),
dp)*ri(3)
705 rn(1:3) = rn(1:3) + ro(1:3)
707 rn = rn/real(nops,
dp)
713 CALL timestop(handle)
725 SUBROUTINE spgr_change_basis(roti, roto, nop, h1, h2)
727 INTEGER,
DIMENSION(:, :, :) :: roti
728 REAL(kind=
dp),
DIMENSION(:, :, :) :: roto
730 REAL(kind=
dp),
DIMENSION(3, 3) :: h1, h2
732 CHARACTER(LEN=*),
PARAMETER :: routinen =
'spgr_change_basis'
734 INTEGER :: handle, ir
735 REAL(kind=
dp),
DIMENSION(3, 3) :: h1ih2, h2ih1, ih1, ih2, r, s
737 CALL timeset(routinen, handle)
741 h2ih1 = matmul(h2, ih1)
742 h1ih2 = matmul(h1, ih2)
745 r(:, :) = roti(:, :, ir)
748 roto(:, :, ir) = r(:, :)
751 CALL timestop(handle)
753 END SUBROUTINE spgr_change_basis
766 TYPE(
spgr_type),
INTENT(IN),
POINTER :: spgr
767 TYPE(
cell_type),
INTENT(IN),
POINTER :: cell
768 REAL(kind=
dp),
DIMENSION(3, 3),
INTENT(INOUT) :: stress
770 CHARACTER(LEN=*),
PARAMETER :: routinen =
'spgr_apply_rotations_stress'
772 INTEGER :: handle, i, ir, j, k, l, nop
773 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:, :, :) :: roto
774 REAL(kind=
dp),
DIMENSION(3, 3) :: hmat1, hmat2, r, stin
776 CALL timeset(routinen, handle)
778 hmat1 = transpose(cell%hmat)
785 nop = spgr%n_operations_subset
787 ALLOCATE (roto(3, 3, nop))
789 CALL spgr_change_basis(spgr%rotations_subset, roto, spgr%n_operations_subset, hmat1, hmat2)
794 r(:, :) = roto(:, :, ir)
799 stress(i, j) = stress(i, j) + (r(k, i)*r(l, j)*stin(k, l))
805 stress = stress/real(nop,
dp)
809 CALL timestop(handle)
822 TYPE(
spgr_type),
INTENT(IN),
POINTER :: spgr
826 IF (spgr%iunit > 0)
THEN
827 WRITE (spgr%iunit,
'(/,T2,A,A)')
"----------------------------------------", &
828 "---------------------------------------"
829 WRITE (spgr%iunit,
"(T2,A,T25,A,T77,A)")
"----",
"SPACE GROUP SYMMETRY INFORMATION",
"----"
830 WRITE (spgr%iunit,
'(T2,A,A)')
"----------------------------------------", &
831 "---------------------------------------"
832 IF (spgr%symlib)
THEN
833 WRITE (spgr%iunit,
'(T2,A,T73,I8)')
"SPGR| SPACE GROUP NUMBER:", &
834 spgr%space_group_number
835 WRITE (spgr%iunit,
'(T2,A,T70,A11)')
"SPGR| INTERNATIONAL SYMBOL:", &
836 trim(adjustr(spgr%international_symbol))
837 WRITE (spgr%iunit,
'(T2,A,T75,A6)')
"SPGR| POINT GROUP SYMBOL:", &
838 trim(adjustr(spgr%pointgroup_symbol))
839 WRITE (spgr%iunit,
'(T2,A,T74,A7)')
"SPGR| SCHOENFLIES SYMBOL:", &
840 trim(adjustr(spgr%schoenflies))
841 WRITE (spgr%iunit,
'(T2,A,T73,I8)')
"SPGR| NUMBER OF SYMMETRY OPERATIONS:", &
843 WRITE (spgr%iunit,
'(T2,A,T73,I8)')
"SPGR| NUMBER OF UNIQUE ROTATIONS:", &
844 spgr%n_operations_subset
845 WRITE (spgr%iunit,
'(T2,A,T73,I8)')
"SPGR| NUMBER OF REDUCED SYMMETRY OPERATIONS:", &
846 spgr%n_reduced_operations
847 WRITE (spgr%iunit,
'(T2,A,T65,I8,I8)')
"SPGR| NUMBER OF PARTICLES AND SYMMETRIZED PARTICLES:", &
848 spgr%nparticle, spgr%nparticle_sym
849 WRITE (spgr%iunit,
'(T2,A,T65,I8,I8)')
"SPGR| NUMBER OF ATOMS AND SYMMETRIZED ATOMS:", &
850 spgr%n_atom, spgr%n_atom_sym
851 WRITE (spgr%iunit,
'(T2,A,T65,I8,I8)')
"SPGR| NUMBER OF CORES AND SYMMETRIZED CORES:", &
852 spgr%n_core, spgr%n_core_sym
853 WRITE (spgr%iunit,
'(T2,A,T65,I8,I8)')
"SPGR| NUMBER OF SHELLS AND SYMMETRIZED SHELLS:", &
854 spgr%n_shell, spgr%n_shell_sym
855 IF (spgr%print_atoms)
THEN
856 WRITE (spgr%iunit, *)
"SPGR| ACTIVE REDUCED SYMMETRY OPERATIONS:", spgr%lop
857 WRITE (spgr%iunit,
'(/,T2,A,A)')
"----------------------------------------", &
858 "---------------------------------------"
859 WRITE (spgr%iunit,
'(T2,A,T34,A,T77,A)')
"----",
"EQUIVALENT ATOMS",
"----"
860 WRITE (spgr%iunit,
'(T2,A,A)')
"----------------------------------------", &
861 "---------------------------------------"
862 DO i = 1, spgr%nparticle
863 DO j = 1, spgr%n_operations
864 WRITE (spgr%iunit,
'(T2,A,T52,I8,I8,I8)')
"SPGR| ATOM | SYMMETRY OPERATION | EQUIVALENT ATOM", &
865 i, j, spgr%eqatom(j, i)
868 WRITE (spgr%iunit,
'(T2,A,A)')
"----------------------------------------", &
869 "---------------------------------------"
870 DO i = 1, spgr%n_operations
871 WRITE (spgr%iunit,
'(T2,A,T46,i4,T51,3I10,/,T51,3I10,/,T51,3I10)') &
872 "SPGR| SYMMETRY OPERATION #:", i, (spgr%rotations(j, :, i), j=1, 3)
873 WRITE (spgr%iunit,
'(T51,3F10.5)') spgr%translations(:, i)
877 WRITE (spgr%iunit,
"(T2,A)")
"SPGLIB for Crystal Symmetry Information determination is not availale"
894 REAL(kind=
dp),
DIMENSION(3, 3),
INTENT(IN) :: stress
895 TYPE(
spgr_type),
INTENT(IN),
POINTER :: spgr
897 REAL(kind=
dp),
DIMENSION(3) :: eigval
898 REAL(kind=
dp),
DIMENSION(3, 3) :: eigvec, stress_tensor
900 stress_tensor(:, :) = stress(:, :)*
pascal*1.0e-9_dp
902 IF (spgr%iunit > 0)
THEN
903 WRITE (unit=spgr%iunit, fmt=
'(/,T2,A)') &
904 'SPGR STRESS| Symmetrized stress tensor [GPa]'
905 WRITE (unit=spgr%iunit, fmt=
'(T2,A,T19,3(19X,A1))') &
906 'SPGR STRESS|',
'x',
'y',
'z'
907 WRITE (unit=spgr%iunit, fmt=
'(T2,A,T26,3(1X,ES19.11))') &
908 'SPGR STRESS| x', stress_tensor(1, 1:3)
909 WRITE (unit=spgr%iunit, fmt=
'(T2,A,T26,3(1X,ES19.11))') &
910 'SPGR STRESS| y', stress_tensor(2, 1:3)
911 WRITE (unit=spgr%iunit, fmt=
'(T2,A,T26,3(1X,ES19.11))') &
912 'SPGR STRESS| z', stress_tensor(3, 1:3)
913 WRITE (unit=spgr%iunit, fmt=
'(T2,A,T66,ES20.11)') &
914 'SPGR STRESS| 1/3 Trace', (stress_tensor(1, 1) + &
915 stress_tensor(2, 2) + &
916 stress_tensor(3, 3))/3.0_dp
917 WRITE (unit=spgr%iunit, fmt=
'(T2,A,T66,ES20.11)') &
918 'SPGR STRESS| Determinant',
det_3x3(stress_tensor(1:3, 1), &
919 stress_tensor(1:3, 2), &
920 stress_tensor(1:3, 3))
922 eigvec(:, :) = 0.0_dp
923 CALL jacobi(stress_tensor, eigval, eigvec)
924 WRITE (unit=spgr%iunit, fmt=
'(/,T2,A)') &
925 'SPGR STRESS| Eigenvectors and eigenvalues of the symmetrized stress tensor [GPa]'
926 WRITE (unit=spgr%iunit, fmt=
'(T2,A,T19,3(1X,I19))') &
927 'SPGR STRESS|', 1, 2, 3
928 WRITE (unit=spgr%iunit, fmt=
'(T2,A,T26,3(1X,ES19.11))') &
929 'SPGR STRESS| Eigenvalues', eigval(1:3)
930 WRITE (unit=spgr%iunit, fmt=
'(T2,A,T26,3(1X,F19.12))') &
931 'SPGR STRESS| x', eigvec(1, 1:3)
932 WRITE (unit=spgr%iunit, fmt=
'(T2,A,T26,3(1X,F19.12))') &
933 'SPGR STRESS| y', eigvec(2, 1:3)
934 WRITE (unit=spgr%iunit, fmt=
'(T2,A,T26,3(1X,F19.12))') &
935 'SPGR STRESS| z', eigvec(3, 1:3)
Define the atomic kind types and their sub types.
subroutine, public get_atomic_kind(atomic_kind, fist_potential, element_symbol, name, mass, kind_number, natom, atom_list, rcov, rvdw, z, qeff, apol, cpol, mm_radius, shell, shell_active, damping)
Get attributes of an atomic kind.
collects all references to literature in CP2K as new algorithms / method are included from literature...
integer, save, public togo2018
Handles all functions related to the CELL.
subroutine, public init_cell(cell, hmat, periodic)
Initialise/readjust a simulation cell after hmat has been changed.
subroutine, public cell_create(cell, hmat, periodic, tag)
allocates and initializes a cell
Handles all functions related to the CELL.
subroutine, public scaled_to_real(r, s, cell)
Transform scaled cell coordinates real coordinates. r=h*s.
subroutine, public real_to_scaled(s, r, cell)
Transform real to scaled cell coordinates. s=h_inv*r.
subroutine, public cell_copy(cell_in, cell_out, tag)
Copy cell variable.
types that represent a subsys, i.e. a part of the system
subroutine, public cp_subsys_get(subsys, ref_count, atomic_kinds, atomic_kind_set, particles, particle_set, local_particles, molecules, molecule_set, molecule_kinds, molecule_kind_set, local_molecules, para_env, colvar_p, shell_particles, core_particles, gci, multipoles, natom, nparticle, ncore, nshell, nkind, atprop, virial, results, cell, cell_ref, use_ref_cell)
returns information about various attributes of the given subsys
contains a functional that calculates the energy and its derivatives for the geometry optimizer
Defines the basic variable types.
integer, parameter, public dp
Collection of simple mathematical functions and subroutines.
subroutine, public jacobi(a, d, v)
Jacobi matrix diagonalization. The eigenvalues are returned in vector d and the eigenvectors are retu...
pure real(kind=dp) function, dimension(3, 3), public inv_3x3(a)
Returns the inverse of the 3 x 3 matrix a.
represent a simple array based list of the given type
Definition of physical constants:
real(kind=dp), parameter, public pascal
Space Group Symmetry Type Module (version 1.0, Ferbruary 12, 2021)
subroutine, public cleanup_spgr_type(spgr)
Cleanup all buffers the SPGR type.
Space Group Symmetry Module (version 1.0, January 16, 2020)
subroutine, public spgr_create(scoor, types, cell, gopt_env, eps_symmetry, pol, ranges, nparticle, n_atom, n_core, n_shell, iunit, print_atoms)
routine creates the space group structure
subroutine, public spgr_find_equivalent_atoms(spgr, scoord)
routine indentifies the equivalent atoms for each rotation matrix.
subroutine, public print_spgr(spgr)
routine prints Space Group Information.
subroutine, public spgr_apply_rotations_stress(spgr, cell, stress)
routine applies the rotation matrices to the stress tensor.
subroutine, public spgr_apply_rotations_coord(spgr, coord)
routine applies the rotation matrices to the coordinates.
subroutine, public identify_space_group(subsys, geo_section, gopt_env, iunit)
routine indentifies the space group and finds rotation matrices.
subroutine, public spgr_apply_rotations_force(spgr, force)
routine applies the rotation matrices to the forces.
subroutine, public spgr_write_stress_tensor(stress, spgr)
Variable precision output of the symmetrized stress tensor.
Interface for SPGLIB symmetry routines.
integer function, public spg_get_international(symbol, lattice, position, types, num_atom, symprec)
...
integer function, public spg_get_multiplicity(lattice, position, types, num_atom, symprec)
...
integer function, public spg_get_pointgroup(symbol, trans_mat, rotations, num_rotations)
...
integer function, public spg_get_symmetry(rotation, translation, max_size, lattice, position, types, num_atom, symprec)
...
integer function, public spg_get_schoenflies(symbol, lattice, position, types, num_atom, symprec)
...
Utilities for string manipulations.
integer function, public strlcpy_c2f(fstring, cstring)
Copy the content of a \0-terminated C-string to a finite-length Fortran string.
Type defining parameters related to the simulation cell.
represents a system: atoms, molecules, their pos,vel,...
calculates the potential energy of a system, and its derivatives
represent a list of objects