(git:e887d2c)
Loading...
Searching...
No Matches
space_groups.F
Go to the documentation of this file.
1!--------------------------------------------------------------------------------------------------!
2! CP2K: A general program to perform molecular dynamics simulations !
3! Copyright 2000-2026 CP2K developers group <https://cp2k.org> !
4! !
5! SPDX-License-Identifier: GPL-2.0-or-later !
6!--------------------------------------------------------------------------------------------------!
7
8! **************************************************************************************************
9!> \brief Space Group Symmetry Module (version 1.0, January 16, 2020)
10!> \par History
11!> Pierre-André Cazade [pcazade] 01.2020 - University of Limerick
12!> \author Pierre-André Cazade (first version)
13! **************************************************************************************************
16 USE bibliography, ONLY: togo2018,&
17 cite_reference
18 USE cell_methods, ONLY: cell_create,&
20 USE cell_types, ONLY: cell_copy,&
21 cell_type,&
26 USE gopt_f_types, ONLY: gopt_f_type
35 USE kinds, ONLY: dp
36 USE mathlib, ONLY: det_3x3,&
37 inv_3x3,&
38 jacobi
40 USE physcon, ONLY: pascal
49#include "../base/base_uses.f90"
50
51 IMPLICIT NONE
52
53 PRIVATE
54
55 CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'space_groups'
56
60
61CONTAINS
62
63! **************************************************************************************************
64!> \brief routine creates the space group structure
65!> \param scoor ...
66!> \param types ...
67!> \param cell ...
68!> \param gopt_env ...
69!> \param eps_symmetry ...
70!> \param pol ...
71!> \param ranges ...
72!> \param nparticle ...
73!> \param n_atom ...
74!> \param n_core ...
75!> \param n_shell ...
76!> \param iunit ...
77!> \param print_atoms ...
78!> \par History
79!> 01.2020 created [pcazade]
80!> \author Pierre-André Cazade (first version)
81! **************************************************************************************************
82 SUBROUTINE spgr_create(scoor, types, cell, gopt_env, eps_symmetry, pol, ranges, &
83 nparticle, n_atom, n_core, n_shell, iunit, print_atoms)
84
85 REAL(kind=dp), DIMENSION(:, :), INTENT(IN) :: scoor
86 INTEGER, DIMENSION(:), INTENT(IN) :: types
87 TYPE(cell_type), INTENT(IN), POINTER :: cell
88 TYPE(gopt_f_type), INTENT(IN), POINTER :: gopt_env
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
95
96 CHARACTER(LEN=*), PARAMETER :: routinen = 'spgr_create'
97#ifdef __SPGLIB
98 CHARACTER(LEN=1000) :: buffer
99 INTEGER :: ierr, nchars, nop, tra_mat(3, 3)
100#endif
101 INTEGER :: handle, i, j, n_sr_rep
102 INTEGER, ALLOCATABLE, DIMENSION(:) :: tmp_types
103 LOGICAL :: spglib
104 REAL(kind=dp), ALLOCATABLE, DIMENSION(:, :) :: tmp_coor
105 TYPE(spgr_type), POINTER :: spgr
106
107 CALL timeset(routinen, handle)
108
109 spgr => gopt_env%spgr
110 cpassert(ASSOCIATED(spgr))
111
112 CALL cleanup_spgr_type(spgr)
113
114 !..total number of particles (atoms plus shells)
115 IF (PRESENT(nparticle)) THEN
116 cpassert(nparticle == SIZE(scoor, 2))
117 spgr%nparticle = nparticle
118 ELSE
119 spgr%nparticle = SIZE(scoor, 2)
120 END IF
121
122 IF (PRESENT(n_atom)) THEN
123 spgr%n_atom = n_atom
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
128 ELSE
129 spgr%n_atom = spgr%nparticle
130 END IF
131
132 IF (PRESENT(n_core)) THEN
133 spgr%n_core = n_core
134 ELSE IF (PRESENT(n_shell)) THEN
135 spgr%n_core = n_shell
136 END IF
137
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
142 END IF
143
144 IF (.NOT. (spgr%nparticle == (spgr%n_atom + spgr%n_shell))) THEN
145 cpabort("spgr_create: nparticle not equal to natom + nshell.")
146 END IF
147
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
152
153 spgr%iunit = iunit
154 spgr%print_atoms = print_atoms
155
156 ! accuracy for symmetry
157 IF (PRESENT(eps_symmetry)) THEN
158 spgr%eps_symmetry = eps_symmetry
159 END IF
160
161 ! vector to test reduced symmetry
162 IF (PRESENT(pol)) THEN
163 spgr%pol(1) = pol(1)
164 spgr%pol(2) = pol(2)
165 spgr%pol(3) = pol(3)
166 END IF
167
168 ALLOCATE (spgr%lat(spgr%nparticle))
169 spgr%lat = .true.
170
171 IF (PRESENT(ranges)) THEN
172 n_sr_rep = SIZE(ranges, 2)
173 DO i = 1, n_sr_rep
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
182 ELSE
183 cpabort("Symmetry exclusion range larger than actual number of particles.")
184 END IF
185 END DO
186 END DO
187 END IF
188
189 ALLOCATE (tmp_coor(3, spgr%n_atom_sym), tmp_types(spgr%n_atom_sym))
190
191 j = 0
192 DO i = 1, spgr%n_atom
193 IF (spgr%lat(i)) THEN
194 j = j + 1
195 tmp_coor(:, j) = scoor(:, i)
196 tmp_types(j) = types(i)
197 END IF
198 END DO
199
200 !..set cell values
201 NULLIFY (spgr%cell_ref)
202 CALL cell_create(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.")
213 CASE (default_cell_md_id)
214 cpabort("SPACE_GROUP_SYMMETRY is not compatible with md.")
215 CASE DEFAULT
216 cpabort("SPACE_GROUP_SYMMETRY invoked with an unknown optimization method.")
217 END SELECT
218 CASE DEFAULT
219 cpabort("SPACE_GROUP_SYMMETRY is not compatible with md.")
220 END SELECT
221
222 ! atom types
223 ALLOCATE (spgr%atype(spgr%nparticle))
224 spgr%atype(1:spgr%nparticle) = types(1:spgr%nparticle)
225
226 spgr%n_operations = 0
227
228#ifdef __SPGLIB
229 spglib = .true.
230 CALL cite_reference(togo2018)
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)
233 buffer = ''
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.")
238 spglib = .false.
239 ELSE
240 nop = spg_get_multiplicity(transpose(cell%hmat), tmp_coor, tmp_types, &
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
246 spgr%lop = .true.
247 ierr = spg_get_symmetry(spgr%rotations, spgr%translations, nop, transpose(cell%hmat), &
248 tmp_coor, tmp_types, spgr%n_atom_sym, eps_symmetry)
249 ! Schoenflies Symbol
250 ierr = spg_get_schoenflies(spgr%schoenflies, transpose(cell%hmat), tmp_coor, tmp_types, &
251 spgr%n_atom_sym, eps_symmetry)
252 buffer = ''
253 nchars = strlcpy_c2f(buffer, spgr%schoenflies)
254 spgr%schoenflies = buffer(1:nchars)
255
256 ! Point Group
257 tra_mat = 0
258 ierr = spg_get_pointgroup(spgr%pointgroup_symbol, tra_mat, &
259 spgr%rotations, spgr%n_operations)
260 buffer = ''
261 nchars = strlcpy_c2f(buffer, spgr%pointgroup_symbol)
262 spgr%pointgroup_symbol = buffer(1:nchars)
263 END IF
264#else
265 cpabort("Symmetry library SPGLIB not available")
266 spglib = .false.
267#endif
268 spgr%symlib = spglib
269
270 DEALLOCATE (tmp_coor, tmp_types)
271
272 CALL timestop(handle)
273
274 END SUBROUTINE spgr_create
275
276! **************************************************************************************************
277!> \brief routine indentifies the space group and finds rotation matrices.
278!> \param subsys ...
279!> \param geo_section ...
280!> \param gopt_env ...
281!> \param iunit ...
282!> \par History
283!> 01.2020 created [pcazade]
284!> \author Pierre-André Cazade (first version)
285!> \note rotation matrices innclude translations and translation symmetry:
286!> it works with supercells as well.
287! **************************************************************************************************
288 SUBROUTINE identify_space_group(subsys, geo_section, gopt_env, iunit)
289
290 TYPE(cp_subsys_type), INTENT(IN), POINTER :: subsys
291 TYPE(section_vals_type), INTENT(IN), POINTER :: geo_section
292 TYPE(gopt_f_type), INTENT(IN), POINTER :: gopt_env
293 INTEGER, INTENT(IN) :: iunit
294
295 CHARACTER(LEN=*), PARAMETER :: routinen = 'identify_space_group'
296
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
306 TYPE(cell_type), POINTER :: cell
307 TYPE(particle_list_type), POINTER :: core_particles, particles, &
308 shell_particles
309 TYPE(spgr_type), POINTER :: spgr
310
311 CALL timeset(routinen, handle)
312
313 n_sr_rep = 0
314 nparticle = 0
315 n_atom = 0
316 n_core = 0
317 n_shell = 0
318
319 NULLIFY (particles)
320 NULLIFY (core_particles)
321 NULLIFY (shell_particles)
322
323 NULLIFY (cell)
324 cell => subsys%cell
325 cpassert(ASSOCIATED(cell))
326
327 CALL cp_subsys_get(subsys, particles=particles, shell_particles=shell_particles, core_particles=core_particles)
328
329 cpassert(ASSOCIATED(particles))
330 n_atom = particles%n_els
331 ! Check if we have other kinds of particles in this subsystem
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
336 ! The same number of shell and core particles is assumed
337 cpassert(n_core == n_shell)
338 ELSE IF (ASSOCIATED(core_particles)) THEN
339 ! This case should not occur at the moment
340 cpabort("Core particles should not be defined without corresponding shell particles.")
341 ELSE
342 n_core = 0
343 n_shell = 0
344 END IF
345
346 nparticle = n_atom + n_shell
347 ALLOCATE (scoord(3, nparticle), atype(nparticle))
348 DO i = 1, n_atom
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))
353 ELSE
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))
359 END IF
360 END DO
361
362 CALL section_vals_val_get(geo_section, "SPGR_PRINT_ATOMS", l_val=print_atoms)
363 CALL section_vals_val_get(geo_section, "EPS_SYMMETRY", r_val=eps_symmetry)
364 CALL section_vals_val_get(geo_section, "SYMM_REDUCTION", r_vals=pol)
365 CALL section_vals_val_get(geo_section, "SYMM_EXCLUDE_RANGE", n_rep_val=n_sr_rep)
366 IF (n_sr_rep > 0) THEN
367 ALLOCATE (ranges(2, n_sr_rep))
368 DO i = 1, n_sr_rep
369 CALL section_vals_val_get(geo_section, "SYMM_EXCLUDE_RANGE", i_rep_val=i, i_vals=tmp)
370 ranges(:, i) = tmp(:)
371 END DO
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)
375 DEALLOCATE (ranges)
376 ELSE
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)
380 END IF
381
382 NULLIFY (spgr)
383 spgr => gopt_env%spgr
384
385 CALL spgr_find_equivalent_atoms(spgr, scoord)
386 CALL spgr_reduce_symm(spgr)
387 CALL spgr_rotations_subset(spgr)
388
389 DEALLOCATE (scoord, atype)
390
391 CALL timestop(handle)
392
393 END SUBROUTINE identify_space_group
394
395! **************************************************************************************************
396!> \brief routine indentifies the equivalent atoms for each rotation matrix.
397!> \param spgr ...
398!> \param scoord ...
399!> \par History
400!> 01.2020 created [pcazade]
401!> \author Pierre-André Cazade (first version)
402! **************************************************************************************************
403 SUBROUTINE spgr_find_equivalent_atoms(spgr, scoord)
404
405 TYPE(spgr_type), INTENT(INOUT), POINTER :: spgr
406 REAL(kind=dp), ALLOCATABLE, DIMENSION(:, :), &
407 INTENT(IN) :: scoord
408
409 CHARACTER(LEN=*), PARAMETER :: routinen = 'spgr_find_equivalent_atoms'
410
411 INTEGER :: handle, i, ia, ib, ir, j, natom, nop, &
412 nshell
413 REAL(kind=dp) :: diff
414 REAL(kind=dp), DIMENSION(3) :: rb, ri, ro, tr
415 REAL(kind=dp), DIMENSION(3, 3) :: rot
416
417 CALL timeset(routinen, handle)
418
419 nop = spgr%n_operations
420 natom = spgr%n_atom
421 nshell = spgr%n_shell
422
423 IF (.NOT. (spgr%nparticle == (natom + nshell))) THEN
424 cpabort("spgr_find_equivalent_atoms: nparticle not equal to natom + nshell.")
425 END IF
426
427 DO ia = 1, spgr%nparticle
428 spgr%eqatom(:, ia) = ia
429 END DO
430
431 !$OMP PARALLEL DO PRIVATE (ia,ib,ir,ri,rb,ro,rot,tr,diff) SHARED (spgr,scoord,natom,nop) DEFAULT(NONE)
432 DO ia = 1, natom
433 IF (.NOT. spgr%lat(ia)) cycle
434 ri(1:3) = scoord(1:3, ia)
435 DO ir = 1, nop
436 rot(1:3, 1:3) = spgr%rotations(1:3, 1:3, ir)
437 tr(1:3) = spgr%translations(1:3, ir)
438 DO ib = 1, natom
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
450 EXIT
451 END IF
452 END DO
453 END DO
454 END DO
455 !$OMP END PARALLEL DO
456
457 !$OMP PARALLEL DO PRIVATE (i,j,ia,ib,ir,ri,rb,ro,rot,tr,diff) SHARED (spgr,scoord,natom,nshell,nop) DEFAULT(NONE)
458 DO i = 1, nshell
459 ia = natom + i
460 IF (.NOT. spgr%lat(ia)) cycle
461 ri(1:3) = scoord(1:3, ia)
462 DO ir = 1, nop
463 rot(1:3, 1:3) = spgr%rotations(1:3, 1:3, ir)
464 tr(1:3) = spgr%translations(1:3, ir)
465 DO j = 1, nshell
466 ib = natom + j
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
478 EXIT
479 END IF
480 END DO
481 END DO
482 END DO
483 !$OMP END PARALLEL DO
484
485 CALL timestop(handle)
486
487 END SUBROUTINE spgr_find_equivalent_atoms
488
489! **************************************************************************************************
490!> \brief routine looks for operations compatible with efield
491!> \param spgr ...
492!> \par History
493!> 01.2020 created [pcazade]
494!> \author Pierre-André Cazade (first version)
495! **************************************************************************************************
496 SUBROUTINE spgr_reduce_symm(spgr)
497
498 TYPE(spgr_type), INTENT(INOUT), POINTER :: spgr
499
500 CHARACTER(LEN=*), PARAMETER :: routinen = 'spgr_reduce_symm'
501
502 INTEGER :: handle, ia, ib, ir, ja, jb, nop, nops, &
503 nparticle
504 REAL(kind=dp), ALLOCATABLE, DIMENSION(:) :: x, xold
505 REAL(kind=dp), DIMENSION(3) :: ri, ro
506 REAL(kind=dp), DIMENSION(3, 3) :: rot
507
508 CALL timeset(routinen, handle)
509
510 nop = spgr%n_operations
511 nparticle = spgr%nparticle
512 ALLOCATE (x(3*nparticle), xold(3*nparticle))
513 x = 0.0_dp
514 DO ia = 1, nparticle
515 ja = 3*(ia - 1)
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)
519 END DO
520 xold(:) = x(:)
521
522 nops = 0
523 DO ir = 1, nop
524 x = 0.d0
525 spgr%lop(ir) = .true.
526 rot(1:3, 1:3) = spgr%rotations(1:3, 1:3, ir)
527 DO ia = 1, nparticle
528 IF (.NOT. spgr%lat(ia)) cycle
529 ja = 3*(ia - 1)
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)
535 END DO
536 DO ia = 1, nparticle
537 IF (.NOT. spgr%lat(ia)) cycle
538 ib = spgr%eqatom(ir, ia)
539 ja = 3*(ia - 1)
540 jb = 3*(ib - 1)
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))
545 END DO
546 IF (spgr%lop(ir)) nops = nops + 1
547 END DO
548
549 spgr%n_reduced_operations = nops
550
551 DEALLOCATE (x, xold)
552 CALL timestop(handle)
553
554 END SUBROUTINE spgr_reduce_symm
555
556! **************************************************************************************************
557!> \brief routine looks for unique rotations
558!> \param spgr ...
559!> \par History
560!> 01.2020 created [pcazade]
561!> \author Pierre-André Cazade (first version)
562! **************************************************************************************************
563
564 SUBROUTINE spgr_rotations_subset(spgr)
565
566 TYPE(spgr_type), INTENT(INOUT), POINTER :: spgr
567
568 CHARACTER(LEN=*), PARAMETER :: routinen = 'spgr_rotations_subset'
569
570 INTEGER :: handle, i, j
571 INTEGER, DIMENSION(3, 3) :: d
572 LOGICAL, ALLOCATABLE, DIMENSION(:) :: mask
573
574 CALL timeset(routinen, handle)
575
576 ALLOCATE (mask(spgr%n_operations))
577 mask = .true.
578
579 DO i = 1, spgr%n_operations
580 IF (.NOT. spgr%lop(i)) mask(i) = .false.
581 END DO
582
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.
589 END DO
590 END DO
591
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
595 END DO
596
597 ALLOCATE (spgr%rotations_subset(3, 3, spgr%n_operations_subset))
598
599 j = 0
600 DO i = 1, spgr%n_operations
601 IF (mask(i)) THEN
602 j = j + 1
603 spgr%rotations_subset(:, :, j) = spgr%rotations(:, :, i)
604 END IF
605 END DO
606
607 DEALLOCATE (mask)
608 CALL timestop(handle)
609
610 END SUBROUTINE spgr_rotations_subset
611
612! **************************************************************************************************
613!> \brief routine applies the rotation matrices to the coordinates.
614!> \param spgr ...
615!> \param coord ...
616!> \par History
617!> 01.2020 created [pcazade]
618!> \author Pierre-André Cazade (first version)
619! **************************************************************************************************
620 SUBROUTINE spgr_apply_rotations_coord(spgr, coord)
621
622 TYPE(spgr_type), INTENT(IN), POINTER :: spgr
623 REAL(kind=dp), DIMENSION(:), INTENT(INOUT) :: coord
624
625 CHARACTER(LEN=*), PARAMETER :: routinen = 'spgr_apply_rotations_coord'
626
627 INTEGER :: handle, ia, ib, ir, ja, jb, nop, nops, &
628 nparticle
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
632
633 CALL timeset(routinen, handle)
634
635 ALLOCATE (cold(SIZE(coord)))
636 cold(:) = coord(:)
637
638 nop = spgr%n_operations
639 nparticle = spgr%nparticle
640 nops = spgr%n_reduced_operations
641
642 !$OMP PARALLEL DO PRIVATE (ia,ib,ja,jb,ir,ri,ro,rf,rn,rot,tr) SHARED (spgr,coord,nparticle,nop,nops) DEFAULT(NONE)
643 DO ia = 1, nparticle
644 IF (.NOT. spgr%lat(ia)) cycle
645 ja = 3*(ia - 1)
646 CALL real_to_scaled(rf(1:3), coord(ja + 1:ja + 3), spgr%cell_ref)
647 rn(1:3) = 0.d0
648 DO ir = 1, nop
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)
653 jb = 3*(ib - 1)
654 CALL real_to_scaled(ri(1:3), coord(jb + 1:jb + 3), spgr%cell_ref)
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)
662 END DO
663 rn = rn/real(nops, dp)
664 CALL scaled_to_real(coord(ja + 1:ja + 3), rn(1:3), spgr%cell_ref)
665 END DO
666 !$OMP END PARALLEL DO
667
668 DEALLOCATE (cold)
669 CALL timestop(handle)
670
671 END SUBROUTINE spgr_apply_rotations_coord
672
673! **************************************************************************************************
674!> \brief routine applies the rotation matrices to the forces.
675!> \param spgr ...
676!> \param force ...
677!> \par History
678!> 01.2020 created [pcazade]
679!> \author Pierre-André Cazade (first version)
680! **************************************************************************************************
681 SUBROUTINE spgr_apply_rotations_force(spgr, force)
682
683 TYPE(spgr_type), INTENT(IN), POINTER :: spgr
684 REAL(kind=dp), DIMENSION(:), INTENT(INOUT) :: force
685
686 CHARACTER(LEN=*), PARAMETER :: routinen = 'spgr_apply_rotations_force'
687
688 INTEGER :: handle, ia, ib, ir, ja, jb, nop, nops, &
689 nparticle
690 REAL(kind=dp), ALLOCATABLE, DIMENSION(:) :: fold
691 REAL(kind=dp), DIMENSION(3) :: ri, rn, ro
692 REAL(kind=dp), DIMENSION(3, 3) :: rot
693
694 CALL timeset(routinen, handle)
695
696 ALLOCATE (fold(SIZE(force)))
697 fold(:) = force(:)
698
699 nop = spgr%n_operations
700 nparticle = spgr%nparticle
701 nops = spgr%n_reduced_operations
702
703 !$OMP PARALLEL DO PRIVATE (ia,ib,ja,jb,ir,ri,ro,rn,rot) SHARED (spgr,force,nparticle,nop,nops) DEFAULT(NONE)
704 DO ia = 1, nparticle
705 IF (.NOT. spgr%lat(ia)) cycle
706 ja = 3*(ia - 1)
707 rn(1:3) = 0.d0
708 DO ir = 1, nop
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)
712 jb = 3*(ib - 1)
713 CALL real_to_scaled(ri(1:3), force(jb + 1:jb + 3), spgr%cell_ref)
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)
718 END DO
719 rn = rn/real(nops, dp)
720 CALL scaled_to_real(force(ja + 1:ja + 3), rn(1:3), spgr%cell_ref)
721 END DO
722 !$OMP END PARALLEL DO
723
724 DEALLOCATE (fold)
725 CALL timestop(handle)
726
727 END SUBROUTINE spgr_apply_rotations_force
728
729! **************************************************************************************************
730!> \brief ...
731!> \param roti ...
732!> \param roto ...
733!> \param nop ...
734!> \param h1 ...
735!> \param h2 ...
736! **************************************************************************************************
737 SUBROUTINE spgr_change_basis(roti, roto, nop, h1, h2)
738
739 INTEGER, DIMENSION(:, :, :) :: roti
740 REAL(kind=dp), DIMENSION(:, :, :) :: roto
741 INTEGER :: nop
742 REAL(kind=dp), DIMENSION(3, 3) :: h1, h2
743
744 CHARACTER(LEN=*), PARAMETER :: routinen = 'spgr_change_basis'
745
746 INTEGER :: handle, ir
747 REAL(kind=dp), DIMENSION(3, 3) :: h1ih2, h2ih1, ih1, ih2, r, s
748
749 CALL timeset(routinen, handle)
750
751 ih1 = inv_3x3(h1)
752 ih2 = inv_3x3(h2)
753 h2ih1 = matmul(h2, ih1)
754 h1ih2 = matmul(h1, ih2)
755
756 DO ir = 1, nop
757 r(:, :) = roti(:, :, ir)
758 s = matmul(h2ih1, r)
759 r = matmul(s, h1ih2)
760 roto(:, :, ir) = r(:, :)
761 END DO
762
763 CALL timestop(handle)
764
765 END SUBROUTINE spgr_change_basis
766
767! **************************************************************************************************
768!> \brief routine applies the rotation matrices to the stress tensor.
769!> \param spgr ...
770!> \param cell ...
771!> \param stress ...
772!> \par History
773!> 01.2020 created [pcazade]
774!> \author Pierre-André Cazade (first version)
775! **************************************************************************************************
776 SUBROUTINE spgr_apply_rotations_stress(spgr, cell, stress)
777
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
781
782 CHARACTER(LEN=*), PARAMETER :: routinen = 'spgr_apply_rotations_stress'
783
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
787
788 CALL timeset(routinen, handle)
789
790 hmat1 = transpose(cell%hmat)
791
792 hmat2 = 0d0
793 hmat2(1, 1) = 1.d0
794 hmat2(2, 2) = 1.d0
795 hmat2(3, 3) = 1.d0
796
797 nop = spgr%n_operations_subset
798
799 ALLOCATE (roto(3, 3, nop))
800
801 CALL spgr_change_basis(spgr%rotations_subset, roto, spgr%n_operations_subset, hmat1, hmat2)
802
803 stin = stress
804 stress = 0.d0
805 DO ir = 1, nop
806 r(:, :) = roto(:, :, ir)
807 DO i = 1, 3
808 DO j = 1, 3
809 DO k = 1, 3
810 DO l = 1, 3
811 stress(i, j) = stress(i, j) + (r(k, i)*r(l, j)*stin(k, l))
812 END DO
813 END DO
814 END DO
815 END DO
816 END DO
817 stress = stress/real(nop, dp)
818
819 DEALLOCATE (roto)
820
821 CALL timestop(handle)
822
823 END SUBROUTINE spgr_apply_rotations_stress
824
825! **************************************************************************************************
826!> \brief routine prints Space Group Information.
827!> \param spgr ...
828!> \par History
829!> 01.2020 created [pcazade]
830!> \author Pierre-André Cazade (first version)
831! **************************************************************************************************
832 SUBROUTINE print_spgr(spgr)
833
834 TYPE(spgr_type), INTENT(IN), POINTER :: spgr
835
836 INTEGER :: i, j
837
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:", &
854 spgr%n_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)
878 END DO
879 END DO
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)
886 END DO
887 END IF
888 ELSE
889 WRITE (spgr%iunit, "(T2,A)") "SPGLIB for Crystal Symmetry Information determination is not availale"
890 END IF
891 END IF
892
893 END SUBROUTINE print_spgr
894
895! **************************************************************************************************
896!> \brief Variable precision output of the symmetrized stress tensor
897!>
898!> \param stress tensor ...
899!> \param spgr ...
900!> \par History
901!> 07.2020 adapted to spgr [pcazade]
902!> \author MK (26.08.2010).
903! **************************************************************************************************
904 SUBROUTINE spgr_write_stress_tensor(stress, spgr)
905
906 REAL(kind=dp), DIMENSION(3, 3), INTENT(IN) :: stress
907 TYPE(spgr_type), INTENT(IN), POINTER :: spgr
908
909 REAL(kind=dp), DIMENSION(3) :: eigval
910 REAL(kind=dp), DIMENSION(3, 3) :: eigvec, stress_tensor
911
912 stress_tensor(:, :) = stress(:, :)*pascal*1.0e-9_dp
913
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))
933 eigval(:) = 0.0_dp
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)
948 END IF
949
950 END SUBROUTINE spgr_write_stress_tensor
951
952END MODULE space_groups
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.
Definition cell_types.F:15
subroutine, public scaled_to_real(r, s, cell)
Transform scaled cell coordinates real coordinates. r=h*s.
Definition cell_types.F:521
subroutine, public real_to_scaled(s, r, cell)
Transform real to scaled cell coordinates. s=h_inv*r.
Definition cell_types.F:491
subroutine, public cell_copy(cell_in, cell_out, tag)
Copy cell variable.
Definition cell_types.F:131
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
collects all constants needed in input so that they can be used without circular dependencies
integer, parameter, public default_cell_geo_opt_id
integer, parameter, public default_cell_method_id
integer, parameter, public default_minimization_method_id
integer, parameter, public default_ts_method_id
integer, parameter, public default_cell_md_id
integer, parameter, public default_cell_direct_id
objects that represent the structure of input sections and the data contained in an input section
subroutine, public section_vals_val_get(section_vals, keyword_name, i_rep_section, i_rep_val, n_rep_val, val, l_val, i_val, r_val, c_val, l_vals, i_vals, r_vals, c_vals, explicit)
returns the requested value
Defines the basic variable types.
Definition kinds.F:23
integer, parameter, public dp
Definition kinds.F:34
Collection of simple mathematical functions and subroutines.
Definition mathlib.F:15
subroutine, public jacobi(a, d, v)
Jacobi matrix diagonalization. The eigenvalues are returned in vector d and the eigenvectors are retu...
Definition mathlib.F:1476
pure real(kind=dp) function, dimension(3, 3), public inv_3x3(a)
Returns the inverse of the 3 x 3 matrix a.
Definition mathlib.F:524
represent a simple array based list of the given type
Definition of physical constants:
Definition physcon.F:68
real(kind=dp), parameter, public pascal
Definition physcon.F:174
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.
Definition spglib_f08.F:119
integer function, public spg_get_international(symbol, lattice, position, types, num_atom, symprec)
...
Definition spglib_f08.F:202
integer function, public spg_get_multiplicity(lattice, position, types, num_atom, symprec)
...
Definition spglib_f08.F:177
integer function, public spg_get_pointgroup(symbol, trans_mat, rotations, num_rotations)
...
Definition spglib_f08.F:256
integer function, public spg_get_symmetry(rotation, translation, max_size, lattice, position, types, num_atom, symprec)
...
Definition spglib_f08.F:147
integer function, public spg_get_schoenflies(symbol, lattice, position, types, num_atom, symprec)
...
Definition spglib_f08.F:230
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.
Definition cell_types.F:60
represents a system: atoms, molecules, their pos,vel,...
calculates the potential energy of a system, and its derivatives