1.. _definitions_and_conventions:
2
3Definitions and conventions
4============================
5
6Information in this page is valid for spglib 1.8.1 or later. The
7definitions of transformation matrix and origin shift were different
8in the previous versions.
9
10.. contents::
11   :depth: 2
12   :local:
13
14
15References
16-----------
17
18Some references about crystallographic definitions and conventions are
19shown below. Though spglib may not follow them fully, it doesn't mean
20spglib doesn't respect them, rather it is due to spglib-author's lack of
21understanding the crystallography ashamedly.
22
23* `International Tables for Crystallography <http://it.iucr.org/>`_.
24* `Bilbao Crystallographic Server <http://www.cryst.ehu.es/>`_. The
25  references of many useful papers are found at
26  http://www.cryst.ehu.es/wiki/index.php/Articles.
27* Ulrich Müller, "Symmetry Relationships between Crystal Structures"
28* E. Parthé, K. Cenzual, and R. E. Gladyshevskii, "Standardization of
29  crystal structure data as an aid to the classification of crystal
30  structure types", Journal of Alloys and Compounds, **197**, 291-301
31  (1993). [`doi2
32  <https://dx.doi.org/10.1016/0925-8388(93)90049-S>`_]
33* E. Parthé and L. M. Gelato, "The ’best’ unit cell for monoclinic
34  structures consistent with b axis unique and cell choice 1
35  of international tables for crystallography (1983)", Acta
36  Cryst. A **41**, 142-151 (1985) [`doi3
37  <https://doi.org/10.1107/S0108767385000289>`_]
38* E. Parthé and L. M. Gelato, "The standardization of inorganic
39  crystal-structure data", Acta Cryst. A
40  **40**, 169-183 (1984) [`doi4
41  <https://doi.org/10.1107/S0108767384000416>`_]
42* S. Hall, "Space-group notation with an explicit origin", Acta
43  Cryst. A **37**, 517-525 (1981) [`doi1
44  <https://doi.org/10.1107/S0567739481001228>`_]
45
46Space group operation and change of basis
47------------------------------------------
48
49Basis vectors :math:`(\mathbf{a}, \mathbf{b}, \mathbf{c})` or :math:`(\mathbf{a}_1, \mathbf{a}_2, \mathbf{a}_3)`
50^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
51
52In spglib, basis vectors are represented by three column vectors:
53
54.. math::
55   :label: basis_vectors
56
57   \mathbf{a}= \begin{pmatrix}
58   a_x \\
59   a_y \\
60   a_z \\
61   \end{pmatrix},
62   \mathbf{b}= \begin{pmatrix}
63   b_x \\
64   b_y \\
65   b_z \\
66   \end{pmatrix},
67   \mathbf{c}= \begin{pmatrix}
68   c_x \\
69   c_y \\
70   c_z \\
71   \end{pmatrix},
72
73in Cartesian coordinates. Depending on the situation,
74:math:`(\mathbf{a}_1, \mathbf{a}_2, \mathbf{a}_3)` is used instead of
75:math:`(\mathbf{a}, \mathbf{b}, \mathbf{c})`.
76
77Atomic point coordinates :math:`\boldsymbol{x}`
78^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
79
80Coordinates of an atomic point :math:`\boldsymbol{x}` are represented
81as three fractional values relative to basis vectors as follows,
82
83.. math::
84   :label: coordinate_triplet
85
86   \boldsymbol{x}= \begin{pmatrix}
87   x_1 \\
88   x_2 \\
89   x_3 \\
90   \end{pmatrix},
91
92where :math:`0 \le x_i < 1`. A position vector :math:`\mathbf{x}` in
93Cartesian coordinates is obtained by
94
95.. math::
96   :label: position_vector_1
97
98   \mathbf{x} = (\mathbf{a}, \mathbf{b}, \mathbf{c}) \boldsymbol{x}.
99
100or
101
102.. math::
103   :label: position_vector_2
104
105   \mathbf{x} = \sum_i x_i \mathbf{a}_i.
106
107Symmetry operation :math:`(\boldsymbol{W}, \boldsymbol{w})`
108^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
109
110A symmetry operation consists of a pair of the rotation part
111:math:`\boldsymbol{W}` and translation part :math:`\boldsymbol{w}`,
112and is represented as :math:`(\boldsymbol{W}, \boldsymbol{w})` in the
113spglib document. The symmetry operation transfers :math:`\boldsymbol{x}` to
114:math:`\tilde{\boldsymbol{x}}` as follows:
115
116.. math::
117   :label: space_group_operation
118
119   \tilde{\boldsymbol{x}} = \boldsymbol{W}\boldsymbol{x} + \boldsymbol{w}.
120
121.. _def_transformation_and_origin_shift:
122
123Transformation matrix :math:`\boldsymbol{P}` and origin shift :math:`\boldsymbol{p}`
124^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
125
126The transformation matrix :math:`\boldsymbol{P}` changes choice of
127basis vectors as follows
128
129.. math::
130   :label: transformation_matrix
131
132   ( \mathbf{a} \; \mathbf{b} \; \mathbf{c} )
133   = ( \mathbf{a}_\mathrm{s} \; \mathbf{b}_\mathrm{s} \;
134   \mathbf{c}_\mathrm{s} )  \boldsymbol{P},
135
136where :math:`( \mathbf{a} \; \mathbf{b} \; \mathbf{c} )` and :math:`(
137\mathbf{a}_\mathrm{s} \; \mathbf{b}_\mathrm{s} \;
138\mathbf{c}_\mathrm{s} )` are the basis vectors of an arbitrary system
139and of a starndardized system, respectively. In general, the
140transformation matrix is not limited for the transformation from the
141standardized system, but can be used in between any systems possibly
142transformed. It has to be emphasized that the transformation matrix
143**doesn't** rotate a crystal in Cartesian coordinates, but just
144changes the choices of basis vectors.
145
146The origin shift :math:`\boldsymbol{p}` gives the vector from the
147origin of the standardized system :math:`\boldsymbol{O}_\mathrm{s}` to
148the origin of the arbitrary system :math:`\boldsymbol{O}`,
149
150.. math::
151   :label: origin_shift
152
153   \boldsymbol{p} = \boldsymbol{O} - \boldsymbol{O}_\mathrm{s}.
154
155Origin shift **doesn't** move a crystal in Cartesian coordinates, but
156just changes the origin to measure the coordinates of atomic points.
157
158
159A change of basis is described by the combination of the
160transformation matrix and the origin shift denoted by
161:math:`(\boldsymbol{P}, \boldsymbol{p})` where first the
162transformation matrix is applied and then origin shift. The points in
163the standardized system :math:`\boldsymbol{x}_\mathrm{s}` and
164arbitrary system :math:`\boldsymbol{x}` are related by
165
166.. math::
167   :label: change_of_basis_1
168
169   \boldsymbol{x}_\mathrm{s} = \boldsymbol{P}\boldsymbol{x} +
170   \boldsymbol{p},
171
172or equivalently,
173
174.. math::
175   :label: change_of_basis_2
176
177   \boldsymbol{x} = \boldsymbol{P}^{-1}\boldsymbol{x}_\mathrm{s} -
178   \boldsymbol{P}^{-1}\boldsymbol{p}.
179
180
181A graphical example is shown below.
182
183.. |cob| image:: change-of-basis.png
184         :width: 20%
185
186|cob|
187
188(click the figure to enlarge)
189
190In this example,
191
192.. math::
193
194   \boldsymbol{P} = \begin{pmatrix}
195   \frac{1}{2} & \frac{1}{2} & 0 \\
196   \frac{\bar{1}}{2} & \frac{1}{2} & 0 \\
197   0 & 0 & 1
198   \end{pmatrix}.
199
200Difference between rotation and transformation matrices
201^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
202
203A rotation matrix rotates (or mirrors, inverts) the crystal body with
204respect to origin. A transformation matrix changes the choice of the
205basis vectors, but does not rotate the crystal body.
206
207A space group operation having no translation part sends an atom to
208another point by
209
210.. math::
211
212   \tilde{\boldsymbol{x}} = \boldsymbol{W}\boldsymbol{x},
213
214where :math:`\tilde{\boldsymbol{x}}` and :math:`\boldsymbol{x}` are
215represented with respect to the same basis vectors :math:`(\mathbf{a},
216\mathbf{b}, \mathbf{c})`. Equivalently the rotation is achieved by
217rotating the basis vectors:
218
219.. math::
220   :label: rotation_matrix
221
222   (\tilde{\mathbf{a}}, \tilde{\mathbf{b}}, \tilde{\mathbf{c}}) =
223   (\mathbf{a}, \mathbf{b}, \mathbf{c}) \boldsymbol{W}
224
225with keeping the representation of the atomic point coordinates
226:math:`\boldsymbol{x}` because
227
228.. math::
229
230   \tilde{\mathbf{x}} = (\mathbf{a}, \mathbf{b}, \mathbf{c}) \tilde{\boldsymbol{x}}
231   = (\mathbf{a}, \mathbf{b}, \mathbf{c}) \boldsymbol{W}
232   \boldsymbol{x}
233   = (\tilde{\mathbf{a}}, \tilde{\mathbf{b}}, \tilde{\mathbf{c}})
234   \boldsymbol{x}.
235
236The transformation matrix changes the choice of the basis vectors as:
237
238.. math::
239
240   (\mathbf{a}', \mathbf{b}', \mathbf{c}') = (\mathbf{a}, \mathbf{b},
241   \mathbf{c}) \boldsymbol{P}.
242
243The atomic position vector is not altered by this transformation, i.e.,
244
245.. math::
246
247   (\mathbf{a}', \mathbf{b}', \mathbf{c}') \boldsymbol{x}' =
248   (\mathbf{a}, \mathbf{b}, \mathbf{c}) \boldsymbol{x}.
249
250However the representation of the atomic point coordinates changes as follows:
251
252.. math::
253
254   \boldsymbol{P} \boldsymbol{x}' = \boldsymbol{x}.
255
256because
257
258.. math::
259
260   (\mathbf{a}, \mathbf{b}, \mathbf{c}) \boldsymbol{P} \boldsymbol{x}'
261   = (\mathbf{a}', \mathbf{b}', \mathbf{c}') \boldsymbol{x}' =
262   (\mathbf{a}, \mathbf{b}, \mathbf{c}) \boldsymbol{x}.
263
264
265.. _def_standardized_unit_cell:
266
267Spglib conventions of standardized unit cell
268---------------------------------------------
269
270The standardization in spglib is achieved by :ref:`a change of basis
271transformation <def_transformation_and_origin_shift>`. If
272:ref:`idealization <def_idealize_cell>` as shown below is further
273applied, the crystal can be rigidly rotated in Cartesian
274coordinates.
275
276Choice of basis vectors
277^^^^^^^^^^^^^^^^^^^^^^^^
278
279Using the APIs ``spg_get_dataset``,
280``spg_get_dataset_with_hall_number``, or ``spg_standardize_cell``, the
281starndardized unit cell is obtained. The "starndardized unit cell" in
282this document means that the (conventional) unit cell structure is
283standardized by the crystal symmetry and lengths of basis
284vectors. This standardization in spglib is not unique, but upto space
285group operations and generators of Euclidean normalizer. Crystals are
286categorized by Hall symbols in 530 different types in terms of 230
287space group types, unique axes, settings, and cell choices. Moreover
288in spglib, lengths of basis vectors are used to choose the order of
289:math:`(\mathbf{a}, \mathbf{b}, \mathbf{c})` if the order can not be
290determined only by the symmetrical conventions.
291
292.. _def_standardized_primitive_cell:
293
294Transformation to the primitive cell
295^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
296
297In the standardized unit cells, there are five different centring
298types available, base centrings of A and C, rhombohedral (R), body centred
299(I), and face centred (F). The transformation is applied to the
300standardized unit cell by
301
302.. math::
303   :label: transformation_to_primitive
304
305   ( \mathbf{a}_\mathrm{p} \; \mathbf{b}_\mathrm{p} \; \mathbf{c}_\mathrm{p} )
306   = ( \mathbf{a}_\mathrm{s} \; \mathbf{b}_\mathrm{s} \;
307   \mathbf{c}_\mathrm{s} )  \boldsymbol{P}_\mathrm{c},
308
309where :math:`\mathbf{a}_\mathrm{p}`, :math:`\mathbf{b}_\mathrm{p}`,
310and :math:`\mathbf{c}_\mathrm{p}` are the basis vectors of the
311primitive cell and :math:`\boldsymbol{P}_\mathrm{c}` is the
312transformation matrix from the standardized unit cell to the primitive
313cell. :math:`\boldsymbol{P}_\mathrm{c}` for centring types are given
314as follows:
315
316.. math::
317
318   \boldsymbol{P}_\mathrm{A} =
319   \begin{pmatrix}
320   1 & 0 & 0 \\
321   0 & \frac{1}{2} & \frac{\bar{1}}{2} \\
322   0 & \frac{1}{2} & \frac{{1}}{2}
323   \end{pmatrix},
324   \boldsymbol{P}_\mathrm{C} =
325   \begin{pmatrix}
326   \frac{1}{2} & \frac{{1}}{2} & 0 \\
327   \frac{\bar{1}}{2} & \frac{1}{2} & 0\\
328   0 & 0 & 1
329   \end{pmatrix},
330   \boldsymbol{P}_\mathrm{R} =
331   \begin{pmatrix}
332   \frac{2}{3} & \frac{\bar{1}}{3} & \frac{\bar{1}}{3} \\
333   \frac{1}{3} & \frac{{1}}{3} & \frac{\bar{2}}{3} \\
334   \frac{1}{3} & \frac{{1}}{3} & \frac{{1}}{3}
335   \end{pmatrix},
336   \boldsymbol{P}_\mathrm{I} =
337   \begin{pmatrix}
338   \frac{\bar{1}}{2} & \frac{{1}}{2} & \frac{{1}}{2} \\
339   \frac{{1}}{2} & \frac{\bar{1}}{2} & \frac{{1}}{2} \\
340   \frac{{1}}{2} & \frac{{1}}{2} & \frac{\bar{1}}{2}
341   \end{pmatrix},
342   \boldsymbol{P}_\mathrm{F} =
343   \begin{pmatrix}
344   0 & \frac{{1}}{2} & \frac{{1}}{2} \\
345   \frac{{1}}{2} & 0 & \frac{{1}}{2} \\
346   \frac{{1}}{2} & \frac{{1}}{2} & 0
347   \end{pmatrix}.
348
349The choice of transformation matrix depends on purposes. These
350transformation matrices above are just the spglib choices and may not
351be the best.
352
353For rhombohedral lattice systems with the H setting (hexagonal
354lattice), :math:`\boldsymbol{P}_\mathrm{R}` is applied to obtain
355primitive basis vectors, but for that with the R setting (rhombohedral
356lattice), no transformation matrix is applied because it is already
357the primitive cell.
358
359.. _def_idealize_cell:
360
361Idealization of unit cell structure
362^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
363
364Spglib allows tolerance parameters to match a slightly distorted unit
365cell structure to a space group type with some higher symmetry. Using
366obtained symmetry operations, the distortion is removed to idealize
367the unit cell structure. The coordinates of atomic points are
368idealized using respective site-symmetries (Grosse-Kunstleve *et
369al*. (2002)). The basis vectors are idealized by forceing them into
370respective lattice shapes as follows. In this treatment, except for
371triclinic crystals, crystals can be rotated in Cartesian coordinates,
372which is the different type of transformation from that of the
373change-of-basis transformation explained above.
374
375Triclinic lattice
376""""""""""""""""""
377
378- Niggli reduced cell is used for choosing :math:`\mathbf{a}, \mathbf{b}, \mathbf{c}`.
379- :math:`\mathbf{a}` is set along :math:`+x` direction of Cartesian coordinates.
380- :math:`\mathbf{b}` is set in :math:`x\text{-}y` plane of Cartesian
381  coordinates so that :math:`\mathbf{a}\times\mathbf{b}` is along
382  :math:`+z` direction of Cartesian coordinates.
383
384Monoclinic lattice
385"""""""""""""""""""
386
387- :math:`b` axis is taken as the unique axis.
388- :math:`\alpha = 90^\circ` and :math:`\gamma = 90^\circ`
389- :math:`90^\circ < \beta < 120^\circ`.
390
391- :math:`\mathbf{a}` is set along :math:`+x` direction of Cartesian coordinates.
392- :math:`\mathbf{b}` is set along :math:`+y` direction of Cartesian coordinates.
393- :math:`\mathbf{c}` is set in :math:`x\text{-}z` plane of Cartesian coordinates.
394
395Orthorhombic lattice
396"""""""""""""""""""""
397
398- :math:`\alpha = \beta = \gamma = 90^\circ`.
399
400- :math:`\mathbf{a}` is set along :math:`+x` direction of Cartesian coordinates.
401- :math:`\mathbf{b}` is set along :math:`+y` direction of Cartesian coordinates.
402- :math:`\mathbf{c}` is set along :math:`+z` direction of Cartesian coordinates.
403
404Tetragonal lattice
405"""""""""""""""""""
406
407- :math:`\alpha = \beta = \gamma = 90^\circ`.
408- :math:`a=b`.
409
410- :math:`\mathbf{a}` is set along :math:`+x` direction of Cartesian coordinates.
411- :math:`\mathbf{b}` is set along :math:`+y` direction of Cartesian coordinates.
412- :math:`\mathbf{c}` is set along :math:`+z` direction of Cartesian coordinates.
413
414Rhombohedral lattice
415"""""""""""""""""""""
416
417- :math:`\alpha = \beta = \gamma`.
418- :math:`a=b=c`.
419
420- Let :math:`\mathbf{a}`, :math:`\mathbf{b}`, and :math:`\mathbf{c}`
421  projected on :math:`x\text{-}y` plane in Cartesian coordinates be
422  :math:`\mathbf{a}_{xy}`, :math:`\mathbf{b}_{xy}`, and
423  :math:`\mathbf{c}_{xy}`, respectively, and their angles be
424  :math:`\alpha_{xy}`, :math:`\beta_{xy}`,
425  :math:`\gamma_{xy}`, respectively.
426- Let :math:`\mathbf{a}`, :math:`\mathbf{b}`, and :math:`\mathbf{c}`
427  projected along :math:`z`-axis in Cartesian coordinates be
428  :math:`\mathbf{a}_{z}`, :math:`\mathbf{b}_{z}`, and
429  :math:`\mathbf{c}_{z}`, respectively.
430
431- :math:`\mathbf{a}_{xy}` is set along the ray :math:`30^\circ`
432  rotated counter-clockwise from the :math:`+x`
433  direction of Cartesian coordinates, and :math:`\mathbf{b}_{xy}` and
434  :math:`\mathbf{c}_{xy}` are placed by angles :math:`120^\circ` and
435  :math:`240^\circ` from :math:`\mathbf{a}_{xy}` counter-clockwise,
436  respectively.
437- :math:`\alpha_{xy} = \beta_{xy} = \gamma_{xy} = 120^\circ`.
438- :math:`a_{xy} = b_{xy} = c_{xy}`.
439- :math:`a_{z} = b_{z} = c_{z}`.
440
441
442Hexagonal lattice
443""""""""""""""""""
444
445- :math:`\alpha = \beta = 90^\circ`.
446- :math:`\gamma = 120^\circ`.
447- :math:`a=b`.
448
449- :math:`\mathbf{a}` is set along :math:`+x` direction of Cartesian coordinates.
450- :math:`\mathbf{b}` is set in :math:`x\text{-}y` plane of Cartesian coordinates.
451- :math:`\mathbf{c}` is set along :math:`+z` direction of Cartesian coordinates.
452
453Cubic lattice
454""""""""""""""
455
456- :math:`\alpha = \beta = \gamma = 90^\circ`.
457- :math:`a=b=c`.
458
459- :math:`\mathbf{a}` is set along :math:`+x` direction of Cartesian coordinates.
460- :math:`\mathbf{b}` is set along :math:`+y` direction of Cartesian coordinates.
461- :math:`\mathbf{c}` is set along :math:`+z` direction of Cartesian coordinates.
462
463Rotation introduced by idealization
464^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
465
466In the idealization step presented above, the input unit cell crystal
467strcuture can be rotated in the Cartesian coordinates.  The rotation
468matrix :math:`\boldsymbol{R}` of this rotation is defined by
469
470.. math::
471   :label: rigid_rotation_matrix
472
473   ( \bar{\mathbf{a}}_\mathrm{s} \;
474   \bar{\mathbf{b}}_\mathrm{s} \; \bar{\mathbf{c}}_\mathrm{s} )
475   = ( \boldsymbol{R} \mathbf{a}_\mathrm{s} \;
476   \boldsymbol{R} \mathbf{b}_\mathrm{s} \; \boldsymbol{R}
477   \mathbf{c}_\mathrm{s} ).
478
479This rotation matrix rotates the standardized
480crystal structure before idealization :math:`( \mathbf{a}_\mathrm{s}
481\; \mathbf{b}_\mathrm{s} \; \mathbf{c}_\mathrm{s} )` to that after
482idealization :math:`( \bar{\mathbf{a}}_\mathrm{s} \;
483\bar{\mathbf{b}}_\mathrm{s} \; \bar{\mathbf{c}}_\mathrm{s} )` in
484Cartesian coordinates of the given input unit cell.
485
486Examples
487--------
488
489Crystallographic choice and rigid rotation
490^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
491
492The following example of a python script gives a crystal structure of
493Br whose space group type is *Cmce*. The basis vectors
494:math:`(\mathbf{a}, \mathbf{b}, \mathbf{c})` are fixed by the symmetry
495crystal in the standardization. The C-centrng determines the c-axis,
496and *m* and *c* operations in *Cmce* fix which directions a- and
497b-axes should be with respect to each other axis. This is the first
498one choice appearing in the list of Hall symbols among 6 different
499choices for this space group type.
500
501::
502
503   import spglib
504
505   # Mind that the a, b, c axes are given in row vectors here,
506   # but the formulation above is given for the column vectors.
507   lattice = [[7.17851431, 0, 0],  # a
508              [0, 3.99943947, 0],  # b
509              [0, 0, 8.57154746]]  # c
510   points = [[0.0, 0.84688439, 0.1203133],
511             [0.0, 0.65311561, 0.6203133],
512             [0.0, 0.34688439, 0.3796867],
513             [0.0, 0.15311561, 0.8796867],
514             [0.5, 0.34688439, 0.1203133],
515             [0.5, 0.15311561, 0.6203133],
516             [0.5, 0.84688439, 0.3796867],
517             [0.5, 0.65311561, 0.8796867]]
518   numbers = [35,] * len(points)
519   cell = (lattice, points, numbers)
520   dataset = spglib.get_symmetry_dataset(cell)
521   print("Space group type: %s (%d)"
522         % (dataset['international'], dataset['number']))
523   print("Transformation matrix:")
524   for x in dataset['transformation_matrix']:
525       print("  %2d %2d %2d" % tuple(x))
526   print("Origin shift: %f %f %f" % tuple(dataset['origin_shift']))
527
528This python script is saved in the file ``example.py``. Then we get
529
530::
531
532   % python example.py
533   Space group type: Cmce (64)
534   Transformation matrix:
535      1  0  0
536      0  1  0
537      0  0  1
538   Origin shift: 0.000000 0.500000 0.500000
539
540No rotation was introduced in the idealization. Next, we swap a- and c-axes.
541
542::
543
544   import spglib
545
546   # Mind that the a, b, c axes are given in row vectors here,
547   # but the formulation above is given for the column vectors.
548   lattice = [[8.57154746, 0, 0],  # a
549              [0, 3.99943947, 0],  # b
550              [0, 0, 7.17851431]]  # c
551   points = [[0.1203133, 0.84688439, 0.0],
552             [0.6203133, 0.65311561, 0.0],
553             [0.3796867, 0.34688439, 0.0],
554             [0.8796867, 0.15311561, 0.0],
555             [0.1203133, 0.34688439, 0.5],
556             [0.6203133, 0.15311561, 0.5],
557             [0.3796867, 0.84688439, 0.5],
558             [0.8796867, 0.65311561, 0.5]]
559   numbers = [35,] * len(points)
560   cell = (lattice, points, numbers)
561   dataset = spglib.get_symmetry_dataset(cell)
562   print("Space group type: %s (%d)"
563         % (dataset['international'], dataset['number']))
564   print("Transformation matrix:")
565   for x in dataset['transformation_matrix']:
566       print("  %2d %2d %2d" % tuple(x))
567   print("Origin shift: %f %f %f" % tuple(dataset['origin_shift']))
568
569By this,
570
571::
572
573   % python spglib-example2.py
574   Space group type: Cmce (64)
575   Transformation matrix:
576      0  0  1
577      0  1  0
578     -1  0  0
579   Origin shift: 0.000000 0.000000 0.000000
580
581We get a non-identity transformation matrix, which wants to transform
582back to the original (above) crystal structure by swapping a- and
583c-axes. The transformation back of the basis vectors is achieved by
584Eq. :eq:`transformation_matrix`. Next, we try to rotate rigidly the
585crystal structure by :math:`45^\circ` around c-axis in Cartesian
586coordinates from the first one::
587
588   import spglib
589
590   # Mind that the a, b, c axes are given in row vectors here,
591   # but the formulation above is given for the column vectors.
592   lattice = [[5.0759761474456697, 5.0759761474456697, 0],  # a
593              [-2.8280307701821314, 2.8280307701821314, 0],  # b
594              [0, 0, 8.57154746]]  # c
595   points = [[0.0, 0.84688439, 0.1203133],
596             [0.0, 0.65311561, 0.6203133],
597             [0.0, 0.34688439, 0.3796867],
598             [0.0, 0.15311561, 0.8796867],
599             [0.5, 0.34688439, 0.1203133],
600             [0.5, 0.15311561, 0.6203133],
601             [0.5, 0.84688439, 0.3796867],
602             [0.5, 0.65311561, 0.8796867]]
603   numbers = [35,] * len(points)
604   cell = (lattice, points, numbers)
605   dataset = spglib.get_symmetry_dataset(cell)
606   print("Space group type: %s (%d)"
607         % (dataset['international'], dataset['number']))
608   print("Transformation matrix:")
609   for x in dataset['transformation_matrix']:
610       print("  %2d %2d %2d" % tuple(x))
611   print("Origin shift: %f %f %f" % tuple(dataset['origin_shift']))
612
613and
614
615::
616
617   % python spglib-example3.py
618   Space group type: Cmce (64)
619   Transformation matrix:
620      1  0  0
621      0  1  0
622      0  0  1
623   Origin shift: 0.000000 0.000000 0.500000
624
625The transformation matrix is kept unchanged even though the crystal
626structure is rotated in Cartesian coordinates. The origin shift is
627different but it changes only the order of atoms, so effectively it
628does nothing.
629
630Transformation to a primitive cell
631^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
632
633There are infinite number of choices of primitive cell. The
634transformation from a primitive cell basis vectors to the other
635primitive cell basis vectors is always done by an integer matrix
636because any lattice points can be generated by the linear combination
637of the three primitive basis vectors.
638
639When we have a non-primitive cell basis vectors as given in the above
640example::
641
642   lattice = [[7.17851431, 0, 0],  # a
643              [0, 3.99943947, 0],  # b
644              [0, 0, 8.57154746]]  # c
645
646This has the C-centring, so it must be transformed to a primitive
647cell. A possible transformation is shown at
648:ref:`def_standardized_primitive_cell`, which is
649:math:`\boldsymbol{P}_\mathrm{C}`. With the following script::
650
651   import numpy as np
652   lattice = [[7.17851431, 0, 0],  # a
653              [0, 3.99943947, 0],  # b
654              [0, 0, 8.57154746]]  # c
655   Pc = [[0.5, 0.5, 0],
656         [-0.5, 0.5, 0],
657         [0, 0, 1]]
658   print(np.dot(np.transpose(lattice), Pc).T)  # given in row vectors
659
660we get the primitive cell basis vectors (shown in row vectors)::
661
662   [[ 3.58925715 -1.99971973  0.        ]
663    [ 3.58925715  1.99971973  0.        ]
664    [ 0.          0.          8.57154746]]
665
666``find_primitive`` gives a primitive cell that is obtained by
667transforming standardized and idealized crystal structure to the
668primitive cell using the transformation matrix. Therefore by this
669script::
670
671   import spglib
672
673   lattice = [[7.17851431, 0, 0],
674              [0, 3.99943947, 0],
675              [0, 0, 8.57154746]]
676   points = [[0.0, 0.84688439, 0.1203133],
677             [0.0, 0.65311561, 0.6203133],
678             [0.0, 0.34688439, 0.3796867],
679             [0.0, 0.15311561, 0.8796867],
680             [0.5, 0.34688439, 0.1203133],
681             [0.5, 0.15311561, 0.6203133],
682             [0.5, 0.84688439, 0.3796867],
683             [0.5, 0.65311561, 0.8796867]]
684   numbers = [8,] * len(points)
685   cell = (lattice, points, numbers)
686
687   primitive_cell = spglib.find_primitive(cell)
688   print(primitive_cell[0])
689
690we get::
691
692   [[ 3.58925715 -1.99971973  0.        ]
693    [ 3.58925715  1.99971973  0.        ]
694    [ 0.          0.          8.57154746]]
695
696This is same as what we manually obtained above.
697Even when the basis vectors are rigidly rotated as::
698
699   lattice = [[5.0759761474456697, 5.0759761474456697, 0],
700              [-2.8280307701821314, 2.8280307701821314, 0],
701              [0, 0, 8.57154746]]
702
703the relationship of a, b, c axes is unchanged. Therefore the same
704transformation matrix to the primitive cell can be used. Then we get::
705
706   [[3.95200346 1.12397269 0.        ]
707    [1.12397269 3.95200346 0.        ]
708    [0.         0.         8.57154746]]
709
710However applying ``find_primitive`` rigidly rotates automatically and
711so the following script doesn't give this basis vectors::
712
713   import spglib
714
715   lattice = [[5.0759761474456697, 5.0759761474456697, 0],
716              [-2.8280307701821314, 2.8280307701821314, 0],
717              [0, 0, 8.57154746]]
718   points = [[0.0, 0.84688439, 0.1203133],
719             [0.0, 0.65311561, 0.6203133],
720             [0.0, 0.34688439, 0.3796867],
721             [0.0, 0.15311561, 0.8796867],
722             [0.5, 0.34688439, 0.1203133],
723             [0.5, 0.15311561, 0.6203133],
724             [0.5, 0.84688439, 0.3796867],
725             [0.5, 0.65311561, 0.8796867]]
726   numbers = [8,] * len(points)
727   cell = (lattice, points, numbers)
728
729   primitive_cell = spglib.find_primitive(cell)
730   print(primitive_cell[0])
731
732but gives those with respect to the idealized ones::
733
734   [[ 3.58925715 -1.99971973  0.        ]
735    [ 3.58925715  1.99971973  0.        ]
736    [ 0.          0.          8.57154746]]
737
738To obtain the rotated primitive cell basis vectors, we can use
739``standardize_cell`` as shown below::
740
741   import spglib
742
743   lattice = [[5.0759761474456697, 5.0759761474456697, 0],
744              [-2.8280307701821314, 2.8280307701821314, 0],
745              [0, 0, 8.57154746]]
746   points = [[0.0, 0.84688439, 0.1203133],
747             [0.0, 0.65311561, 0.6203133],
748             [0.0, 0.34688439, 0.3796867],
749             [0.0, 0.15311561, 0.8796867],
750             [0.5, 0.34688439, 0.1203133],
751             [0.5, 0.15311561, 0.6203133],
752             [0.5, 0.84688439, 0.3796867],
753             [0.5, 0.65311561, 0.8796867]]
754   numbers = [8,] * len(points)
755   cell = (lattice, points, numbers)
756   primitive_cell = spglib.standardize_cell(cell, to_primitive=1, no_idealize=1)
757   print(primitive_cell[0])
758
759then we get::
760
761   [[3.95200346 1.12397269 0.        ]
762    [1.12397269 3.95200346 0.        ]
763    [0.         0.         8.57154746]]
764
765which is equivalent to that we get manually. However using
766``standardize_cell``, distortion is not removed for the distorted
767crystal structure.
768
769Computing rigid rotation introduced by idealization
770^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
771
772In spglib, rigid rotation is purposely introduced in the idealization
773step though this is unlikely as a crystallographic operation.
774
775The crystal structure in the following script is the same as shown
776above, which is the one :math:`45^\circ` rotated around c-axis::
777
778    import spglib
779
780    # Mind that the a, b, c axes are given in row vectors here,
781    # but the formulation above is given for the column vectors.
782    lattice = [[5.0759761474456697, 5.0759761474456697, 0],  # a
783               [-2.8280307701821314, 2.8280307701821314, 0],  # b
784               [0, 0, 8.57154746]]  # c
785    points = [[0.0, 0.84688439, 0.1203133],
786              [0.0, 0.65311561, 0.6203133],
787              [0.0, 0.34688439, 0.3796867],
788              [0.0, 0.15311561, 0.8796867],
789              [0.5, 0.34688439, 0.1203133],
790              [0.5, 0.15311561, 0.6203133],
791              [0.5, 0.84688439, 0.3796867],
792              [0.5, 0.65311561, 0.8796867]]
793    numbers = [35,] * len(points)
794    cell = (lattice, points, numbers)
795    dataset = spglib.get_symmetry_dataset(cell)
796    print("Space group type: %s (%d)"
797          % (dataset['international'], dataset['number']))
798    print("Transformation matrix:")
799    for x in dataset['transformation_matrix']:
800        print("  %2d %2d %2d" % tuple(x))
801    print("std_lattice_after_idealization:")
802    print(dataset['std_lattice'])
803
804we get
805
806::
807
808   Space group type: Cmce (64)
809   Transformation matrix:
810      1  0  0
811      0  1  0
812      0  0  1
813   std_lattice_after_idealization:
814   [[7.17851431 0.         0.        ]
815    [0.         3.99943947 0.        ]
816    [0.         0.         8.57154746]]
817
818From Eq. :eq:`transformation_matrix`, the standardized basis vectors
819**before** idealization :math:`( \mathbf{a}_\mathrm{s} \; \mathbf{b}_\mathrm{s}
820\; \mathbf{c}_\mathrm{s} )` is obtained by (after ``import numpy as np``)
821
822::
823
824   std_lattice_before_idealization = np.dot(
825       np.transpose(lattice),
826       np.linalg.inv(dataset['transformation_matrix'])).T
827   print(std_lattice_before_idealization)
828
829which is (in row vectors)
830
831::
832
833   [[ 5.07597615  5.07597615  0.        ]
834    [-2.82803077  2.82803077  0.        ]
835    [ 0.          0.          8.57154746]]
836
837This is different from the standardized basis vectors **after**
838idealization :math:`( \bar{\mathbf{a}}_\mathrm{s} \;
839\bar{\mathbf{b}}_\mathrm{s} \; \bar{\mathbf{c}}_\mathrm{s} )`.  Unless
840this crystal strucutre is distorted from the crystal structure that
841has the ideal symmetry, this means that the crystal was rotated
842rigidly in the idealization step by
843
844.. math::
845
846   ( \bar{\mathbf{a}}_\mathrm{s} \;
847   \bar{\mathbf{b}}_\mathrm{s} \; \bar{\mathbf{c}}_\mathrm{s} )
848   = ( \boldsymbol{R} \mathbf{a}_\mathrm{s} \;
849   \boldsymbol{R} \mathbf{b}_\mathrm{s} \; \boldsymbol{R}
850   \mathbf{c}_\mathrm{s} ).
851
852where :math:`\boldsymbol{R}` is the rotation
853matrix. This is computed by
854
855::
856
857   R = np.dot(dataset['std_lattice'].T,
858              np.linalg.inv(std_lattice_before_idealization.T))
859   print(R)
860
861and we get
862
863::
864
865   [[ 0.70710678  0.70710678  0.        ]
866    [-0.70710678  0.70710678  0.        ]
867    [ 0.          0.          1.        ]]
868
869This equals to
870
871.. math::
872
873   \begin{pmatrix}
874   \cos\theta & -\sin\theta & 0 \\
875   \sin\theta & \cos\theta & 0 \\
876   0 & 0 & 1
877   \end{pmatrix},
878
879with :math:`\theta = -\pi/4` and :math:`\det(\boldsymbol{R})=1` when
880no distortion. ``dataset['std_rotation_matrix'])`` gives
881approximately the same result::
882
883   [[ 0.70710678  0.70710678  0.        ]
884    [-0.70710678  0.70710678  0.        ]
885    [ 0.          0.          1.        ]]
886
887In summary,
888
889.. math::
890
891
892   ( \bar{\mathbf{a}}_\mathrm{s} \;
893   \bar{\mathbf{b}}_\mathrm{s} \; \bar{\mathbf{c}}_\mathrm{s} )  \boldsymbol{P}
894   = ( \boldsymbol{R} \mathbf{a} \; \boldsymbol{R} \mathbf{b} \;
895   \boldsymbol{R} \mathbf{c} ).
896
897The atomic point coordinates in :math:`( \bar{\mathbf{a}}_\mathrm{s}
898\; \bar{\mathbf{b}}_\mathrm{s} \; \bar{\mathbf{c}}_\mathrm{s} )`
899are simply obtained by Eq. :eq:`change_of_basis_1` since the
900rotation doesn't affect them.
901