(git:ccc2433)
spglib_f08.F
Go to the documentation of this file.
1 !--------------------------------------------------------------------------------------------------!
2 ! CP2K: A general program to perform molecular dynamics simulations !
3 ! Copyright 2000-2024 CP2K developers group <https://cp2k.org> !
4 ! !
5 ! SPDX-License-Identifier: GPL-2.0-or-later !
6 !--------------------------------------------------------------------------------------------------!
7 
8 ! **************************************************************************************************
9 !> \brief Interface for SPGLIB symmetry routines
10 !> \par History
11 !> \author jgh
12 ! **************************************************************************************************
13 #ifdef __SPGLIB
14 MODULE spglib_f08
15 
16  USE iso_c_binding, ONLY: c_char,&
17  c_double,&
18  c_int
19 
20  IMPLICIT NONE
21 
22  PRIVATE
23 
24  INTERFACE
25 
26  FUNCTION spg_get_symmetry(rotation, translation, max_size, lattice, &
27  position, types, num_atom, symprec) bind(c)
28  import c_int, c_double
29  INTEGER(c_int), INTENT(inout) :: rotation(3, 3, *)
30  REAL(c_double), INTENT(inout) :: translation(3, *)
31  INTEGER(c_int), INTENT(in), value :: max_size
32  REAL(c_double), INTENT(in) :: lattice(3, 3), position(3, *)
33  INTEGER(c_int), INTENT(in) :: types(*)
34  INTEGER(c_int), INTENT(in), value :: num_atom
35  REAL(c_double), INTENT(in), value :: symprec
36  INTEGER(c_int) :: spg_get_symmetry
37  END FUNCTION spg_get_symmetry
38 
39  FUNCTION spg_get_multiplicity(lattice, position, types, num_atom, symprec) bind(c)
40  import c_int, c_double
41  REAL(c_double), INTENT(in) :: lattice(3, 3), position(3, *)
42  INTEGER(c_int), INTENT(in) :: types(*)
43  INTEGER(c_int), INTENT(in), value :: num_atom
44  REAL(c_double), INTENT(in), value :: symprec
45  INTEGER(c_int) :: spg_get_multiplicity
46  END FUNCTION spg_get_multiplicity
47 
48  FUNCTION spg_get_international(symbol, lattice, position, types, num_atom, symprec) bind(c)
49  import c_char, c_int, c_double
50  CHARACTER(kind=c_char), INTENT(out) :: symbol(11)
51  REAL(c_double), INTENT(in) :: lattice(3, 3), position(3, *)
52  INTEGER(c_int), INTENT(in) :: types(*)
53  INTEGER(c_int), INTENT(in), value :: num_atom
54  REAL(c_double), INTENT(in), value :: symprec
55  INTEGER(c_int) :: spg_get_international ! the number corresponding to 'symbol'. 0 on failure
56  END FUNCTION spg_get_international
57 
58  FUNCTION spg_get_schoenflies(symbol, lattice, position, types, num_atom, symprec) bind(c)
59  import c_char, c_int, c_double
60  CHARACTER(kind=c_char), INTENT(out) :: symbol(7)
61  REAL(c_double), INTENT(in) :: lattice(3, 3), position(3, *)
62  INTEGER(c_int), INTENT(in) :: types(*)
63  INTEGER(c_int), INTENT(in), value :: num_atom
64  REAL(c_double), INTENT(in), value :: symprec
65  INTEGER(c_int) :: spg_get_schoenflies ! the number corresponding to 'symbol'. 0 on failure
66  END FUNCTION spg_get_schoenflies
67 
68  FUNCTION spg_get_pointgroup(symbol, trans_mat, rotations, num_rotations) bind(c)
69  import c_char, c_int, c_double
70  CHARACTER(kind=c_char) :: symbol(6)
71  INTEGER(c_int), INTENT(inout) :: trans_mat(3, 3)
72  INTEGER(c_int), INTENT(in) :: rotations(3, 3, *)
73  INTEGER(c_int), INTENT(in), value :: num_rotations
74  INTEGER(c_int) :: spg_get_pointgroup
75  END FUNCTION spg_get_pointgroup
76 
77  FUNCTION spg_get_ir_reciprocal_mesh(grid_point, map, mesh, &
78  is_shift, is_time_reversal, lattice, position, types, num_atom, symprec) bind(c)
79  import c_int, c_double
80 ! Beware the map refers to positions starting at 0
81  INTEGER(c_int), INTENT(out) :: grid_point(3, *), map(*) ! size is product(mesh)
82  INTEGER(c_int), INTENT(in) :: mesh(3), is_shift(3)
83  INTEGER(c_int), INTENT(in), value :: is_time_reversal
84  REAL(c_double), INTENT(in) :: lattice(3, 3), position(3, *)
85  INTEGER(c_int), INTENT(in) :: types(*)
86  INTEGER(c_int), INTENT(in), value :: num_atom
87  REAL(c_double), INTENT(in), value :: symprec
88  INTEGER(c_int) :: spg_get_ir_reciprocal_mesh ! the number of points in the reduced mesh
89  END FUNCTION spg_get_ir_reciprocal_mesh
90 
91  FUNCTION spg_get_major_version() bind(c)
92  import c_int
93  INTEGER(c_int) :: spg_get_major_version
94  END FUNCTION spg_get_major_version
95 
96  FUNCTION spg_get_minor_version() bind(c)
97  import c_int
98  INTEGER(c_int) :: spg_get_minor_version
99  END FUNCTION spg_get_minor_version
100 
101  FUNCTION spg_get_micro_version() bind(c)
102  import c_int
103  INTEGER(c_int) :: spg_get_micro_version
104  END FUNCTION spg_get_micro_version
105 
106  END INTERFACE
107 
111 
112 END MODULE spglib_f08
113 #else
114 ! **************************************************************************************************
115 !> \brief This is a stub for the Interface for SPGLIB symmetry routines
116 !> \par History
117 !> \author jgh
118 ! **************************************************************************************************
120 
121  USE kinds, ONLY: dp
122 
123 #include "./base/base_uses.f90"
124 
125  IMPLICIT NONE
126  PRIVATE
130 
131 CONTAINS
132 
133 ! **************************************************************************************************
134 !> \brief ...
135 !> \param rotation ...
136 !> \param translation ...
137 !> \param max_size ...
138 !> \param lattice ...
139 !> \param position ...
140 !> \param types ...
141 !> \param num_atom ...
142 !> \param symprec ...
143 !> \return ...
144 ! **************************************************************************************************
145  FUNCTION spg_get_symmetry(rotation, translation, max_size, lattice, &
146  position, types, num_atom, symprec)
147  INTEGER, INTENT(inout) :: rotation(:, :, :)
148  REAL(kind=dp), INTENT(inout) :: translation(:, :)
149  INTEGER, INTENT(in) :: max_size
150  REAL(kind=dp), INTENT(in) :: lattice(:, :), position(:, :)
151  INTEGER, INTENT(in) :: types(:), num_atom
152  REAL(kind=dp), INTENT(in) :: symprec
153  INTEGER :: spg_get_symmetry
154 
155  mark_used(rotation)
156  mark_used(translation)
157  mark_used(max_size)
158  mark_used(lattice)
159  mark_used(position)
160  mark_used(types)
161  mark_used(num_atom)
162  mark_used(symprec)
163  spg_get_symmetry = 0
164  cpabort("Requires linking to the SPGLIB library.")
165  END FUNCTION spg_get_symmetry
166 
167 ! **************************************************************************************************
168 !> \brief ...
169 !> \param lattice ...
170 !> \param position ...
171 !> \param types ...
172 !> \param num_atom ...
173 !> \param symprec ...
174 !> \return ...
175 ! **************************************************************************************************
176  FUNCTION spg_get_multiplicity(lattice, position, types, num_atom, symprec)
177  REAL(kind=dp), INTENT(in) :: lattice(:, :), position(:, :)
178  INTEGER, INTENT(in) :: types(:), num_atom
179  REAL(kind=dp), INTENT(in) :: symprec
180  INTEGER :: spg_get_multiplicity
181 
182  mark_used(lattice)
183  mark_used(position)
184  mark_used(types)
185  mark_used(num_atom)
186  mark_used(symprec)
188  cpabort("Requires linking to the SPGLIB library.")
189  END FUNCTION spg_get_multiplicity
190 
191 ! **************************************************************************************************
192 !> \brief ...
193 !> \param symbol ...
194 !> \param lattice ...
195 !> \param position ...
196 !> \param types ...
197 !> \param num_atom ...
198 !> \param symprec ...
199 !> \return ...
200 ! **************************************************************************************************
201  FUNCTION spg_get_international(symbol, lattice, position, types, num_atom, symprec)
202  CHARACTER, INTENT(out) :: symbol(11)
203  REAL(kind=dp), INTENT(in) :: lattice(:, :), position(:, :)
204  INTEGER, INTENT(in) :: types(:), num_atom
205  REAL(kind=dp), INTENT(in) :: symprec
206  INTEGER :: spg_get_international
207 
208 ! the number corresponding to 'symbol'. 0 on failure
209  mark_used(lattice)
210  mark_used(position)
211  mark_used(types)
212  mark_used(num_atom)
213  mark_used(symprec)
215  symbol = " "
216  cpabort("Requires linking to the SPGLIB library.")
217  END FUNCTION spg_get_international
218 
219 ! **************************************************************************************************
220 !> \brief ...
221 !> \param symbol ...
222 !> \param lattice ...
223 !> \param position ...
224 !> \param types ...
225 !> \param num_atom ...
226 !> \param symprec ...
227 !> \return ...
228 ! **************************************************************************************************
229  FUNCTION spg_get_schoenflies(symbol, lattice, position, types, num_atom, symprec)
230  CHARACTER, INTENT(out) :: symbol(7)
231  REAL(kind=dp), INTENT(in) :: lattice(:, :), position(:, :)
232  INTEGER, INTENT(in) :: types(:), num_atom
233  REAL(kind=dp), INTENT(in) :: symprec
234  INTEGER :: spg_get_schoenflies
235 
236 ! the number corresponding to 'symbol'. 0 on failure
237  mark_used(lattice)
238  mark_used(position)
239  mark_used(types)
240  mark_used(num_atom)
241  mark_used(symprec)
243  symbol = " "
244  cpabort("Requires linking to the SPGLIB library.")
245  END FUNCTION spg_get_schoenflies
246 
247 ! **************************************************************************************************
248 !> \brief ...
249 !> \param symbol ...
250 !> \param trans_mat ...
251 !> \param rotations ...
252 !> \param num_rotations ...
253 !> \return ...
254 ! **************************************************************************************************
255  FUNCTION spg_get_pointgroup(symbol, trans_mat, rotations, num_rotations)
256  CHARACTER :: symbol(6)
257  INTEGER, INTENT(inout) :: trans_mat(:, :)
258  INTEGER, INTENT(in) :: rotations(:, :, :), num_rotations
259  INTEGER :: spg_get_pointgroup
260 
261  mark_used(trans_mat)
262  mark_used(rotations)
263  mark_used(num_rotations)
265  symbol = " "
266  cpabort("Requires linking to the SPGLIB library.")
267  END FUNCTION spg_get_pointgroup
268 
269 ! **************************************************************************************************
270 !> \brief ...
271 !> \param grid_point ...
272 !> \param map ...
273 !> \param mesh ...
274 !> \param is_shift ...
275 !> \param is_time_reversal ...
276 !> \param lattice ...
277 !> \param position ...
278 !> \param types ...
279 !> \param num_atom ...
280 !> \param symprec ...
281 !> \return ...
282 ! **************************************************************************************************
283  FUNCTION spg_get_ir_reciprocal_mesh(grid_point, map, mesh, &
284  is_shift, is_time_reversal, lattice, position, types, num_atom, symprec)
285  INTEGER, INTENT(out) :: grid_point(:, :), map(:)
286  INTEGER, INTENT(in) :: mesh(:), is_shift(:), is_time_reversal
287  REAL(kind=dp), INTENT(in) :: lattice(:, :), position(:, :)
288  INTEGER, INTENT(in) :: types(:), num_atom
289  REAL(kind=dp), INTENT(in) :: symprec
290  INTEGER :: spg_get_ir_reciprocal_mesh
291 
292 ! size is product(mesh)
293 
294 ! the number of points in the reduced mesh
295  mark_used(grid_point)
296  mark_used(map)
297  mark_used(mesh)
298  mark_used(is_shift)
299  mark_used(is_time_reversal)
300  mark_used(lattice)
301  mark_used(position)
302  mark_used(types)
303  mark_used(num_atom)
304  mark_used(symprec)
306  cpabort("Requires linking to the SPGLIB library.")
307  END FUNCTION spg_get_ir_reciprocal_mesh
308 
309 ! **************************************************************************************************
310 !> \brief ...
311 !> \return ...
312 ! **************************************************************************************************
314  INTEGER :: spg_get_major_version
315 
317  END FUNCTION spg_get_major_version
318 
319 ! **************************************************************************************************
320 !> \brief ...
321 !> \return ...
322 ! **************************************************************************************************
324  INTEGER :: spg_get_minor_version
325 
327  END FUNCTION spg_get_minor_version
328 
329 ! **************************************************************************************************
330 !> \brief ...
331 !> \return ...
332 ! **************************************************************************************************
334  INTEGER :: spg_get_micro_version
335 
337  END FUNCTION spg_get_micro_version
338 
339 END MODULE spglib_f08
340 
341 #endif
Defines the basic variable types.
Definition: kinds.F:23
integer, parameter, public dp
Definition: kinds.F:34
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_micro_version()
...
Definition: spglib_f08.F:334
integer function, public spg_get_minor_version()
...
Definition: spglib_f08.F:324
integer function, public spg_get_ir_reciprocal_mesh(grid_point, map, mesh, is_shift, is_time_reversal, lattice, position, types, num_atom, symprec)
...
Definition: spglib_f08.F:285
integer function, public spg_get_major_version()
...
Definition: spglib_f08.F:314
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