Coverage for gpaw/spherical_harmonics.py: 99%
116 statements
« prev ^ index » next coverage.py v7.7.1, created at 2025-07-20 00:19 +0000
« prev ^ index » next coverage.py v7.7.1, created at 2025-07-20 00:19 +0000
1r"""
2Real-valued spherical harmonics
5=== === === =======
6 L l m
7=== === === =======
8 0 0 0 1
9 1 1 -1 y
10 2 1 0 z
11 3 1 1 x
12 4 2 -2 xy
13 5 2 -1 yz
14 6 2 0 3z2-r2
15 7 2 1 zx
16 8 2 2 x2-y2
17=== === === =======
19For a more complete list, see c/bmgs/sharmonic.py
20"""
22import numpy as np
24from math import pi
25from collections import defaultdict
26from gpaw import debug, GPAW_NO_C_EXTENSION
28__all__ = ['Y', 'YL', 'nablarlYL', 'Yl']
30if GPAW_NO_C_EXTENSION:
31 Yl = None
32else:
33 from gpaw.cgpaw import spherical_harmonics
34 if debug:
35 def Yl(l: int, R_v: np.ndarray, rlY_m: np.ndarray) -> None:
36 assert R_v.dtype == float and R_v.shape == (3,)
37 assert rlY_m.dtype == float and rlY_m.shape == (2 * l + 1,)
38 spherical_harmonics(l, R_v, rlY_m)
39 else:
40 Yl = spherical_harmonics
42names = [['1'],
43 ['y', 'z', 'x'],
44 ['xy', 'yz', '3z2-r2', 'zx', 'x2-y2'],
45 ['3x2y-y3', 'xyz', '4yz2-y3-x2y', '2z3-3x2z-3y2z', '4xz2-x3-xy2',
46 'x2z-y2z', 'x3-3xy2']]
49def Y(L, x, y, z):
50 result = 0.0
51 for c, n in YL[L]:
52 result += c * x**n[0] * y**n[1] * z**n[2]
53 return result
56def Yarr(L_M, R_Av):
57 """
58 Calculate spherical harmonics L_M at positions R_Av, where
59 A is some array like index.
60 """
61 Y_MA = np.zeros((len(L_M), *R_Av.shape[:-1]))
62 for M, L in enumerate(L_M):
63 for c, n in YL[L]: # could be vectorized further
64 Y_MA[M] += c * np.prod(np.power(R_Av, n), axis=-1)
65 return Y_MA
68def nablarlYL(L, R):
69 """Calculate the gradient of a real solid spherical harmonic."""
70 x, y, z = R
71 dYdx = dYdy = dYdz = 0.0
72 terms = YL[L]
73 # The 'abs' avoids error in case powx == 0
74 for N, (powx, powy, powz) in terms:
75 dYdx += N * powx * x**abs(powx - 1) * y**powy * z**powz
76 dYdy += N * powy * x**powx * y**abs(powy - 1) * z**powz
77 dYdz += N * powz * x**powx * y**powy * z**abs(powz - 1)
78 return dYdx, dYdy, dYdz
81g = [1.0]
82for l in range(11):
83 g.append(g[-1] * (l + 0.5))
86def gam(n0, n1, n2):
87 h0 = n0 // 2
88 h1 = n1 // 2
89 h2 = n2 // 2
90 if 2 * h0 != n0 or 2 * h1 != n1 or 2 * h2 != n2:
91 return 0.0
92 return 2.0 * pi * g[h0] * g[h1] * g[h2] / g[1 + h0 + h1 + h2]
95def Y0(l, m):
96 """Sympy version of spherical harmonics."""
97 from fractions import Fraction
98 from sympy import assoc_legendre, sqrt, simplify, factorial as fac, I, pi
99 from sympy.abc import x, y, z
100 c = sqrt((2 * l + 1) * fac(l - m) / fac(l + m) / 4 / pi)
101 if m > 0:
102 return simplify(c * (x + I * y)**m / (1 - z**2)**Fraction(m, 2) *
103 assoc_legendre(l, m, z))
104 return simplify(c * (x - I * y)**(-m) / (1 - z**2)**Fraction(-m, 2) *
105 assoc_legendre(l, m, z))
108def S(l, m):
109 """Sympy version of real valued spherical harmonics."""
110 from sympy import I, Number, sqrt
111 if m > 0:
112 return (Y0(l, m) + (-1)**m * Y0(l, -m)) / sqrt(2) * (-1)**m
113 if m < 0:
114 return -(Y0(l, m) - Number(-1)**m * Y0(l, -m)) / (sqrt(2) * I)
115 return Y0(l, m)
118def poly_coeffs(l, m):
119 """Sympy coefficients for polynomiunm in x, y and z."""
120 from sympy import Poly
121 from sympy.abc import x, y, z
122 Y = S(l, m)
123 coeffs = {}
124 for nx, coef in enumerate(reversed(Poly(Y, x).all_coeffs())):
125 for ny, coef in enumerate(reversed(Poly(coef, y).all_coeffs())):
126 for nz, coef in enumerate(reversed(Poly(coef, z).all_coeffs())):
127 if coef:
128 coeffs[(nx, ny, nz)] = coef
129 return coeffs
132def fix_exponents(coeffs, l):
133 """Make sure exponents add up to l."""
134 from sympy import Number
135 new = defaultdict(lambda: Number(0))
136 for (nx, ny, nz), coef in coeffs.items():
137 if nx + ny + nz == l:
138 new[(nx, ny, nz)] += coef
139 else:
140 new[(nx + 2, ny, nz)] += coef
141 new[(nx, ny + 2, nz)] += coef
142 new[(nx, ny, nz + 2)] += coef
144 new = {nxyz: coef for nxyz, coef in new.items() if coef}
146 if not all(sum(nxyz) == l for nxyz in new):
147 new = fix_exponents(new, l)
149 return new
152def print_YL_table_code():
153 """Generate YL table using sympy.
155 This will generate slightly more accurate numbers, but we will not update
156 right now because then we would also have to update
157 c/bmgs/spherical_harminics.c.
158 """
159 print('# Computer generated table - do not touch!')
160 print('YL = [')
161 print(' # s, l=0:')
162 print(f' [({(4 * pi)**-0.5}, (0, 0, 0))],')
163 for l in range(1, 9):
164 s = 'spdfghijk'[l]
165 print(f' # {s}, l={l}:')
166 for m in range(-l, l + 1):
167 e = poly_coeffs(l, m)
168 e = fix_exponents(e, l)
169 if l**2 + m + l < len(YL):
170 assert len(e) == len(YL[l**2 + m + l])
171 for c0, n in YL[l**2 + m + l]:
172 c = e[n].evalf()
173 assert abs(c0 - c) < 1e-10
174 terms = []
175 for n, en in e.items():
176 c = float(en)
177 terms.append(f'({c!r}, {n})')
178 print(' [' + ',\n '.join(terms) + '],')
179 print(']')
182def write_c_code(l: int) -> None:
183 print(f' else if (l == {l})')
184 print(' {')
185 for m in range(2 * l + 1):
186 terms = []
187 for c, n in YL[l**2 + m]:
188 terms.append(f'{c!r} * ' + '*'.join('x' * n[0] +
189 'y' * n[1] +
190 'z' * n[2]))
191 print(f' Y_m[{m}] = ' + ' + '.join(terms) + ';')
192 print(' }')
195# Computer generated table - do not touch!
196# The numbers match those in c/bmgs/spherical_harmonics.c and were
197# originally generated with c/bmgs/sharmonic.py (old Python 2 code).
198YL = [
199 # s, l=0:
200 [(0.28209479177387814, (0, 0, 0))],
201 # p, l=1:
202 [(0.4886025119029199, (0, 1, 0))],
203 [(0.4886025119029199, (0, 0, 1))],
204 [(0.4886025119029199, (1, 0, 0))],
205 # d, l=2:
206 [(1.0925484305920792, (1, 1, 0))],
207 [(1.0925484305920792, (0, 1, 1))],
208 [(0.6307831305050401, (0, 0, 2)),
209 (-0.31539156525252005, (0, 2, 0)),
210 (-0.31539156525252005, (2, 0, 0))],
211 [(1.0925484305920792, (1, 0, 1))],
212 [(0.5462742152960396, (2, 0, 0)),
213 (-0.5462742152960396, (0, 2, 0))],
214 # f, l=3:
215 [(-0.5900435899266435, (0, 3, 0)),
216 (1.7701307697799304, (2, 1, 0))],
217 [(2.890611442640554, (1, 1, 1))],
218 [(-0.4570457994644658, (0, 3, 0)),
219 (1.828183197857863, (0, 1, 2)),
220 (-0.4570457994644658, (2, 1, 0))],
221 [(0.7463526651802308, (0, 0, 3)),
222 (-1.1195289977703462, (2, 0, 1)),
223 (-1.1195289977703462, (0, 2, 1))],
224 [(1.828183197857863, (1, 0, 2)),
225 (-0.4570457994644658, (3, 0, 0)),
226 (-0.4570457994644658, (1, 2, 0))],
227 [(1.445305721320277, (2, 0, 1)),
228 (-1.445305721320277, (0, 2, 1))],
229 [(0.5900435899266435, (3, 0, 0)),
230 (-1.7701307697799304, (1, 2, 0))],
231 # g, l=4:
232 [(2.5033429417967046, (3, 1, 0)),
233 (-2.5033429417967046, (1, 3, 0))],
234 [(-1.7701307697799307, (0, 3, 1)),
235 (5.310392309339792, (2, 1, 1))],
236 [(-0.9461746957575601, (3, 1, 0)),
237 (-0.9461746957575601, (1, 3, 0)),
238 (5.6770481745453605, (1, 1, 2))],
239 [(-2.0071396306718676, (2, 1, 1)),
240 (2.676186174229157, (0, 1, 3)),
241 (-2.0071396306718676, (0, 3, 1))],
242 [(0.6347132814912259, (2, 2, 0)),
243 (-2.5388531259649034, (2, 0, 2)),
244 (0.31735664074561293, (0, 4, 0)),
245 (-2.5388531259649034, (0, 2, 2)),
246 (0.31735664074561293, (4, 0, 0)),
247 (0.8462843753216345, (0, 0, 4))],
248 [(2.676186174229157, (1, 0, 3)),
249 (-2.0071396306718676, (3, 0, 1)),
250 (-2.0071396306718676, (1, 2, 1))],
251 [(2.8385240872726802, (2, 0, 2)),
252 (0.47308734787878004, (0, 4, 0)),
253 (-0.47308734787878004, (4, 0, 0)),
254 (-2.8385240872726802, (0, 2, 2))],
255 [(1.7701307697799307, (3, 0, 1)),
256 (-5.310392309339792, (1, 2, 1))],
257 [(-3.755014412695057, (2, 2, 0)),
258 (0.6258357354491761, (0, 4, 0)),
259 (0.6258357354491761, (4, 0, 0))],
260 # h, l=5:
261 [(-6.5638205684017015, (2, 3, 0)),
262 (3.2819102842008507, (4, 1, 0)),
263 (0.6563820568401701, (0, 5, 0))],
264 [(8.302649259524165, (3, 1, 1)),
265 (-8.302649259524165, (1, 3, 1))],
266 [(-3.913906395482003, (0, 3, 2)),
267 (0.4892382994352504, (0, 5, 0)),
268 (-1.467714898305751, (4, 1, 0)),
269 (-0.9784765988705008, (2, 3, 0)),
270 (11.741719186446009, (2, 1, 2))],
271 [(-4.793536784973324, (3, 1, 1)),
272 (-4.793536784973324, (1, 3, 1)),
273 (9.587073569946648, (1, 1, 3))],
274 [(-5.435359814348363, (0, 3, 2)),
275 (0.9058933023913939, (2, 3, 0)),
276 (-5.435359814348363, (2, 1, 2)),
277 (3.6235732095655755, (0, 1, 4)),
278 (0.45294665119569694, (4, 1, 0)),
279 (0.45294665119569694, (0, 5, 0))],
280 [(3.508509673602708, (2, 2, 1)),
281 (-4.678012898136944, (0, 2, 3)),
282 (1.754254836801354, (0, 4, 1)),
283 (-4.678012898136944, (2, 0, 3)),
284 (1.754254836801354, (4, 0, 1)),
285 (0.9356025796273888, (0, 0, 5))],
286 [(-5.435359814348363, (3, 0, 2)),
287 (3.6235732095655755, (1, 0, 4)),
288 (0.45294665119569694, (5, 0, 0)),
289 (0.9058933023913939, (3, 2, 0)),
290 (-5.435359814348363, (1, 2, 2)),
291 (0.45294665119569694, (1, 4, 0))],
292 [(-2.396768392486662, (4, 0, 1)),
293 (2.396768392486662, (0, 4, 1)),
294 (4.793536784973324, (2, 0, 3)),
295 (-4.793536784973324, (0, 2, 3))],
296 [(3.913906395482003, (3, 0, 2)),
297 (-0.4892382994352504, (5, 0, 0)),
298 (0.9784765988705008, (3, 2, 0)),
299 (-11.741719186446009, (1, 2, 2)),
300 (1.467714898305751, (1, 4, 0))],
301 [(2.075662314881041, (4, 0, 1)),
302 (-12.453973889286246, (2, 2, 1)),
303 (2.075662314881041, (0, 4, 1))],
304 [(-6.5638205684017015, (3, 2, 0)),
305 (0.6563820568401701, (5, 0, 0)),
306 (3.2819102842008507, (1, 4, 0))],
307 # i, l=6:
308 [(4.099104631151485, (5, 1, 0)),
309 (-13.663682103838287, (3, 3, 0)),
310 (4.099104631151485, (1, 5, 0))],
311 [(11.83309581115876, (4, 1, 1)),
312 (-23.66619162231752, (2, 3, 1)),
313 (2.366619162231752, (0, 5, 1))],
314 [(20.182596029148968, (3, 1, 2)),
315 (-2.0182596029148967, (5, 1, 0)),
316 (2.0182596029148967, (1, 5, 0)),
317 (-20.182596029148968, (1, 3, 2))],
318 [(-7.369642076119388, (0, 3, 3)),
319 (-5.527231557089541, (2, 3, 1)),
320 (2.7636157785447706, (0, 5, 1)),
321 (22.108926228358165, (2, 1, 3)),
322 (-8.29084733563431, (4, 1, 1))],
323 [(-14.739284152238776, (3, 1, 2)),
324 (14.739284152238776, (1, 1, 4)),
325 (1.842410519029847, (3, 3, 0)),
326 (0.9212052595149235, (5, 1, 0)),
327 (-14.739284152238776, (1, 3, 2)),
328 (0.9212052595149235, (1, 5, 0))],
329 [(2.9131068125936572, (0, 5, 1)),
330 (-11.652427250374629, (0, 3, 3)),
331 (5.8262136251873144, (2, 3, 1)),
332 (-11.652427250374629, (2, 1, 3)),
333 (2.9131068125936572, (4, 1, 1)),
334 (4.660970900149851, (0, 1, 5))],
335 [(5.721228204086558, (4, 0, 2)),
336 (-7.628304272115411, (0, 2, 4)),
337 (-0.9535380340144264, (2, 4, 0)),
338 (1.0171072362820548, (0, 0, 6)),
339 (-0.9535380340144264, (4, 2, 0)),
340 (5.721228204086558, (0, 4, 2)),
341 (-0.3178460113381421, (0, 6, 0)),
342 (-7.628304272115411, (2, 0, 4)),
343 (-0.3178460113381421, (6, 0, 0)),
344 (11.442456408173117, (2, 2, 2))],
345 [(-11.652427250374629, (3, 0, 3)),
346 (4.660970900149851, (1, 0, 5)),
347 (2.9131068125936572, (5, 0, 1)),
348 (5.8262136251873144, (3, 2, 1)),
349 (-11.652427250374629, (1, 2, 3)),
350 (2.9131068125936572, (1, 4, 1))],
351 [(7.369642076119388, (2, 0, 4)),
352 (-7.369642076119388, (0, 2, 4)),
353 (-0.46060262975746175, (2, 4, 0)),
354 (-7.369642076119388, (4, 0, 2)),
355 (0.46060262975746175, (4, 2, 0)),
356 (-0.46060262975746175, (0, 6, 0)),
357 (0.46060262975746175, (6, 0, 0)),
358 (7.369642076119388, (0, 4, 2))],
359 [(7.369642076119388, (3, 0, 3)),
360 (-2.7636157785447706, (5, 0, 1)),
361 (5.527231557089541, (3, 2, 1)),
362 (-22.108926228358165, (1, 2, 3)),
363 (8.29084733563431, (1, 4, 1))],
364 [(2.522824503643621, (4, 2, 0)),
365 (5.045649007287242, (0, 4, 2)),
366 (-30.273894043723452, (2, 2, 2)),
367 (-0.5045649007287242, (0, 6, 0)),
368 (-0.5045649007287242, (6, 0, 0)),
369 (5.045649007287242, (4, 0, 2)),
370 (2.522824503643621, (2, 4, 0))],
371 [(2.366619162231752, (5, 0, 1)),
372 (-23.66619162231752, (3, 2, 1)),
373 (11.83309581115876, (1, 4, 1))],
374 [(-10.247761577878714, (4, 2, 0)),
375 (-0.6831841051919143, (0, 6, 0)),
376 (0.6831841051919143, (6, 0, 0)),
377 (10.247761577878714, (2, 4, 0))],
378 # j, l=7:
379 [(14.850417383016522, (2, 5, 0)),
380 (4.950139127672174, (6, 1, 0)),
381 (-24.75069563836087, (4, 3, 0)),
382 (-0.7071627325245963, (0, 7, 0))],
383 [(-52.91921323603801, (3, 3, 1)),
384 (15.875763970811402, (1, 5, 1)),
385 (15.875763970811402, (5, 1, 1))],
386 [(-2.5945778936013015, (6, 1, 0)),
387 (2.5945778936013015, (4, 3, 0)),
388 (-62.26986944643124, (2, 3, 2)),
389 (4.670240208482342, (2, 5, 0)),
390 (6.226986944643123, (0, 5, 2)),
391 (31.13493472321562, (4, 1, 2)),
392 (-0.5189155787202603, (0, 7, 0))],
393 [(41.513246297620825, (3, 1, 3)),
394 (12.453973889286246, (1, 5, 1)),
395 (-41.513246297620825, (1, 3, 3)),
396 (-12.453973889286246, (5, 1, 1))],
397 [(-18.775072063475285, (2, 3, 2)),
398 (-0.4693768015868821, (0, 7, 0)),
399 (0.4693768015868821, (2, 5, 0)),
400 (2.3468840079344107, (4, 3, 0)),
401 (-12.516714708983523, (0, 3, 4)),
402 (37.55014412695057, (2, 1, 4)),
403 (1.4081304047606462, (6, 1, 0)),
404 (9.387536031737643, (0, 5, 2)),
405 (-28.162608095212928, (4, 1, 2))],
406 [(13.27598077334948, (3, 3, 1)),
407 (6.63799038667474, (1, 5, 1)),
408 (-35.402615395598616, (3, 1, 3)),
409 (21.24156923735917, (1, 1, 5)),
410 (-35.402615395598616, (1, 3, 3)),
411 (6.63799038667474, (5, 1, 1))],
412 [(-0.4516580379125865, (0, 7, 0)),
413 (10.839792909902076, (0, 5, 2)),
414 (-1.3549741137377596, (2, 5, 0)),
415 (-1.3549741137377596, (4, 3, 0)),
416 (-21.679585819804153, (0, 3, 4)),
417 (-21.679585819804153, (2, 1, 4)),
418 (5.781222885281108, (0, 1, 6)),
419 (-0.4516580379125865, (6, 1, 0)),
420 (21.679585819804153, (2, 3, 2)),
421 (10.839792909902076, (4, 1, 2))],
422 [(-11.471758521216831, (2, 0, 5)),
423 (1.0925484305920792, (0, 0, 7)),
424 (-11.471758521216831, (0, 2, 5)),
425 (28.67939630304208, (2, 2, 3)),
426 (-2.3899496919201733, (6, 0, 1)),
427 (-7.16984907576052, (4, 2, 1)),
428 (14.33969815152104, (4, 0, 3)),
429 (-2.3899496919201733, (0, 6, 1)),
430 (-7.16984907576052, (2, 4, 1)),
431 (14.33969815152104, (0, 4, 3))],
432 [(10.839792909902076, (1, 4, 2)),
433 (-0.4516580379125865, (7, 0, 0)),
434 (21.679585819804153, (3, 2, 2)),
435 (-1.3549741137377596, (5, 2, 0)),
436 (-0.4516580379125865, (1, 6, 0)),
437 (-21.679585819804153, (3, 0, 4)),
438 (-1.3549741137377596, (3, 4, 0)),
439 (5.781222885281108, (1, 0, 6)),
440 (-21.679585819804153, (1, 2, 4)),
441 (10.839792909902076, (5, 0, 2))],
442 [(10.620784618679584, (2, 0, 5)),
443 (-10.620784618679584, (0, 2, 5)),
444 (3.31899519333737, (6, 0, 1)),
445 (3.31899519333737, (4, 2, 1)),
446 (-17.701307697799308, (4, 0, 3)),
447 (-3.31899519333737, (0, 6, 1)),
448 (-3.31899519333737, (2, 4, 1)),
449 (17.701307697799308, (0, 4, 3))],
450 [(-1.4081304047606462, (1, 6, 0)),
451 (0.4693768015868821, (7, 0, 0)),
452 (18.775072063475285, (3, 2, 2)),
453 (-0.4693768015868821, (5, 2, 0)),
454 (12.516714708983523, (3, 0, 4)),
455 (-2.3468840079344107, (3, 4, 0)),
456 (28.162608095212928, (1, 4, 2)),
457 (-37.55014412695057, (1, 2, 4)),
458 (-9.387536031737643, (5, 0, 2))],
459 [(10.378311574405206, (4, 0, 3)),
460 (-3.1134934723215615, (0, 6, 1)),
461 (15.56746736160781, (2, 4, 1)),
462 (-62.26986944643124, (2, 2, 3)),
463 (10.378311574405206, (0, 4, 3)),
464 (-3.1134934723215615, (6, 0, 1)),
465 (15.56746736160781, (4, 2, 1))],
466 [(-2.5945778936013015, (1, 6, 0)),
467 (-62.26986944643124, (3, 2, 2)),
468 (-0.5189155787202603, (7, 0, 0)),
469 (31.13493472321562, (1, 4, 2)),
470 (2.5945778936013015, (3, 4, 0)),
471 (6.226986944643123, (5, 0, 2)),
472 (4.670240208482342, (5, 2, 0))],
473 [(2.6459606618019005, (6, 0, 1)),
474 (-2.6459606618019005, (0, 6, 1)),
475 (-39.68940992702851, (4, 2, 1)),
476 (39.68940992702851, (2, 4, 1))],
477 [(0.7071627325245963, (7, 0, 0)),
478 (-14.850417383016522, (5, 2, 0)),
479 (24.75069563836087, (3, 4, 0)),
480 (-4.950139127672174, (1, 6, 0))],
481 # k, l=8:
482 [(-5.831413281398639, (1, 7, 0)),
483 (40.81989296979047, (3, 5, 0)),
484 (-40.81989296979047, (5, 3, 0)),
485 (5.831413281398639, (7, 1, 0))],
486 [(-2.9157066406993195, (0, 7, 1)),
487 (61.22983945468571, (2, 5, 1)),
488 (-102.04973242447618, (4, 3, 1)),
489 (20.409946484895237, (6, 1, 1))],
490 [(7.452658724833596, (3, 5, 0)),
491 (-3.1939965963572554, (1, 7, 0)),
492 (44.715952349001576, (1, 5, 2)),
493 (7.452658724833596, (5, 3, 0)),
494 (-149.05317449667191, (3, 3, 2)),
495 (-3.1939965963572554, (7, 1, 0)),
496 (44.715952349001576, (5, 1, 2))],
497 [(31.04919559888297, (2, 5, 1)),
498 (-3.449910622098108, (0, 7, 1)),
499 (13.799642488392433, (0, 5, 3)),
500 (17.24955311049054, (4, 3, 1)),
501 (-137.99642488392433, (2, 3, 3)),
502 (-17.24955311049054, (6, 1, 1)),
503 (68.99821244196217, (4, 1, 3))],
504 [(-1.9136660990373229, (3, 5, 0)),
505 (-1.9136660990373229, (1, 7, 0)),
506 (45.927986376895745, (1, 5, 2)),
507 (-76.54664396149292, (1, 3, 4)),
508 (1.9136660990373229, (7, 1, 0)),
509 (1.9136660990373229, (5, 3, 0)),
510 (-45.927986376895745, (5, 1, 2)),
511 (76.54664396149292, (3, 1, 4))],
512 [(18.528992329433162, (4, 3, 1)),
513 (3.705798465886632, (2, 5, 1)),
514 (-49.41064621182176, (2, 3, 3)),
515 (-3.705798465886632, (0, 7, 1)),
516 (24.70532310591088, (0, 5, 3)),
517 (-19.764258484728707, (0, 3, 5)),
518 (11.117395397659896, (6, 1, 1)),
519 (-74.11596931773265, (4, 1, 3)),
520 (59.29277545418611, (2, 1, 5))],
521 [(-0.912304516869819, (7, 1, 0)),
522 (-2.7369135506094566, (5, 3, 0)),
523 (27.36913550609457, (5, 1, 2)),
524 (-2.7369135506094566, (3, 5, 0)),
525 (54.73827101218914, (3, 3, 2)),
526 (-72.98436134958551, (3, 1, 4)),
527 (-0.912304516869819, (1, 7, 0)),
528 (27.36913550609457, (1, 5, 2)),
529 (-72.98436134958551, (1, 3, 4)),
530 (29.193744539834206, (1, 1, 6))],
531 [(-3.8164436064572986, (6, 1, 1)),
532 (-11.449330819371895, (4, 3, 1)),
533 (30.53154885165839, (4, 1, 3)),
534 (-11.449330819371895, (2, 5, 1)),
535 (61.06309770331678, (2, 3, 3)),
536 (-36.63785862199006, (2, 1, 5)),
537 (-3.8164436064572986, (0, 7, 1)),
538 (30.53154885165839, (0, 5, 3)),
539 (-36.63785862199006, (0, 3, 5)),
540 (6.978639737521917, (0, 1, 7))],
541 [(0.31803696720477487, (8, 0, 0)),
542 (1.2721478688190995, (6, 2, 0)),
543 (-10.177182950552796, (6, 0, 2)),
544 (1.9082218032286493, (4, 4, 0)),
545 (-30.53154885165839, (4, 2, 2)),
546 (30.53154885165839, (4, 0, 4)),
547 (1.2721478688190995, (2, 6, 0)),
548 (-30.53154885165839, (2, 4, 2)),
549 (61.06309770331678, (2, 2, 4)),
550 (-16.283492720884475, (2, 0, 6)),
551 (0.31803696720477487, (0, 8, 0)),
552 (-10.177182950552796, (0, 6, 2)),
553 (30.53154885165839, (0, 4, 4)),
554 (-16.283492720884475, (0, 2, 6)),
555 (1.1631066229203195, (0, 0, 8))],
556 [(-3.8164436064572986, (7, 0, 1)),
557 (-11.449330819371895, (5, 2, 1)),
558 (30.53154885165839, (5, 0, 3)),
559 (-11.449330819371895, (3, 4, 1)),
560 (61.06309770331678, (3, 2, 3)),
561 (-36.63785862199006, (3, 0, 5)),
562 (-3.8164436064572986, (1, 6, 1)),
563 (30.53154885165839, (1, 4, 3)),
564 (-36.63785862199006, (1, 2, 5)),
565 (6.978639737521917, (1, 0, 7))],
566 [(0.912304516869819, (2, 6, 0)),
567 (-13.684567753047284, (2, 4, 2)),
568 (0.4561522584349095, (0, 8, 0)),
569 (-13.684567753047284, (0, 6, 2)),
570 (36.492180674792756, (0, 4, 4)),
571 (-14.596872269917103, (0, 2, 6)),
572 (-0.4561522584349095, (8, 0, 0)),
573 (-0.912304516869819, (6, 2, 0)),
574 (13.684567753047284, (6, 0, 2)),
575 (13.684567753047284, (4, 2, 2)),
576 (-36.492180674792756, (4, 0, 4)),
577 (14.596872269917103, (2, 0, 6))],
578 [(-3.705798465886632, (5, 2, 1)),
579 (-18.528992329433162, (3, 4, 1)),
580 (49.41064621182176, (3, 2, 3)),
581 (-11.117395397659896, (1, 6, 1)),
582 (74.11596931773265, (1, 4, 3)),
583 (-59.29277545418611, (1, 2, 5)),
584 (3.705798465886632, (7, 0, 1)),
585 (-24.70532310591088, (5, 0, 3)),
586 (19.764258484728707, (3, 0, 5))],
587 [(-4.784165247593307, (4, 4, 0)),
588 (-1.9136660990373229, (2, 6, 0)),
589 (57.40998297111968, (2, 4, 2)),
590 (0.4784165247593307, (0, 8, 0)),
591 (-11.481996594223936, (0, 6, 2)),
592 (19.13666099037323, (0, 4, 4)),
593 (-1.9136660990373229, (6, 2, 0)),
594 (57.40998297111968, (4, 2, 2)),
595 (-114.81996594223936, (2, 2, 4)),
596 (0.4784165247593307, (8, 0, 0)),
597 (-11.481996594223936, (6, 0, 2)),
598 (19.13666099037323, (4, 0, 4))],
599 [(17.24955311049054, (3, 4, 1)),
600 (-17.24955311049054, (1, 6, 1)),
601 (68.99821244196217, (1, 4, 3)),
602 (31.04919559888297, (5, 2, 1)),
603 (-137.99642488392433, (3, 2, 3)),
604 (-3.449910622098108, (7, 0, 1)),
605 (13.799642488392433, (5, 0, 3))],
606 [(-7.452658724833596, (2, 6, 0)),
607 (0.5323327660595426, (0, 8, 0)),
608 (-7.452658724833596, (0, 6, 2)),
609 (111.78988087250394, (2, 4, 2)),
610 (7.452658724833596, (6, 2, 0)),
611 (-111.78988087250394, (4, 2, 2)),
612 (-0.5323327660595426, (8, 0, 0)),
613 (7.452658724833596, (6, 0, 2))],
614 [(-20.409946484895237, (1, 6, 1)),
615 (102.04973242447618, (3, 4, 1)),
616 (-61.22983945468571, (5, 2, 1)),
617 (2.9157066406993195, (7, 0, 1))],
618 [(0.7289266601748299, (0, 8, 0)),
619 (-20.409946484895237, (2, 6, 0)),
620 (51.02486621223809, (4, 4, 0)),
621 (-20.409946484895237, (6, 2, 0)),
622 (0.7289266601748299, (8, 0, 0))],
623]