File: py_bindings.h

package info (click to toggle)
spglib 2.7.0-1
  • links: PTS, VCS
  • area: main
  • in suites: sid
  • size: 14,180 kB
  • sloc: ansic: 125,066; python: 7,717; cpp: 2,197; f90: 2,143; ruby: 792; makefile: 22; sh: 18
file content (208 lines) | stat: -rw-r--r-- 8,621 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
// Copyright (C) 2025 Spglib team
// SPDX-License-Identifier: BSD-3-Clause

#pragma once

#include <optional>

#include <pybind11/numpy.h>
#include <pybind11/pybind11.h>
#include <pybind11/stl.h>

namespace py = pybind11;
namespace spglib {
using array_double =
    py::array_t<double, py::array::c_style | py::array::forcecast>;
using array_int = py::array_t<int, py::array::c_style | py::array::forcecast>;
using array_uintp =
    py::array_t<uintptr_t, py::array::c_style | py::array::forcecast>;
using array_size_t =
    py::array_t<size_t, py::array::c_style | py::array::forcecast>;

// Specialized classes wrapping the numpy-like python arrays incorporating
// sanity checks like ndim and shape consistency checks.

class Lattice {
    array_double array;

   public:
    explicit Lattice(array_double &&_array);
    double (*data())[3];
    [[nodiscard]] double const (*data() const)[3];
};

class Rotations {
    array_int array;

   public:
    int const n_operations;
    explicit Rotations(array_int &&_array);
    int (*data())[3][3];
    [[nodiscard]] int const (*data() const)[3][3];
};

class Translations {
    array_double array;

   public:
    int const n_operations;
    explicit Translations(array_double &&_array);
    double (*data())[3];
    [[nodiscard]] double const (*data() const)[3];
};

class Symmetries {
   public:
    Rotations const rotations;
    Translations const translations;
    int const n_operations;
    Symmetries(Rotations &&_rotations, Translations &&_translations);
};

class Positions {
    array_double array;

   public:
    int const n_atoms;
    explicit Positions(array_double &&_array);
    double (*data())[3];
    [[nodiscard]] double const (*data() const)[3];
};

class AtomTypes {
    array_int array;

   public:
    int const n_atoms;
    explicit AtomTypes(array_int &&_array);
    int *data();
    [[nodiscard]]
    int const *data() const;
};

class Magmoms {
    array_double array;

   public:
    int const n_atoms;
    explicit Magmoms(array_double &&_array);
    double *data();
    [[nodiscard]]
    double const *data() const;
};

class Atoms {
   public:
    Positions const positions;
    AtomTypes const types;
    int const n_atoms;
    Atoms(Positions &&_positions, AtomTypes &&_types);
};

class SpglibError : public std::exception {
    std::string msg;

