Octopus
spglib_f08.F90
Go to the documentation of this file.
1module spglib_f08
2
3 use iso_c_binding, only: c_char, c_int, c_double, c_ptr, c_null_char, c_f_pointer, c_associated
4
5 implicit none
6
7 private
8
9 enum, bind(C)
10 enumerator :: SPGLIB_SUCCESS = 0
11 enumerator :: SPGERR_SPACEGROUP_SEARCH_FAILED
12 enumerator :: SPGERR_CELL_STANDARDIZATION_FAILED
13 enumerator :: SPGERR_SYMMETRY_OPERATION_SEARCH_FAILED
14 enumerator :: SPGERR_ATOMS_TOO_CLOSE
15 enumerator :: SPGERR_POINTGROUP_NOT_FOUND
19 enumerator :: spgerr_none
20 end enum
23 integer(c_int) :: number
24 character(len=11) :: international_short
25 character(len=20) :: international_full
26 character(len=32) :: international
27 character(len=7) :: schoenflies
28 integer(c_int) :: hall_number
29 character(len=17) :: hall_symbol
30 character(len=6) :: choice
31 character(len=6) :: pointgroup_international
32 character(len=4) :: pointgroup_schoenflies
33 integer(c_int) :: arithmetic_crystal_class_number
34 character(len=7) :: arithmetic_crystal_class_symbol
38 integer(c_int) :: spacegroup_number
39 integer(c_int) :: hall_number
40 character(len=11) :: international_symbol
41 character(len=17) :: hall_symbol
42 character(len=6) :: choice
43 real(c_double) :: transformation_matrix(3, 3)
44 real(c_double) :: origin_shift(3)
45 integer(c_int) :: n_operations
46 integer(c_int), allocatable :: rotations(:, :, :)
47 real(c_double), allocatable :: translations(:, :)
48 integer(c_int) :: n_atoms
49 integer(c_int), allocatable :: wyckoffs(:)
50 character(len=7), allocatable :: site_symmetry_symbols(:)
51 integer(c_int), allocatable :: equivalent_atoms(:) !Beware mapping refers to positions starting at 0
52 integer(c_int), allocatable :: crystallographic_orbits(:) !Beware mapping refers to positions starting at 0
53 real(c_double) :: primitive_lattice(3, 3)
54 integer(c_int), allocatable :: mapping_to_primitive(:) !Beware mapping refers to positions starting at 0
55 integer(c_int) :: n_std_atoms
56 real(c_double) :: std_lattice(3, 3)
57 integer(c_int), allocatable :: std_types(:)
58 real(c_double), allocatable :: std_positions(:, :)
59 real(c_double) :: std_rotation_matrix(3, 3)
60 integer(c_int), allocatable :: std_mapping_to_primitive(:) !Beware mapping refers to positions starting at 0
61 character(len=6) :: pointgroup_symbol
62 integer(kind(SPGLIB_SUCCESS)) :: spglib_error
65 interface
67 function spg_get_symmetry(rotation, translation, max_size, lattice, &
68 & position, types, num_atom, symprec) bind(c)
69 import c_int, c_double
70 integer(c_int), intent(out) :: rotation(3, 3, *)
71 real(c_double), intent(out) :: translation(3, *)
72 integer(c_int), intent(in), value :: max_size
73 real(c_double), intent(in) :: lattice(3, 3), position(3, *)
74 integer(c_int), intent(in) :: types(*)
75 integer(c_int), intent(in), value :: num_atom
76 real(c_double), intent(in), value :: symprec
77 integer(c_int) :: spg_get_symmetry
78 end function spg_get_symmetry
79
80 function spgat_get_symmetry(rotation, translation, max_size, lattice, &
81 & position, types, num_atom, symprec, angle_tolerance) bind(c)
82 import c_int, c_double
83 integer(c_int), intent(out) :: rotation(3, 3, *)
84 real(c_double), intent(out) :: translation(3, *)
85 integer(c_int), intent(in), value :: max_size
86 real(c_double), intent(in) :: lattice(3, 3), position(3, *)
87 integer(c_int), intent(in) :: types(*)
88 integer(c_int), intent(in), value :: num_atom
89 real(c_double), intent(in), value :: symprec, angle_tolerance
90 integer(c_int) :: spgat_get_symmetry
91 end function spgat_get_symmetry
92
93 function spg_get_symmetry_with_collinear_spin(rotation, translation, &
94 & equivalent_atoms, max_size, lattice, position, types, spins, &
95 & num_atom, symprec) bind(c)
96 import c_int, c_double
97 integer(c_int), intent(out) :: rotation(3, 3, *)
98 real(c_double), intent(out) :: translation(3, *)
99 integer(c_int), intent(out) :: equivalent_atoms(*)
100 integer(c_int), intent(in), value :: max_size
101 real(c_double), intent(in) :: lattice(3, 3), position(3, *)
102 integer(c_int), intent(in) :: types(*)
103 real(c_double), intent(in) :: spins(*)
104 integer(c_int), intent(in), value :: num_atom
105 real(c_double), intent(in), value :: symprec
108
109 function spgat_get_symmetry_with_collinear_spin(rotation, translation, &
110 & equivalent_atoms, max_size, lattice, position, types, spins, &
111 & num_atom, symprec, angle_tolerance) bind(c)
112 import c_int, c_double
113 integer(c_int), intent(out) :: rotation(3, 3, *)
114 real(c_double), intent(out) :: translation(3, *)
115 integer(c_int), intent(out) :: equivalent_atoms(*)
116 integer(c_int), intent(in), value :: max_size
117 real(c_double), intent(in) :: lattice(3, 3), position(3, *)
118 integer(c_int), intent(in) :: types(*)
119 real(c_double), intent(in) :: spins(*)
120 integer(c_int), intent(in), value :: num_atom
121 real(c_double), intent(in), value :: symprec, angle_tolerance
124
125 function spg_get_symmetry_with_site_tensors(rotation, translation, &
126 & equivalent_atoms, primitive_lattice, spin_flips, max_size, lattice, &
127 & position, types, tensors, tensor_rank, num_atom, with_time_reversal, &
128 & is_axial, symprec) bind(c)
129 import c_int, c_double
130 integer(c_int), intent(out) :: rotation(3, 3, *)
131 real(c_double), intent(out) :: translation(3, *)
132 integer(c_int), intent(out) :: equivalent_atoms(*)
133 real(c_double), intent(out) :: primitive_lattice(3, 3)
134 integer(c_int), intent(out) :: spin_flips(*)
135 integer(c_int), intent(in), value :: max_size
136 real(c_double), intent(in) :: lattice(3, 3), position(3, *)
137 integer(c_int), intent(in) :: types(*)
138 real(c_double), intent(in) :: tensors(*)
139 integer(c_int), intent(in), value :: tensor_rank
140 integer(c_int), intent(in), value :: num_atom
141 integer(c_int), intent(in), value :: with_time_reversal
142 integer(c_int), intent(in), value :: is_axial
143 real(c_double), intent(in), value :: symprec
144 integer(c_int) :: spg_get_symmetry_with_site_tensors
146
147 function spgat_get_symmetry_with_site_tensors(rotation, translation, &
148 & equivalent_atoms, primitive_lattice, spin_flips, max_size, lattice, &
149 & position, types, tensors, tensor_rank, num_atom, with_time_reversal, &
150 & is_axial, symprec, angle_tolerance) bind(c)
151 import c_int, c_double
152 integer(c_int), intent(out) :: rotation(3, 3, *)
153 real(c_double), intent(out) :: translation(3, *)
154 integer(c_int), intent(out) :: equivalent_atoms(*)
155 real(c_double), intent(out) :: primitive_lattice(3, 3)
156 integer(c_int), intent(out) :: spin_flips(*)
157 integer(c_int), intent(in), value :: max_size
158 real(c_double), intent(in) :: lattice(3, 3), position(3, *)
159 integer(c_int), intent(in) :: types(*)
160 real(c_double), intent(in) :: tensors(*)
161 integer(c_int), intent(in), value :: tensor_rank
162 integer(c_int), intent(in), value :: num_atom
163 integer(c_int), intent(in), value :: with_time_reversal
164 integer(c_int), intent(in), value :: is_axial
165 real(c_double), intent(in), value :: symprec, angle_tolerance
168
169 function spg_get_multiplicity(lattice, position, types, num_atom, &
170 & symprec) bind(c)
171 import c_int, c_double
172 real(c_double), intent(in) :: lattice(3, 3), position(3, *)
173 integer(c_int), intent(in) :: types(*)
174 integer(c_int), intent(in), value :: num_atom
175 real(c_double), intent(in), value :: symprec
176 integer(c_int) :: spg_get_multiplicity
177 end function spg_get_multiplicity
178
179 function spgat_get_multiplicity(lattice, position, types, num_atom, &
180 & symprec, angle_tolerance) bind(c)
181 import c_int, c_double
182 real(c_double), intent(in) :: lattice(3, 3), position(3, *)
183 integer(c_int), intent(in) :: types(*)
184 integer(c_int), intent(in), value :: num_atom
185 real(c_double), intent(in), value :: symprec, angle_tolerance
186 integer(c_int) :: spgat_get_multiplicity
187 end function spgat_get_multiplicity
188
189 function spg_delaunay_reduce(lattice, symprec) bind(c)
190 import c_int, c_double
191 real(c_double), intent(inout) :: lattice(3, 3)
192 real(c_double), intent(in), value :: symprec
193 integer(c_int) :: spg_delaunay_reduce
194 end function spg_delaunay_reduce
196 function spg_niggli_reduce(lattice, symprec) bind(c)
197 import c_int, c_double
198 real(c_double), intent(inout) :: lattice(3, 3)
199 real(c_double), intent(in), value :: symprec
200 integer(c_int) :: spg_niggli_reduce
201 end function spg_niggli_reduce
203 function spg_find_primitive(lattice, position, types, num_atom, symprec) &
204 & bind(c)
205 import c_int, c_double
206 real(c_double), intent(inout) :: lattice(3, 3), position(3, *)
207 integer(c_int), intent(inout) :: types(*)
208 integer(c_int), intent(in), value :: num_atom
209 real(c_double), intent(in), value :: symprec
210 integer(c_int) :: spg_find_primitive
211 end function spg_find_primitive
212
213 function spgat_find_primitive(lattice, position, types, num_atom, symprec, &
214 & angle_tolerance) bind(c)
215 import c_int, c_double
216 real(c_double), intent(inout) :: lattice(3, 3), position(3, *)
217 integer(c_int), intent(inout) :: types(*)
218 integer(c_int), intent(in), value :: num_atom
219 real(c_double), intent(in), value :: symprec, angle_tolerance
220 integer(c_int) :: spgat_find_primitive
221 end function spgat_find_primitive
222
223 function spg_get_international(symbol, lattice, position, types, num_atom, &
224 & symprec) bind(c)
225 import c_char, c_int, c_double
226 character(kind=c_char), intent(out) :: symbol(11)
227 real(c_double), intent(in) :: lattice(3, 3), position(3, *)
228 integer(c_int), intent(in) :: types(*)
229 integer(c_int), intent(in), value :: num_atom
230 real(c_double), intent(in), value :: symprec
231 integer(c_int) :: spg_get_international ! the number corresponding to 'symbol'. 0 on failure
232 end function spg_get_international
233
234 function spgat_get_international(symbol, lattice, position, types, num_atom, &
235 & symprec, angle_tolerance) bind(c)
236 import c_char, c_int, c_double
237 character(kind=c_char), intent(out) :: symbol(11)
238 real(c_double), intent(in) :: lattice(3, 3), position(3, *)
239 integer(c_int), intent(in) :: types(*)
240 integer(c_int), intent(in), value :: num_atom
241 real(c_double), intent(in), value :: symprec, angle_tolerance
242 integer(c_int) :: spgat_get_international ! the number corresponding to 'symbol'. 0 on failure
243 end function spgat_get_international
244
245 function spg_get_schoenflies(symbol, lattice, position, types, num_atom, &
246 & symprec) bind(c)
247 import c_char, c_int, c_double
248 character(kind=c_char), intent(out) :: symbol(7)
249 real(c_double), intent(in) :: lattice(3, 3), position(3, *)
250 integer(c_int), intent(in) :: types(*)
251 integer(c_int), intent(in), value :: num_atom
252 real(c_double), intent(in), value :: symprec
253 integer(c_int) :: spg_get_schoenflies ! the number corresponding to 'symbol'. 0 on failure
254 end function spg_get_schoenflies
255
256 function spgat_get_schoenflies(symbol, lattice, position, types, num_atom, &
257 & symprec, angle_tolerance) bind(c)
258 import c_char, c_int, c_double
259 character(kind=c_char), intent(out) :: symbol(7)
260 real(c_double), intent(in) :: lattice(3, 3), position(3, *)
261 integer(c_int), intent(in) :: types(*)
262 integer(c_int), intent(in), value :: num_atom
263 real(c_double), intent(in), value :: symprec, angle_tolerance
264 integer(c_int) :: spgat_get_schoenflies ! the number corresponding to 'symbol'. 0 on failure
265 end function spgat_get_schoenflies
266
267 function spg_get_pointgroup(symbol, trans_mat, rotations, num_rotations) &
268 & bind(c)
269 import c_char, c_int, c_double
270 character(kind=c_char) :: symbol(6)
271 integer(c_int), intent(out) :: trans_mat(3, 3)
272 integer(c_int), intent(in) :: rotations(3, 3, *)
273 integer(c_int), intent(in), value :: num_rotations
274 integer(c_int) :: spg_get_pointgroup
275 end function spg_get_pointgroup
276
277 function spg_refine_cell(lattice, position, types, num_atom, symprec) bind(c)
278 import c_int, c_double
279 real(c_double), intent(inout) :: lattice(3, 3), position(3, *)
280 integer(c_int), intent(inout) :: types(*)
281 integer(c_int), intent(in), value :: num_atom
282 real(c_double), intent(in), value :: symprec
283 integer(c_int) :: spg_refine_cell
284 end function spg_refine_cell
285
286 function spgat_refine_cell(lattice, position, types, num_atom, symprec, &
287 & angle_tolerance) bind(c)
288 import c_int, c_double
289 real(c_double), intent(inout) :: lattice(3, 3), position(3, *)
290 integer(c_int), intent(inout) :: types(*)
291 integer(c_int), intent(in), value :: num_atom
292 real(c_double), intent(in), value :: symprec, angle_tolerance
293 integer(c_int) :: spgat_refine_cell
294 end function spgat_refine_cell
295
296 function spg_standardize_cell(lattice, position, types, num_atom, &
297 & to_primitive, no_idealize, symprec) bind(c)
298 import c_int, c_double
299 real(c_double), intent(inout) :: lattice(3, 3), position(3, *)
300 integer(c_int), intent(inout) :: types(*)
301 integer(c_int), intent(in), value :: num_atom
302 integer(c_int), intent(in), value :: to_primitive, no_idealize
303 real(c_double), intent(in), value :: symprec
304 integer(c_int) :: spg_standardize_cell
305 end function spg_standardize_cell
306
307 function spgat_standardize_cell(lattice, position, types, num_atom, &
308 & to_primitive, no_idealize, symprec, angle_tolerance) bind(c)
309 import c_int, c_double
310 real(c_double), intent(inout) :: lattice(3, 3), position(3, *)
311 integer(c_int), intent(inout) :: types(*)
312 integer(c_int), intent(in), value :: num_atom
313 integer(c_int), intent(in), value :: to_primitive, no_idealize
314 real(c_double), intent(in), value :: symprec, angle_tolerance
315 integer(c_int) :: spgat_standardize_cell
316 end function spgat_standardize_cell
317
318 function spg_get_ir_reciprocal_mesh(grid_point, map, mesh, is_shift, &
319 & is_time_reversal, lattice, position, types, num_atom, symprec) bind(c)
320 import c_int, c_double
321 ! Beware the map refers to positions starting at 0
322 integer(c_int), intent(out) :: grid_point(3, *), map(*) ! size is product(mesh)
323 integer(c_int), intent(in) :: mesh(3), is_shift(3)
324 integer(c_int), intent(in), value :: is_time_reversal
325 real(c_double), intent(in) :: lattice(3, 3), position(3, *)
326 integer(c_int), intent(in) :: types(*)
327 integer(c_int), intent(in), value :: num_atom
328 real(c_double), intent(in), value :: symprec
329 integer(c_int) :: spg_get_ir_reciprocal_mesh ! the number of points in the reduced mesh
330 end function spg_get_ir_reciprocal_mesh
331
332 function spg_get_stabilized_reciprocal_mesh(grid_point, map, mesh, is_shift, &
333 & is_time_reversal, lattice, num_rot, rotations, num_q, qpoints) bind(c)
334 import c_int, c_double
335 ! Beware the map refers to positions starting at 0
336 integer(c_int), intent(out) :: grid_point(3, *), map(*)
337 integer(c_int), intent(in) :: mesh(3)
338 integer(c_int), intent(in) :: is_shift(3)
339 integer(c_int), intent(in), value :: is_time_reversal
340 real(c_double), intent(in) :: lattice(3, 3)
341 integer(c_int), intent(in), value :: num_rot
342 integer(c_int), intent(in) :: rotations(3, 3, *)
343 integer(c_int), intent(in), value :: num_q
344 real(c_double), intent(in) :: qpoints(3, *)
345 integer(c_int) :: spg_get_stabilized_reciprocal_mesh
347
348 function spg_get_error_code() bind(c, name='spg_get_error_code')
349 import spglib_success
350 integer(kind(SPGLIB_SUCCESS)) :: spg_get_error_code
351 end function spg_get_error_code
352
353 end interface
355 public :: spglibdataset, spg_get_dataset, &
372
373contains
374
375 function spg_get_error_message(spglib_error)
376 integer(kind(SPGLIB_SUCCESS)) :: spglib_error
377 character(len=32) :: spg_get_error_message
378
379 character, pointer, dimension(:) :: message
380 type(c_ptr) :: message_ptr
382 integer :: i
383
384 interface
385
386 function spg_get_error_message_c(spglib_error) &
387 & bind(c, name='spg_get_error_message')
388 import c_ptr, spglib_success
389
390 integer(kind(SPGLIB_SUCCESS)), value :: spglib_error
391 type(c_ptr) :: spg_get_error_message_c
392
393 end function spg_get_error_message_c
394
395 end interface
396
397 message_ptr = spg_get_error_message_c(spglib_error)
398 call c_f_pointer(message_ptr, message, [len(spg_get_error_message)])
399
401 do i = 1, len(spg_get_error_message)
402 if (message(i) == c_null_char) exit
403 spg_get_error_message(i:i) = message(i)
404 end do
405
406 end function spg_get_error_message
407
408 function spg_get_spacegroup_type(hall_number) result(spgtype)
409 integer(c_int), intent(in) :: hall_number
410 type(spglibspacegrouptype) :: spgtype
411
412 type, bind(c) :: spglibspacegrouptype_c
413 integer(c_int) :: number
414 character(kind=c_char) :: international_short(11)
415 character(kind=c_char) :: international_full(20)
416 character(kind=c_char) :: international(32)
417 character(kind=c_char) :: schoenflies(7)
418 integer(c_int) :: hall_number
419 character(kind=c_char) :: hall_symbol(17)
420 character(kind=c_char) :: choice(6)
421 character(kind=c_char) :: pointgroup_international(6)
422 character(kind=c_char) :: pointgroup_schoenflies(4)
423 integer(c_int) :: arithmetic_crystal_class_number
424 character(kind=c_char) :: arithmetic_crystal_class_symbol(7)
425 end type spglibspacegrouptype_c
426
427 interface
428
429 function spg_get_spacegroup_type_c(hall_number) &
430 & bind(c, name='spg_get_spacegroup_type')
431 import c_int, spglibspacegrouptype_c
432 integer(c_int), intent(in), value :: hall_number
433 type(spglibspacegrouptype_c) :: spg_get_spacegroup_type_c
434 end function spg_get_spacegroup_type_c
435
436 end interface
437
438 type(spglibspacegrouptype_c):: spgtype_c
439 integer :: i
440
441 spgtype_c = spg_get_spacegroup_type_c(hall_number)
442
443 spgtype%number = spgtype_c%number
444 spgtype%hall_number = spgtype_c%hall_number
445 spgtype%arithmetic_crystal_class_number = &
446 & spgtype_c%arithmetic_crystal_class_number
447
448 do i = 1, size(spgtype_c%international_short)
449 if (spgtype_c%international_short(i) == c_null_char) then
450 spgtype%international_short(i:) = ' '
451 exit
452 end if
453 spgtype%international_short(i:i) = spgtype_c%international_short(i)
454 end do
455
456 do i = 1, size(spgtype_c%international_full)
457 if (spgtype_c%international_full(i) == c_null_char) then
458 spgtype%international_full(i:) = ' '
459 exit
460 end if
461 spgtype%international_full(i:i) = spgtype_c%international_full(i)
462 end do
463
464 do i = 1, size(spgtype_c%international)
465 if (spgtype_c%international(i) == c_null_char) then
466 spgtype%international(i:) = ' '
467 exit
468 end if
469 spgtype%international(i:i) = spgtype_c%international(i)
470 end do
471
472 do i = 1, size(spgtype_c%schoenflies)
473 if (spgtype_c%schoenflies(i) == c_null_char) then
474 spgtype%schoenflies(i:) = ' '
475 exit
476 end if
477 spgtype%schoenflies(i:i) = spgtype_c%schoenflies(i)
478 end do
479
480 do i = 1, size(spgtype_c%hall_symbol)
481 if (spgtype_c%hall_symbol(i) == c_null_char) then
482 spgtype%hall_symbol(i:) = ' '
483 exit
484 end if
485 spgtype%hall_symbol(i:i) = spgtype_c%hall_symbol(i)
486 end do
487
488 do i = 1, size(spgtype_c%choice)
489 if (spgtype_c%choice(i) == c_null_char) then
490 spgtype%choice(i:) = ' '
491 exit
492 end if
493 spgtype%choice(i:i) = spgtype_c%choice(i)
494 end do
495
496 do i = 1, size(spgtype_c%pointgroup_international)
497 if (spgtype_c%pointgroup_international(i) == c_null_char) then
498 spgtype%pointgroup_international(i:) = ' '
499 exit
500 end if
501 spgtype%pointgroup_international(i:i) = &
502 & spgtype_c%pointgroup_international(i)
503 end do
504
505 do i = 1, size(spgtype_c%pointgroup_schoenflies)
506 if (spgtype_c%pointgroup_schoenflies(i) == c_null_char) then
507 spgtype%pointgroup_schoenflies(i:) = ' '
508 exit
509 end if
510 spgtype%pointgroup_schoenflies(i:i) = &
511 & spgtype_c%pointgroup_schoenflies(i)
512 end do
513
514 do i = 1, size(spgtype_c%arithmetic_crystal_class_symbol)
515 if (spgtype_c%arithmetic_crystal_class_symbol(i) == c_null_char) then
516 spgtype%arithmetic_crystal_class_symbol(i:) = ' '
517 exit
518 end if
519 spgtype%arithmetic_crystal_class_symbol(i:i) = &
520 & spgtype_c%arithmetic_crystal_class_symbol(i)
521 end do
522
523 end function spg_get_spacegroup_type
524
525 function spg_get_dataset(lattice, position, types, num_atom, symprec) result(dset)
526
527 real(c_double), intent(in) :: lattice(3, 3)
528 real(c_double), intent(in) :: position(3, *)
529 integer(c_int), intent(in) :: types(*)
530 integer(c_int), intent(in), value :: num_atom
531 real(c_double), intent(in), value :: symprec
532 type(spglibdataset) :: dset
533
534 type, bind(c) :: spglibdataset_c
535 integer(c_int) :: spacegroup_number
536 integer(c_int) :: hall_number
537 character(kind=c_char) :: international_symbol(11)
538 character(kind=c_char) :: hall_symbol(17)
539 character(kind=c_char) :: choice(6)
540 real(c_double) :: transformation_matrix(3, 3)
541 real(c_double) :: origin_shift(3)
542 integer(c_int) :: n_operations
543 type(c_ptr) :: rotations
544 type(c_ptr) :: translations
545 integer(c_int) :: n_atoms
546 type(c_ptr) :: wyckoffs
547 type(c_ptr) :: site_symmetry_symbols
548 type(c_ptr) :: equivalent_atoms
549 type(c_ptr) :: crystallographic_orbits
550 real(c_double) :: primitive_lattice(3, 3)
551 type(c_ptr) :: mapping_to_primitive
552 integer(c_int) :: n_std_atoms
553 real(c_double) :: std_lattice(3, 3)
554 type(c_ptr) :: std_types
555 type(c_ptr) :: std_positions
556 real(c_double) :: std_rotation_matrix(3, 3)
557 type(c_ptr) :: std_mapping_to_primitive
558 character(kind=c_char) :: pointgroup_symbol(6)
559 end type spglibdataset_c
560
561 interface
562 function spg_get_dataset_c(lattice, position, types, num_atom, symprec) &
563 & bind(c, name='spg_get_dataset')
564 import c_int, c_double, c_ptr
565 real(c_double), intent(in) :: lattice(3, 3)
566 real(c_double), intent(in) :: position(3, *)
567 integer(c_int), intent(in) :: types(*)
568 integer(c_int), intent(in), value :: num_atom
569 real(c_double), intent(in), value :: symprec
570 type(c_ptr) :: spg_get_dataset_c
571 end function spg_get_dataset_c
572
573 subroutine spg_free_dataset_c(dataset) bind(c, name='spg_free_dataset')
574 import spglibdataset_c
575 type(spglibdataset_c), intent(inout) :: dataset
576 end subroutine spg_free_dataset_c
577
578 end interface
579
580 type(spglibdataset_c), pointer :: dset_c
581 type(c_ptr) :: dataset_ptr_c
582 integer(c_int) :: n_operations, n_atoms, n_std_atoms
583 integer :: i
584 real(c_double), pointer :: translations(:, :), std_positions(:, :)
585 integer(c_int), pointer :: rotations(:, :, :), wyckoffs(:), equivalent_atoms(:)
586 integer(c_int), pointer :: crystallographic_orbits(:), std_types(:)
587
588 dataset_ptr_c = spg_get_dataset_c(lattice, position, types, num_atom, symprec)
589
590 if (c_associated(dataset_ptr_c)) then
591
592 dset%spglib_error = spglib_success
593
594 call c_f_pointer(dataset_ptr_c, dset_c)
595
596 dset%spacegroup_number = dset_c%spacegroup_number
597 dset%hall_number = dset_c%hall_number
598 dset%transformation_matrix = dset_c%transformation_matrix
599 dset%origin_shift = dset_c%origin_shift
600 dset%n_operations = dset_c%n_operations
601 dset%n_atoms = dset_c%n_atoms
602 dset%n_std_atoms = dset_c%n_std_atoms
603 dset%std_lattice = dset_c%std_lattice
604 dset%std_rotation_matrix = dset_c%std_rotation_matrix
605 dset%primitive_lattice = dset_c%primitive_lattice
606
607 ! Copy C strings to Fortran characters, converting C NULL to Fortran space padded strings
608 do i = 1, size(dset_c%international_symbol)
609 if (dset_c%international_symbol(i) == c_null_char) then
610 dset%international_symbol(i:) = ' '
611 exit
612 end if
613 dset%international_symbol(i:i) = dset_c%international_symbol(i)
614 end do
615
616 do i = 1, size(dset_c%hall_symbol)
617 if (dset_c%hall_symbol(i) == c_null_char) then
618 dset%hall_symbol(i:) = ' '
619 exit
620 end if
621 dset%hall_symbol(i:i) = dset_c%hall_symbol(i)
622 end do
623
624 do i = 1, size(dset_c%choice)
625 if (dset_c%choice(i) == c_null_char) then
626 dset%choice(i:) = ' '
627 exit
628 end if
629 dset%choice(i:i) = dset_c%choice(i)
630 end do
631
632 do i = 1, size(dset_c%pointgroup_symbol)
633 if (dset_c%pointgroup_symbol(i) == c_null_char) then
634 dset%pointgroup_symbol(i:) = ' '
635 exit
636 end if
637 dset%pointgroup_symbol(i:i) = dset_c%pointgroup_symbol(i)
638 end do
639
640 n_operations = dset_c%n_operations
641 n_atoms = dset_c%n_atoms
642 n_std_atoms = dset_c%n_std_atoms
643
644 call c_f_pointer(dset_c%rotations, rotations, shape=[3, 3, n_operations])
645 call c_f_pointer(dset_c%translations, translations, shape=[3, n_operations])
646 call c_f_pointer(dset_c%wyckoffs, wyckoffs, shape=[n_atoms])
647 call c_f_pointer(dset_c%equivalent_atoms, equivalent_atoms, shape=[n_atoms])
648 call c_f_pointer(dset_c%crystallographic_orbits, crystallographic_orbits, shape=[n_atoms])
649 call c_f_pointer(dset_c%std_types, std_types, shape=[n_std_atoms])
650 call c_f_pointer(dset_c%std_positions, std_positions, shape=[3, n_std_atoms])
651
652 allocate (dset%rotations(3, 3, n_operations))
653 allocate (dset%translations(3, n_operations))
654 allocate (dset%wyckoffs(n_atoms))
655 allocate (dset%equivalent_atoms(n_atoms))
656 allocate (dset%crystallographic_orbits(n_atoms))
657 allocate (dset%std_types(n_std_atoms))
658 allocate (dset%std_positions(3, n_std_atoms))
659
660 dset%rotations = rotations
661 dset%translations = translations
662 dset%wyckoffs = wyckoffs
663 dset%equivalent_atoms = equivalent_atoms
664 dset%crystallographic_orbits = crystallographic_orbits
665 dset%std_types = std_types
666 dset%std_positions = std_positions
667
668 call spg_free_dataset_c(dset_c)
669
670 else
671
672 dset%spglib_error = spg_get_error_code()
673
674 dset%spacegroup_number = 0
675 dset%hall_number = 0
676 dset%international_symbol = ' '
677 dset%hall_symbol = ' '
678 dset%choice = ' '
679 dset%transformation_matrix = 0.0_c_double
680 dset%origin_shift = 0.0_c_double
681 dset%n_operations = 0
682 dset%n_atoms = 0
683 dset%n_std_atoms = 0
684 dset%std_lattice = 0.0_c_double
685 dset%primitive_lattice = 0.0_c_double
686 dset%pointgroup_symbol = ' '
687
688 end if
689
690 end function spg_get_dataset
691
692end module spglib_f08
type(spglibdataset) function, public spg_get_dataset(lattice, position, types, num_atom, symprec)
Definition: spglib_f08.F90:532
type(spglibspacegrouptype) function, public spg_get_spacegroup_type(hall_number)
Definition: spglib_f08.F90:415
@ spgerr_niggli_failed
Definition: spglib_f08.F90:22
@ spgerr_delaunay_failed
Definition: spglib_f08.F90:23
@ spgerr_array_size_shortage
Definition: spglib_f08.F90:24
character(len=32) function, public spg_get_error_message(spglib_error)
Definition: spglib_f08.F90:382