49#include "../base/base_uses.f90"
55 CHARACTER(len=*),
PARAMETER,
PRIVATE :: moduleN =
'space_groups'
82 SUBROUTINE spgr_create(scoor, types, cell, gopt_env, eps_symmetry, pol, ranges, &
83 nparticle, n_atom, n_core, n_shell, iunit, print_atoms)
85 REAL(kind=
dp),
DIMENSION(:, :),
INTENT(IN) :: scoor
86 INTEGER,
DIMENSION(:),
INTENT(IN) :: types
87 TYPE(
cell_type),
INTENT(IN),
POINTER :: cell
89 REAL(kind=
dp),
INTENT(IN),
OPTIONAL :: eps_symmetry
90 REAL(kind=
dp),
DIMENSION(3),
INTENT(IN),
OPTIONAL :: pol
91 INTEGER,
DIMENSION(:, :),
INTENT(IN),
OPTIONAL :: ranges
92 INTEGER,
INTENT(IN),
OPTIONAL :: nparticle, n_atom, n_core, n_shell
93 INTEGER,
INTENT(IN) :: iunit
94 LOGICAL,
INTENT(IN) :: print_atoms
96 CHARACTER(LEN=*),
PARAMETER :: routinen =
'spgr_create'
98 CHARACTER(LEN=1000) :: buffer
99 INTEGER :: ierr, nchars, nop, tra_mat(3, 3)
101 INTEGER :: handle, i, j, n_sr_rep
102 INTEGER,
ALLOCATABLE,
DIMENSION(:) :: tmp_types
104 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:, :) :: tmp_coor
107 CALL timeset(routinen, handle)
109 spgr => gopt_env%spgr
110 cpassert(
ASSOCIATED(spgr))
115 IF (
PRESENT(nparticle))
THEN
116 cpassert(nparticle ==
SIZE(scoor, 2))
117 spgr%nparticle = nparticle
119 spgr%nparticle =
SIZE(scoor, 2)
122 IF (
PRESENT(n_atom))
THEN
124 ELSE IF (
PRESENT(n_core))
THEN
125 spgr%n_atom = spgr%nparticle - n_core
126 ELSE IF (
PRESENT(n_shell))
THEN
127 spgr%n_atom = spgr%nparticle - n_shell
129 spgr%n_atom = spgr%nparticle
132 IF (
PRESENT(n_core))
THEN
134 ELSE IF (
PRESENT(n_shell))
THEN
135 spgr%n_core = n_shell
138 IF (
PRESENT(n_shell))
THEN
139 spgr%n_shell = n_shell
140 ELSE IF (
PRESENT(n_core))
THEN
141 spgr%n_shell = n_core
144 IF (.NOT. (spgr%nparticle == (spgr%n_atom + spgr%n_shell)))
THEN
145 cpabort(
"spgr_create: nparticle not equal to natom + nshell.")
148 spgr%nparticle_sym = spgr%nparticle
149 spgr%n_atom_sym = spgr%n_atom
150 spgr%n_core_sym = spgr%n_core
151 spgr%n_shell_sym = spgr%n_shell
154 spgr%print_atoms = print_atoms
157 IF (
PRESENT(eps_symmetry))
THEN
158 spgr%eps_symmetry = eps_symmetry
162 IF (
PRESENT(pol))
THEN
168 ALLOCATE (spgr%lat(spgr%nparticle))
171 IF (
PRESENT(ranges))
THEN
172 n_sr_rep =
SIZE(ranges, 2)
174 DO j = ranges(1, i), ranges(2, i)
175 spgr%lat(j) = .false.
176 spgr%nparticle_sym = spgr%nparticle_sym - 1
177 IF (j <= spgr%n_atom)
THEN
178 spgr%n_atom_sym = spgr%n_atom_sym - 1
179 ELSE IF (j > spgr%n_atom .AND. j <= spgr%nparticle)
THEN
180 spgr%n_core_sym = spgr%n_core_sym - 1
181 spgr%n_shell_sym = spgr%n_shell_sym - 1
183 cpabort(
"Symmetry exclusion range larger than actual number of particles.")
189 ALLOCATE (tmp_coor(3, spgr%n_atom_sym), tmp_types(spgr%n_atom_sym))
192 DO i = 1, spgr%n_atom
193 IF (spgr%lat(i))
THEN
195 tmp_coor(:, j) = scoor(:, i)
196 tmp_types(j) = types(i)
201 NULLIFY (spgr%cell_ref)
203 CALL cell_copy(cell, spgr%cell_ref, tag=
"CELL_OPT_REF")
204 SELECT CASE (gopt_env%type_id)
206 CALL init_cell(spgr%cell_ref, hmat=cell%hmat)
208 SELECT CASE (gopt_env%cell_method_id)
210 CALL init_cell(spgr%cell_ref, hmat=gopt_env%h_ref)
212 cpabort(
"SPACE_GROUP_SYMMETRY should not be invoked during the cell step.")
214 cpabort(
"SPACE_GROUP_SYMMETRY is not compatible with md.")
216 cpabort(
"SPACE_GROUP_SYMMETRY invoked with an unknown optimization method.")
219 cpabort(
"SPACE_GROUP_SYMMETRY is not compatible with md.")
223 ALLOCATE (spgr%atype(spgr%nparticle))
224 spgr%atype(1:spgr%nparticle) = types(1:spgr%nparticle)
226 spgr%n_operations = 0
231 spgr%space_group_number =
spg_get_international(spgr%international_symbol, transpose(cell%hmat), tmp_coor, tmp_types, &
232 spgr%n_atom_sym, eps_symmetry)
234 nchars =
strlcpy_c2f(buffer, spgr%international_symbol)
235 spgr%international_symbol = buffer(1:nchars)
236 IF (spgr%space_group_number == 0)
THEN
237 cpabort(
"Symmetry Library SPGLIB failed, most likely due a problem with the coordinates.")
241 spgr%n_atom_sym, eps_symmetry)
242 ALLOCATE (spgr%rotations(3, 3, nop), spgr%translations(3, nop))
243 ALLOCATE (spgr%eqatom(nop, spgr%nparticle))
244 ALLOCATE (spgr%lop(nop))
245 spgr%n_operations = nop
247 ierr =
spg_get_symmetry(spgr%rotations, spgr%translations, nop, transpose(cell%hmat), &
248 tmp_coor, tmp_types, spgr%n_atom_sym, eps_symmetry)
251 spgr%n_atom_sym, eps_symmetry)
254 spgr%schoenflies = buffer(1:nchars)
259 spgr%rotations, spgr%n_operations)
261 nchars =
strlcpy_c2f(buffer, spgr%pointgroup_symbol)
262 spgr%pointgroup_symbol = buffer(1:nchars)
265 cpabort(
"Symmetry library SPGLIB not available")
270 DEALLOCATE (tmp_coor, tmp_types)
272 CALL timestop(handle)
293 INTEGER,
INTENT(IN) :: iunit
295 CHARACTER(LEN=*),
PARAMETER :: routinen =
'identify_space_group'
297 INTEGER :: handle, i, k, n_atom, n_core, n_shell, &
298 n_sr_rep, nparticle, shell_index
299 INTEGER,
ALLOCATABLE,
DIMENSION(:) :: atype
300 INTEGER,
ALLOCATABLE,
DIMENSION(:, :) :: ranges
301 INTEGER,
DIMENSION(:),
POINTER :: tmp
302 LOGICAL :: print_atoms
303 REAL(kind=
dp) :: eps_symmetry
304 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:, :) :: scoord
305 REAL(kind=
dp),
DIMENSION(:),
POINTER :: pol
311 CALL timeset(routinen, handle)
320 NULLIFY (core_particles)
321 NULLIFY (shell_particles)
325 cpassert(
ASSOCIATED(cell))
327 CALL cp_subsys_get(subsys, particles=particles, shell_particles=shell_particles, core_particles=core_particles)
329 cpassert(
ASSOCIATED(particles))
330 n_atom = particles%n_els
332 IF (
ASSOCIATED(shell_particles))
THEN
333 n_shell = shell_particles%n_els
334 cpassert(
ASSOCIATED(core_particles))
335 n_core = subsys%core_particles%n_els
337 cpassert(n_core == n_shell)
338 ELSE IF (
ASSOCIATED(core_particles))
THEN
340 cpabort(
"Core particles should not be defined without corresponding shell particles.")
346 nparticle = n_atom + n_shell
347 ALLOCATE (scoord(3, nparticle), atype(nparticle))
349 shell_index = particles%els(i)%shell_index
350 IF (shell_index == 0)
THEN
351 CALL real_to_scaled(scoord(1:3, i), particles%els(i)%r(1:3), cell)
352 CALL get_atomic_kind(atomic_kind=particles%els(i)%atomic_kind, kind_number=atype(i))
354 CALL real_to_scaled(scoord(1:3, i), core_particles%els(shell_index)%r(1:3), cell)
355 CALL get_atomic_kind(atomic_kind=core_particles%els(shell_index)%atomic_kind, kind_number=atype(i))
356 k = n_atom + shell_index
357 CALL real_to_scaled(scoord(1:3, k), shell_particles%els(shell_index)%r(1:3), cell)
358 CALL get_atomic_kind(atomic_kind=shell_particles%els(shell_index)%atomic_kind, kind_number=atype(k))
366 IF (n_sr_rep > 0)
THEN
367 ALLOCATE (ranges(2, n_sr_rep))
370 ranges(:, i) = tmp(:)
372 CALL spgr_create(scoord, atype, cell, gopt_env, eps_symmetry=eps_symmetry, pol=pol(1:3), &
373 ranges=ranges, nparticle=nparticle, n_atom=n_atom, &
374 n_core=n_core, n_shell=n_shell, iunit=iunit, print_atoms=print_atoms)
377 CALL spgr_create(scoord, atype, cell, gopt_env, eps_symmetry=eps_symmetry, pol=pol(1:3), &
378 nparticle=nparticle, n_atom=n_atom, &
379 n_core=n_core, n_shell=n_shell, iunit=iunit, print_atoms=print_atoms)
383 spgr => gopt_env%spgr
386 CALL spgr_reduce_symm(spgr)
387 CALL spgr_rotations_subset(spgr)
389 DEALLOCATE (scoord, atype)
391 CALL timestop(handle)
405 TYPE(
spgr_type),
INTENT(INOUT),
POINTER :: spgr
406 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:, :), &
409 CHARACTER(LEN=*),
PARAMETER :: routinen =
'spgr_find_equivalent_atoms'
411 INTEGER :: handle, i, ia, ib, ir, j, natom, nop, &
413 REAL(kind=
dp) :: diff
414 REAL(kind=
dp),
DIMENSION(3) :: rb, ri, ro, tr
415 REAL(kind=
dp),
DIMENSION(3, 3) :: rot
417 CALL timeset(routinen, handle)
419 nop = spgr%n_operations
421 nshell = spgr%n_shell
423 IF (.NOT. (spgr%nparticle == (natom + nshell)))
THEN
424 cpabort(
"spgr_find_equivalent_atoms: nparticle not equal to natom + nshell.")
427 DO ia = 1, spgr%nparticle
428 spgr%eqatom(:, ia) = ia
433 IF (.NOT. spgr%lat(ia)) cycle
434 ri(1:3) = scoord(1:3, ia)
436 rot(1:3, 1:3) = spgr%rotations(1:3, 1:3, ir)
437 tr(1:3) = spgr%translations(1:3, ir)
439 IF (.NOT. spgr%lat(ib)) cycle
440 rb(1:3) = scoord(1:3, ib)
441 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)
442 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)
443 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)
444 ro(1) = ro(1) - real(nint(ro(1) - ri(1)),
dp)
445 ro(2) = ro(2) - real(nint(ro(2) - ri(2)),
dp)
446 ro(3) = ro(3) - real(nint(ro(3) - ri(3)),
dp)
447 diff = norm2(ri(:) - ro(:))
448 IF ((diff < spgr%eps_symmetry) .AND. (spgr%atype(ia) == spgr%atype(ib)))
THEN
449 spgr%eqatom(ir, ia) = ib
460 IF (.NOT. spgr%lat(ia)) cycle
461 ri(1:3) = scoord(1:3, ia)
463 rot(1:3, 1:3) = spgr%rotations(1:3, 1:3, ir)
464 tr(1:3) = spgr%translations(1:3, ir)
467 IF (.NOT. spgr%lat(ib)) cycle
468 rb(1:3) = scoord(1:3, ib)
469 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)
470 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)
471 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)
472 ro(1) = ro(1) - real(nint(ro(1) - ri(1)),
dp)
473 ro(2) = ro(2) - real(nint(ro(2) - ri(2)),
dp)
474 ro(3) = ro(3) - real(nint(ro(3) - ri(3)),
dp)
475 diff = norm2(ri(:) - ro(:))
476 IF ((diff < spgr%eps_symmetry) .AND. (spgr%atype(ia) == spgr%atype(ib)))
THEN
477 spgr%eqatom(ir, ia) = ib
485 CALL timestop(handle)
496 SUBROUTINE spgr_reduce_symm(spgr)
498 TYPE(
spgr_type),
INTENT(INOUT),
POINTER :: spgr
500 CHARACTER(LEN=*),
PARAMETER :: routinen =
'spgr_reduce_symm'
502 INTEGER :: handle, ia, ib, ir, ja, jb, nop, nops, &
504 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:) :: x, xold
505 REAL(kind=
dp),
DIMENSION(3) :: ri, ro
506 REAL(kind=
dp),
DIMENSION(3, 3) :: rot
508 CALL timeset(routinen, handle)
510 nop = spgr%n_operations
511 nparticle = spgr%nparticle
512 ALLOCATE (x(3*nparticle), xold(3*nparticle))
516 x(ja + 1) = x(ja + 1) + spgr%pol(1)
517 x(ja + 2) = x(ja + 2) + spgr%pol(2)
518 x(ja + 3) = x(ja + 3) + spgr%pol(3)
525 spgr%lop(ir) = .true.
526 rot(1:3, 1:3) = spgr%rotations(1:3, 1:3, ir)
528 IF (.NOT. spgr%lat(ia)) cycle
530 ri(1:3) = xold(ja + 1:ja + 3)
531 ro(1) = real(rot(1, 1),
dp)*ri(1) + real(rot(2, 1),
dp)*ri(2) + real(rot(3, 1),
dp)*ri(3)
532 ro(2) = real(rot(1, 2),
dp)*ri(1) + real(rot(2, 2),
dp)*ri(2) + real(rot(3, 2),
dp)*ri(3)
533 ro(3) = real(rot(1, 3),
dp)*ri(1) + real(rot(2, 3),
dp)*ri(2) + real(rot(3, 3),
dp)*ri(3)
534 x(ja + 1:ja + 3) = ro(1:3)
537 IF (.NOT. spgr%lat(ia)) cycle
538 ib = spgr%eqatom(ir, ia)
541 ro = x(jb + 1:jb + 3) - xold(ja + 1:ja + 3)
542 spgr%lop(ir) = (spgr%lop(ir) .AND. (abs(ro(1)) < spgr%eps_symmetry) &
543 .AND. (abs(ro(2)) < spgr%eps_symmetry) &
544 .AND. (abs(ro(3)) < spgr%eps_symmetry))
546 IF (spgr%lop(ir)) nops = nops + 1
549 spgr%n_reduced_operations = nops
552 CALL timestop(handle)
554 END SUBROUTINE spgr_reduce_symm
564 SUBROUTINE spgr_rotations_subset(spgr)
566 TYPE(
spgr_type),
INTENT(INOUT),
POINTER :: spgr
568 CHARACTER(LEN=*),
PARAMETER :: routinen =
'spgr_rotations_subset'
570 INTEGER :: handle, i, j
571 INTEGER,
DIMENSION(3, 3) :: d
572 LOGICAL,
ALLOCATABLE,
DIMENSION(:) :: mask
574 CALL timeset(routinen, handle)
576 ALLOCATE (mask(spgr%n_operations))
579 DO i = 1, spgr%n_operations
580 IF (.NOT. spgr%lop(i)) mask(i) = .false.
583 DO i = 1, spgr%n_operations - 1
584 IF (.NOT. mask(i)) cycle
585 DO j = i + 1, spgr%n_operations
586 IF (.NOT. mask(j)) cycle
587 d(:, :) = spgr%rotations(:, :, j) - spgr%rotations(:, :, i)
588 IF (sum(abs(d)) == 0) mask(j) = .false.
592 spgr%n_operations_subset = 0
593 DO i = 1, spgr%n_operations
594 IF (mask(i)) spgr%n_operations_subset = spgr%n_operations_subset + 1
597 ALLOCATE (spgr%rotations_subset(3, 3, spgr%n_operations_subset))
600 DO i = 1, spgr%n_operations
603 spgr%rotations_subset(:, :, j) = spgr%rotations(:, :, i)
608 CALL timestop(handle)
610 END SUBROUTINE spgr_rotations_subset
622 TYPE(
spgr_type),
INTENT(IN),
POINTER :: spgr
623 REAL(kind=
dp),
DIMENSION(:),
INTENT(INOUT) :: coord
625 CHARACTER(LEN=*),
PARAMETER :: routinen =
'spgr_apply_rotations_coord'
627 INTEGER :: handle, ia, ib, ir, ja, jb, nop, nops, &
629 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:) :: cold
630 REAL(kind=
dp),
DIMENSION(3) :: rf, ri, rn, ro, tr
631 REAL(kind=
dp),
DIMENSION(3, 3) :: rot
633 CALL timeset(routinen, handle)
635 ALLOCATE (cold(
SIZE(coord)))
638 nop = spgr%n_operations
639 nparticle = spgr%nparticle
640 nops = spgr%n_reduced_operations
644 IF (.NOT. spgr%lat(ia)) cycle
649 IF (.NOT. spgr%lop(ir)) cycle
650 ib = spgr%eqatom(ir, ia)
651 rot(1:3, 1:3) = spgr%rotations(1:3, 1:3, ir)
652 tr(1:3) = spgr%translations(1:3, ir)
655 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)
656 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)
657 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)
658 ro(1) = ro(1) - real(nint(ro(1) - rf(1)),
dp)
659 ro(2) = ro(2) - real(nint(ro(2) - rf(2)),
dp)
660 ro(3) = ro(3) - real(nint(ro(3) - rf(3)),
dp)
661 rn(1:3) = rn(1:3) + ro(1:3)
663 rn = rn/real(nops,
dp)
669 CALL timestop(handle)
683 TYPE(
spgr_type),
INTENT(IN),
POINTER :: spgr
684 REAL(kind=
dp),
DIMENSION(:),
INTENT(INOUT) :: force
686 CHARACTER(LEN=*),
PARAMETER :: routinen =
'spgr_apply_rotations_force'
688 INTEGER :: handle, ia, ib, ir, ja, jb, nop, nops, &
690 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:) :: fold
691 REAL(kind=
dp),
DIMENSION(3) :: ri, rn, ro
692 REAL(kind=
dp),
DIMENSION(3, 3) :: rot
694 CALL timeset(routinen, handle)
696 ALLOCATE (fold(
SIZE(force)))
699 nop = spgr%n_operations
700 nparticle = spgr%nparticle
701 nops = spgr%n_reduced_operations
705 IF (.NOT. spgr%lat(ia)) cycle
709 IF (.NOT. spgr%lop(ir)) cycle
710 ib = spgr%eqatom(ir, ia)
711 rot(1:3, 1:3) = spgr%rotations(1:3, 1:3, ir)
714 ro(1) = real(rot(1, 1),
dp)*ri(1) + real(rot(2, 1),
dp)*ri(2) + real(rot(3, 1),
dp)*ri(3)
715 ro(2) = real(rot(1, 2),
dp)*ri(1) + real(rot(2, 2),
dp)*ri(2) + real(rot(3, 2),
dp)*ri(3)
716 ro(3) = real(rot(1, 3),
dp)*ri(1) + real(rot(2, 3),
dp)*ri(2) + real(rot(3, 3),
dp)*ri(3)
717 rn(1:3) = rn(1:3) + ro(1:3)
719 rn = rn/real(nops,
dp)
725 CALL timestop(handle)
737 SUBROUTINE spgr_change_basis(roti, roto, nop, h1, h2)
739 INTEGER,
DIMENSION(:, :, :) :: roti
740 REAL(kind=
dp),
DIMENSION(:, :, :) :: roto
742 REAL(kind=
dp),
DIMENSION(3, 3) :: h1, h2
744 CHARACTER(LEN=*),
PARAMETER :: routinen =
'spgr_change_basis'
746 INTEGER :: handle, ir
747 REAL(kind=
dp),
DIMENSION(3, 3) :: h1ih2, h2ih1, ih1, ih2, r, s
749 CALL timeset(routinen, handle)
753 h2ih1 = matmul(h2, ih1)
754 h1ih2 = matmul(h1, ih2)
757 r(:, :) = roti(:, :, ir)
760 roto(:, :, ir) = r(:, :)
763 CALL timestop(handle)
765 END SUBROUTINE spgr_change_basis
778 TYPE(
spgr_type),
INTENT(IN),
POINTER :: spgr
779 TYPE(
cell_type),
INTENT(IN),
POINTER :: cell
780 REAL(kind=
dp),
DIMENSION(3, 3),
INTENT(INOUT) :: stress
782 CHARACTER(LEN=*),
PARAMETER :: routinen =
'spgr_apply_rotations_stress'
784 INTEGER :: handle, i, ir, j, k, l, nop
785 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:, :, :) :: roto
786 REAL(kind=
dp),
DIMENSION(3, 3) :: hmat1, hmat2, r, stin
788 CALL timeset(routinen, handle)
790 hmat1 = transpose(cell%hmat)
797 nop = spgr%n_operations_subset
799 ALLOCATE (roto(3, 3, nop))
801 CALL spgr_change_basis(spgr%rotations_subset, roto, spgr%n_operations_subset, hmat1, hmat2)
806 r(:, :) = roto(:, :, ir)
811 stress(i, j) = stress(i, j) + (r(k, i)*r(l, j)*stin(k, l))
817 stress = stress/real(nop,
dp)
821 CALL timestop(handle)
834 TYPE(
spgr_type),
INTENT(IN),
POINTER :: spgr
838 IF (spgr%iunit > 0)
THEN
839 WRITE (spgr%iunit,
'(/,T2,A,A)')
"----------------------------------------", &
840 "---------------------------------------"
841 WRITE (spgr%iunit,
"(T2,A,T25,A,T77,A)")
"----",
"SPACE GROUP SYMMETRY INFORMATION",
"----"
842 WRITE (spgr%iunit,
'(T2,A,A)')
"----------------------------------------", &
843 "---------------------------------------"
844 IF (spgr%symlib)
THEN
845 WRITE (spgr%iunit,
'(T2,A,T73,I8)')
"SPGR| SPACE GROUP NUMBER:", &
846 spgr%space_group_number
847 WRITE (spgr%iunit,
'(T2,A,T70,A11)')
"SPGR| INTERNATIONAL SYMBOL:", &
848 trim(adjustr(spgr%international_symbol))
849 WRITE (spgr%iunit,
'(T2,A,T75,A6)')
"SPGR| POINT GROUP SYMBOL:", &
850 trim(adjustr(spgr%pointgroup_symbol))
851 WRITE (spgr%iunit,
'(T2,A,T74,A7)')
"SPGR| SCHOENFLIES SYMBOL:", &
852 trim(adjustr(spgr%schoenflies))
853 WRITE (spgr%iunit,
'(T2,A,T73,I8)')
"SPGR| NUMBER OF SYMMETRY OPERATIONS:", &
855 WRITE (spgr%iunit,
'(T2,A,T73,I8)')
"SPGR| NUMBER OF UNIQUE ROTATIONS:", &
856 spgr%n_operations_subset
857 WRITE (spgr%iunit,
'(T2,A,T73,I8)')
"SPGR| NUMBER OF REDUCED SYMMETRY OPERATIONS:", &
858 spgr%n_reduced_operations
859 WRITE (spgr%iunit,
'(T2,A,T65,I8,I8)')
"SPGR| NUMBER OF PARTICLES AND SYMMETRIZED PARTICLES:", &
860 spgr%nparticle, spgr%nparticle_sym
861 WRITE (spgr%iunit,
'(T2,A,T65,I8,I8)')
"SPGR| NUMBER OF ATOMS AND SYMMETRIZED ATOMS:", &
862 spgr%n_atom, spgr%n_atom_sym
863 WRITE (spgr%iunit,
'(T2,A,T65,I8,I8)')
"SPGR| NUMBER OF CORES AND SYMMETRIZED CORES:", &
864 spgr%n_core, spgr%n_core_sym
865 WRITE (spgr%iunit,
'(T2,A,T65,I8,I8)')
"SPGR| NUMBER OF SHELLS AND SYMMETRIZED SHELLS:", &
866 spgr%n_shell, spgr%n_shell_sym
867 IF (spgr%print_atoms)
THEN
868 WRITE (spgr%iunit, *)
"SPGR| ACTIVE REDUCED SYMMETRY OPERATIONS:", spgr%lop
869 WRITE (spgr%iunit,
'(/,T2,A,A)')
"----------------------------------------", &
870 "---------------------------------------"
871 WRITE (spgr%iunit,
'(T2,A,T34,A,T77,A)')
"----",
"EQUIVALENT ATOMS",
"----"
872 WRITE (spgr%iunit,
'(T2,A,A)')
"----------------------------------------", &
873 "---------------------------------------"
874 DO i = 1, spgr%nparticle
875 DO j = 1, spgr%n_operations
876 WRITE (spgr%iunit,
'(T2,A,T52,I8,I8,I8)')
"SPGR| ATOM | SYMMETRY OPERATION | EQUIVALENT ATOM", &
877 i, j, spgr%eqatom(j, i)
880 WRITE (spgr%iunit,
'(T2,A,A)')
"----------------------------------------", &
881 "---------------------------------------"
882 DO i = 1, spgr%n_operations
883 WRITE (spgr%iunit,
'(T2,A,T46,i4,T51,3I10,/,T51,3I10,/,T51,3I10)') &
884 "SPGR| SYMMETRY OPERATION #:", i, (spgr%rotations(j, :, i), j=1, 3)
885 WRITE (spgr%iunit,
'(T51,3F10.5)') spgr%translations(:, i)
889 WRITE (spgr%iunit,
"(T2,A)")
"SPGLIB for Crystal Symmetry Information determination is not availale"
906 REAL(kind=
dp),
DIMENSION(3, 3),
INTENT(IN) :: stress
907 TYPE(
spgr_type),
INTENT(IN),
POINTER :: spgr
909 REAL(kind=
dp),
DIMENSION(3) :: eigval
910 REAL(kind=
dp),
DIMENSION(3, 3) :: eigvec, stress_tensor
912 stress_tensor(:, :) = stress(:, :)*
pascal*1.0e-9_dp
914 IF (spgr%iunit > 0)
THEN
915 WRITE (unit=spgr%iunit, fmt=
'(/,T2,A)') &
916 'SPGR STRESS| Symmetrized stress tensor [GPa]'
917 WRITE (unit=spgr%iunit, fmt=
'(T2,A,T19,3(19X,A1))') &
918 'SPGR STRESS|',
'x',
'y',
'z'
919 WRITE (unit=spgr%iunit, fmt=
'(T2,A,T26,3(1X,ES19.11))') &
920 'SPGR STRESS| x', stress_tensor(1, 1:3)
921 WRITE (unit=spgr%iunit, fmt=
'(T2,A,T26,3(1X,ES19.11))') &
922 'SPGR STRESS| y', stress_tensor(2, 1:3)
923 WRITE (unit=spgr%iunit, fmt=
'(T2,A,T26,3(1X,ES19.11))') &
924 'SPGR STRESS| z', stress_tensor(3, 1:3)
925 WRITE (unit=spgr%iunit, fmt=
'(T2,A,T66,ES20.11)') &
926 'SPGR STRESS| 1/3 Trace', (stress_tensor(1, 1) + &
927 stress_tensor(2, 2) + &
928 stress_tensor(3, 3))/3.0_dp
929 WRITE (unit=spgr%iunit, fmt=
'(T2,A,T66,ES20.11)') &
930 'SPGR STRESS| Determinant',
det_3x3(stress_tensor(1:3, 1), &
931 stress_tensor(1:3, 2), &
932 stress_tensor(1:3, 3))
934 eigvec(:, :) = 0.0_dp
935 CALL jacobi(stress_tensor, eigval, eigvec)
936 WRITE (unit=spgr%iunit, fmt=
'(/,T2,A)') &
937 'SPGR STRESS| Eigenvectors and eigenvalues of the symmetrized stress tensor [GPa]'
938 WRITE (unit=spgr%iunit, fmt=
'(T2,A,T19,3(1X,I19))') &
939 'SPGR STRESS|', 1, 2, 3
940 WRITE (unit=spgr%iunit, fmt=
'(T2,A,T26,3(1X,ES19.11))') &
941 'SPGR STRESS| Eigenvalues', eigval(1:3)
942 WRITE (unit=spgr%iunit, fmt=
'(T2,A,T26,3(1X,F19.12))') &
943 'SPGR STRESS| x', eigvec(1, 1:3)
944 WRITE (unit=spgr%iunit, fmt=
'(T2,A,T26,3(1X,F19.12))') &
945 'SPGR STRESS| y', eigvec(2, 1:3)
946 WRITE (unit=spgr%iunit, fmt=
'(T2,A,T26,3(1X,F19.12))') &
947 '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