   public:
    SpglibError(std::string_view _msg);
    char const *what() const noexcept override;
};

py::tuple version_tuple();
py::str version_string();
py::str version_full();
py::str commit();
py::dict dataset(Lattice const &lattice, Positions const &positions,
                 AtomTypes const &atom_types, py::int_ hall_number,
                 py::float_ symprec, py::float_ angle_tolerance);
py::dict layer_dataset(Lattice const &lattice, Positions const &positions,
                       AtomTypes const &atom_types, py::int_ aperiodic_dir,
                       py::float_ symprec);
py::dict magnetic_dataset(Lattice const &lattice, Positions const &positions,
                          AtomTypes const &atom_types, array_double magmoms,
                          py::int_ tensor_rank, py::bool_ is_axial,
                          py::float_ symprec, py::float_ angle_tolerance,
                          py::float_ mag_symprec);
py::dict spacegroup_type(py::int_ hall_number);
py::dict spacegroup_type_from_symmetry(Rotations const &rotations,
                                       Translations const &translations,
                                       Lattice const &lattice,
                                       py::float_ symprec);
py::dict magnetic_spacegroup_type(py::int_ uni_number);
py::dict magnetic_spacegroup_type_from_symmetry(
    Rotations const &rotations, Translations const &translations,
    array_int time_reversals, Lattice const &lattice, py::float_ symprec);
py::int_ symmetry_from_database(Rotations &rotations,
                                Translations &translations,
                                py::int_ hall_number);
py::int_ magnetic_symmetry_from_database(Rotations &rotations,
                                         Translations &translations,
                                         array_int time_reversals,
                                         py::int_ uni_number,
                                         py::int_ hall_number);
py::tuple pointgroup(array_int rotations);
py::int_ standardize_cell(Lattice &lattice, Positions &positions,
                          array_int atom_types, py::int_ num_atom,
                          py::int_ to_primative, py::int_ no_idealize,
                          py::float_ symprec, py::float_ angle_tolerance);
py::int_ refine_cell(Lattice &lattice, Positions &positions,
                     AtomTypes &atom_types, py::int_ num_atom,
                     py::float_ symprec, py::float_ angle_tolerance);
py::int_ symmetry(Rotations &rotations, Translations &translations,
                  Lattice const &lattice, Positions const &positions,
                  AtomTypes const &atom_types, py::float_ symprec,
                  py::float_ angle_tolerance);
py::int_ symmetry_with_collinear_spin(
    Rotations &rotations, Translations &translations, array_int equiv_atoms,
    Lattice const &lattice, Positions const &positions,
    AtomTypes const &atom_types, array_double magmoms, py::float_ symprec,
    py::float_ angle_tolerance);
py::int_ symmetry_with_site_tensors(
    Rotations &rotations, Translations &translations, array_int equiv_atoms,
    Lattice &primitive_lattice, array_int spin_flips, Lattice const &lattice,
    Positions const &positions, AtomTypes const &atom_types,
    array_double tensors, py::int_ with_time_reversal, py::int_ is_axial,
    py::float_ symprec, py::float_ angle_tolerance, py::float_ mag_symprec);
py::int_ primitive(Lattice &lattice, Positions &positions,
                   AtomTypes &atom_types, py::float_ symprec,
                   py::float_ angle_tolerance);
py::int_ grid_point_from_address(array_int grid_address, array_int mesh);
py::int_ ir_reciprocal_mesh(array_int grid_address,
                            array_int grid_mapping_table, array_int mesh,
                            array_int is_shift, py::int_ is_time_reversal,
                            Lattice const &lattice, Positions const &positions,
                            AtomTypes const &atom_types, py::float_ symprec);
py::int_ ir_reciprocal_mesh(array_int grid_address,
                            array_size_t grid_mapping_table, array_int mesh,
                            array_int is_shift, py::int_ is_time_reversal,
                            Lattice const &lattice, Positions const &positions,
                            AtomTypes const &atom_types, py::float_ symprec);
py::int_ stabilized_reciprocal_mesh(array_int grid_address,
                                    array_int grid_mapping_table,
                                    array_int mesh, array_int is_shift,
                                    py::int_ is_time_reversal,
                                    Rotations const &rotations,
                                    array_double qpoints);
py::int_ stabilized_reciprocal_mesh(array_int grid_address,
                                    array_size_t grid_mapping_table,
                                    array_int mesh, array_int is_shift,
                                    py::int_ is_time_reversal,
                                    Rotations const &rotations,
                                    array_double qpoints);
void grid_points_by_rotations(array_size_t rot_grid_points,
                              array_int address_orig,
                              Rotations const &rot_reciprocal, array_int mesh,
                              array_int is_shift);
void BZ_grid_points_by_rotations(array_size_t rot_grid_points,
                                 array_int address_orig,
                                 Rotations const &rot_reciprocal,
                                 array_int mesh, array_int is_shift,
                                 array_size_t bz_map);
py::int_ BZ_grid_address(array_int bz_grid_address, array_size_t bz_map,
                         array_int grid_address, array_int mesh,
                         Lattice const &reciprocal_lattice, array_int is_shift);
py::int_ delaunay_reduce(Lattice &lattice, py::float_ symprec);
py::int_ niggli_reduce(Lattice &lattice, py::float_ eps);
py::int_ hall_number_from_symmetry(Rotations const &rotations,
                                   Translations const &translations,
                                   py::float_ symprec);
}  // namespace spglib