Coverage for /builds/kinetik161/ase/ase/phasediagram.py: 66.57%
362 statements
« prev ^ index » next coverage.py v7.2.7, created at 2023-12-10 11:04 +0000
« prev ^ index » next coverage.py v7.2.7, created at 2023-12-10 11:04 +0000
1import fractions
2import functools
3import re
4from collections import OrderedDict
5from typing import Dict, List, Tuple
7import numpy as np
8from scipy.spatial import ConvexHull
10import ase.units as units
11from ase.formula import Formula
13_solvated: List[Tuple[str, Dict[str, int], float, bool, float]] = []
16def parse_formula(formula):
17 aq = formula.endswith('(aq)')
18 if aq:
19 formula = formula[:-4]
20 charge = formula.count('+') - formula.count('-')
21 if charge:
22 formula = formula.rstrip('+-')
23 count = Formula(formula).count()
24 return count, charge, aq
27def float2str(x):
28 f = fractions.Fraction(x).limit_denominator(100)
29 n = f.numerator
30 d = f.denominator
31 if abs(n / d - f) > 1e-6:
32 return f'{f:.3f}'
33 if d == 0:
34 return '0'
35 if f.denominator == 1:
36 return str(n)
37 return f'{f.numerator}/{f.denominator}'
40def solvated(symbols):
41 """Extract solvation energies from database.
43 symbols: str
44 Extract only those molecules that contain the chemical elements
45 given by the symbols string (plus water and H+).
47 Data from:
49 Johnson JW, Oelkers EH, Helgeson HC (1992)
50 Comput Geosci 18(7):899.
51 doi:10.1016/0098-3004(92)90029-Q
53 and:
55 Pourbaix M (1966)
56 Atlas of electrochemical equilibria in aqueous solutions.
57 No. v. 1 in Atlas of Electrochemical Equilibria in Aqueous Solutions.
58 Pergamon Press, New York.
60 Returns list of (name, energy) tuples.
61 """
63 if isinstance(symbols, str):
64 symbols = Formula(symbols).count().keys()
65 if len(_solvated) == 0:
66 for line in _aqueous.splitlines():
67 energy, formula = line.split(',')
68 name = formula + '(aq)'
69 count, charge, aq = parse_formula(name)
70 energy = float(energy) * 0.001 * units.kcal / units.mol
71 _solvated.append((name, count, charge, aq, energy))
72 references = []
73 for name, count, charge, aq, energy in _solvated:
74 for symbol in count:
75 if symbol not in 'HO' and symbol not in symbols:
76 break
77 else:
78 references.append((name, energy))
79 return references
82def bisect(A, X, Y, f):
83 a = []
84 for i in [0, -1]:
85 for j in [0, -1]:
86 if A[i, j] == -1:
87 A[i, j] = f(X[i], Y[j])
88 a.append(A[i, j])
90 if np.ptp(a) == 0:
91 A[:] = a[0]
92 return
93 if a[0] == a[1]:
94 A[0] = a[0]
95 if a[1] == a[3]:
96 A[:, -1] = a[1]
97 if a[3] == a[2]:
98 A[-1] = a[3]
99 if a[2] == a[0]:
100 A[:, 0] = a[2]
101 if not (A == -1).any():
102 return
103 i = len(X) // 2
104 j = len(Y) // 2
105 bisect(A[:i + 1, :j + 1], X[:i + 1], Y[:j + 1], f)
106 bisect(A[:i + 1, j:], X[:i + 1], Y[j:], f)
107 bisect(A[i:, :j + 1], X[i:], Y[:j + 1], f)
108 bisect(A[i:, j:], X[i:], Y[j:], f)
111def print_results(results):
112 total_energy = 0.0
113 print('reference coefficient energy')
114 print('------------------------------------')
115 for name, coef, energy in results:
116 total_energy += coef * energy
117 if abs(coef) < 1e-7:
118 continue
119 print(f'{name:14}{float2str(coef):>10}{energy:12.3f}')
120 print('------------------------------------')
121 print(f'Total energy: {total_energy:22.3f}')
122 print('------------------------------------')
125class Pourbaix:
126 def __init__(self, references, formula=None, T=300.0, **kwargs):
127 """Pourbaix object.
129 references: list of (name, energy) tuples
130 Examples of names: ZnO2, H+(aq), H2O(aq), Zn++(aq), ...
131 formula: str
132 Stoichiometry. Example: ``'ZnO'``. Can also be given as
133 keyword arguments: ``Pourbaix(refs, Zn=1, O=1)``.
134 T: float
135 Temperature in Kelvin.
136 """
138 if formula:
139 assert not kwargs
140 kwargs = parse_formula(formula)[0]
142 if 'O' not in kwargs:
143 kwargs['O'] = 0
144 if 'H' not in kwargs:
145 kwargs['H'] = 0
147 self.kT = units.kB * T
148 self.references = []
149 for name, energy in references:
150 if name == 'O':
151 continue
152 count, charge, aq = parse_formula(name)
153 if all(symbol in kwargs for symbol in count):
154 self.references.append((count, charge, aq, energy, name))
156 self.references.append(({}, -1, False, 0.0, 'e-')) # an electron
158 self.count = kwargs
160 self.N = {'e-': 0}
161 for symbol in kwargs:
162 if symbol not in self.N:
163 self.N[symbol] = len(self.N)
165 def decompose(self, U, pH, verbose=True, concentration=1e-6):
166 """Decompose material.
168 U: float
169 Potential in V.
170 pH: float
171 pH value.
172 verbose: bool
173 Default is True.
174 concentration: float
175 Concentration of solvated references.
177 Returns optimal coefficients and energy:
179 >>> from ase.phasediagram import Pourbaix, solvated
180 >>> refs = solvated('CoO') + [
181 ... ('Co', 0.0),
182 ... ('CoO', -2.509),
183 ... ('Co3O4', -9.402)]
184 >>> pb = Pourbaix(refs, Co=3, O=4)
185 >>> coefs, energy = pb.decompose(U=1.5, pH=0,
186 ... concentration=1e-6,
187 ... verbose=True)
188 0 HCoO2-(aq) -3.974
189 1 CoO2--(aq) -3.098
190 2 H2O(aq) -2.458
191 3 CoOH+(aq) -2.787
192 4 CoO(aq) -2.265
193 5 CoOH++(aq) -1.355
194 6 Co++(aq) -0.921
195 7 H+(aq) 0.000
196 8 Co+++(aq) 1.030
197 9 Co 0.000
198 10 CoO -2.509
199 11 Co3O4 -9.402
200 12 e- -1.500
201 reference coefficient energy
202 ------------------------------------
203 H2O(aq) 4 -2.458
204 Co++(aq) 3 -0.921
205 H+(aq) -8 0.000
206 e- -2 -1.500
207 ------------------------------------
208 Total energy: -9.596
209 ------------------------------------
210 """
212 alpha = np.log(10) * self.kT
213 entropy = -np.log(concentration) * self.kT
215 # We want to minimize np.dot(energies, x) under the constraints:
216 #
217 # np.dot(x, eq2) == eq1
218 #
219 # with bounds[i,0] <= x[i] <= bounds[i, 1].
220 #
221 # First two equations are charge and number of hydrogens, and
222 # the rest are the remaining species.
224 eq1 = [0] + list(self.count.values())
225 eq2 = []
226 energies = []
227 bounds = []
228 names = []
229 for count, charge, aq, energy, name in self.references:
230 eq = np.zeros(len(self.N))
231 eq[0] = charge
232 for symbol, n in count.items():
233 eq[self.N[symbol]] = n
234 eq2.append(eq)
235 if name in ['H2O(aq)', 'H+(aq)', 'e-']:
236 bounds.append((-np.inf, np.inf))
237 if name == 'e-':
238 energy = -U
239 elif name == 'H+(aq)':
240 energy = -pH * alpha
241 else:
242 bounds.append((0, np.inf))
243 if aq:
244 energy -= entropy
245 if verbose:
246 print('{:<5}{:10}{:10.3f}'.format(len(energies),
247 name, energy))
248 energies.append(energy)
249 names.append(name)
251 from scipy.optimize import linprog
253 result = linprog(c=energies,
254 A_eq=np.transpose(eq2),
255 b_eq=eq1,
256 bounds=bounds)
258 if verbose:
259 print_results(zip(names, result.x, energies))
261 return result.x, result.fun
263 def diagram(self, U, pH, plot=True, show=False, ax=None):
264 """Calculate Pourbaix diagram.
266 U: list of float
267 Potentials in V.
268 pH: list of float
269 pH values.
270 plot: bool
271 Create plot.
272 show: bool
273 Open graphical window and show plot.
274 ax: matplotlib axes object
275 When creating plot, plot onto the given axes object.
276 If none given, plot onto the current one.
277 """
278 a = np.empty((len(U), len(pH)), int)
279 a[:] = -1
280 colors = {}
281 f = functools.partial(self.colorfunction, colors=colors)
282 bisect(a, U, pH, f)
283 compositions = [None] * len(colors)
284 names = [ref[-1] for ref in self.references]
285 for indices, color in colors.items():
286 compositions[color] = ' + '.join(names[i] for i in indices
287 if names[i] not in
288 ['H2O(aq)', 'H+(aq)', 'e-'])
289 text = []
290 for i, name in enumerate(compositions):
291 b = (a == i)
292 x = np.dot(b.sum(1), U) / b.sum()
293 y = np.dot(b.sum(0), pH) / b.sum()
294 name = re.sub(r'(\S)([+-]+)', r'\1$^{\2}$', name)
295 name = re.sub(r'(\d+)', r'$_{\1}$', name)
296 text.append((x, y, name))
298 if plot:
299 import matplotlib.cm as cm
300 import matplotlib.pyplot as plt
301 if ax is None:
302 ax = plt.gca()
304 # rasterized pcolormesh has a bug which leaves a tiny
305 # white border. Unrasterized pcolormesh produces
306 # unreasonably large files. Avoid this by using the more
307 # general imshow.
308 ax.imshow(a, cmap=cm.Accent,
309 extent=[min(pH), max(pH), min(U), max(U)],
310 origin='lower',
311 aspect='auto')
313 for x, y, name in text:
314 ax.text(y, x, name, horizontalalignment='center')
315 ax.set_xlabel('pH')
316 ax.set_ylabel('potential [V]')
317 ax.set_xlim(min(pH), max(pH))
318 ax.set_ylim(min(U), max(U))
319 if show:
320 plt.show()
322 return a, compositions, text
324 def colorfunction(self, U, pH, colors):
325 coefs, energy = self.decompose(U, pH, verbose=False)
326 indices = tuple(sorted(np.where(abs(coefs) > 1e-3)[0]))
327 color = colors.get(indices)
328 if color is None:
329 color = len(colors)
330 colors[indices] = color
331 return color
334class PhaseDiagram:
335 def __init__(self, references, filter='', verbose=True):
336 """Phase-diagram.
338 references: list of (name, energy) tuples
339 List of references. The energy must be the total energy and not
340 energy per atom. The names can also be dicts like
341 ``{'Zn': 1, 'O': 2}`` which would be equivalent to ``'ZnO2'``.
342 filter: str or list of str
343 Use only those references that match the given filter.
344 Example: ``filter='ZnO'`` will select those that
345 contain zinc or oxygen.
346 verbose: bool
347 Write information.
348 """
350 if not references:
351 raise ValueError("You must provide a non-empty list of references"
352 " for the phase diagram! "
353 "You have provided '{}'".format(references))
354 filter = parse_formula(filter)[0]
356 self.verbose = verbose
358 self.species = OrderedDict()
359 self.references = []
360 for name, energy in references:
361 if isinstance(name, str):
362 count = parse_formula(name)[0]
363 else:
364 count = name
366 if filter and any(symbol not in filter for symbol in count):
367 continue
369 if not isinstance(name, str):
370 name = Formula.from_dict(count).format('metal')
372 natoms = 0
373 for symbol, n in count.items():
374 natoms += n
375 if symbol not in self.species:
376 self.species[symbol] = len(self.species)
377 self.references.append((count, energy, name, natoms))
379 ns = len(self.species)
380 self.symbols = [None] * ns
381 for symbol, id in self.species.items():
382 self.symbols[id] = symbol
384 if verbose:
385 print('Species:', ', '.join(self.symbols))
386 print('References:', len(self.references))
387 for i, (count, energy, name, natoms) in enumerate(self.references):
388 print(f'{i:<5}{name:10}{energy:10.3f}')
390 self.points = np.zeros((len(self.references), ns + 1))
391 for s, (count, energy, name, natoms) in enumerate(self.references):
392 for symbol, n in count.items():
393 self.points[s, self.species[symbol]] = n / natoms
394 self.points[s, -1] = energy / natoms
396 if len(self.points) == ns:
397 # Simple case that qhull would choke on:
398 self.simplices = np.arange(ns).reshape((1, ns))
399 self.hull = np.ones(ns, bool)
400 else:
401 hull = ConvexHull(self.points[:, 1:])
403 # Find relevant simplices:
404 ok = hull.equations[:, -2] < 0
405 self.simplices = hull.simplices[ok]
407 # Create a mask for those points that are on the convex hull:
408 self.hull = np.zeros(len(self.points), bool)
409 for simplex in self.simplices:
410 self.hull[simplex] = True
412 if verbose:
413 print('Simplices:', len(self.simplices))
415 def decompose(self, formula=None, **kwargs):
416 """Find the combination of the references with the lowest energy.
418 formula: str
419 Stoichiometry. Example: ``'ZnO'``. Can also be given as
420 keyword arguments: ``decompose(Zn=1, O=1)``.
422 Example::
424 pd = PhaseDiagram(...)
425 pd.decompose(Zn=1, O=3)
427 Returns energy, indices of references and coefficients."""
429 if formula:
430 assert not kwargs
431 kwargs = parse_formula(formula)[0]
433 point = np.zeros(len(self.species))
434 N = 0
435 for symbol, n in kwargs.items():
436 point[self.species[symbol]] = n
437 N += n
439 # Find coordinates within each simplex:
440 X = self.points[self.simplices, 1:-1] - point[1:] / N
442 # Find the simplex with positive coordinates that sum to
443 # less than one:
444 eps = 1e-15
445 for i, Y in enumerate(X):
446 try:
447 x = np.linalg.solve((Y[1:] - Y[:1]).T, -Y[0])
448 except np.linalg.linalg.LinAlgError:
449 continue
450 if (x > -eps).all() and x.sum() < 1 + eps:
451 break
452 else:
453 assert False, X
455 indices = self.simplices[i]
456 points = self.points[indices]
458 scaledcoefs = [1 - x.sum()]
459 scaledcoefs.extend(x)
461 energy = N * np.dot(scaledcoefs, points[:, -1])
463 coefs = []
464 results = []
465 for coef, s in zip(scaledcoefs, indices):
466 count, e, name, natoms = self.references[s]
467 coef *= N / natoms
468 coefs.append(coef)
469 results.append((name, coef, e))
471 if self.verbose:
472 print_results(results)
474 return energy, indices, np.array(coefs)
476 def plot(self, ax=None, dims=None, show=False, **plotkwargs):
477 """Make 2-d or 3-d plot of datapoints and convex hull.
479 Default is 2-d for 2- and 3-component diagrams and 3-d for a
480 4-component diagram.
481 """
482 import matplotlib.pyplot as plt
484 N = len(self.species)
486 if dims is None:
487 if N <= 3:
488 dims = 2
489 else:
490 dims = 3
492 if ax is None:
493 projection = None
494 if dims == 3:
495 projection = '3d'
496 from mpl_toolkits.mplot3d import Axes3D
497 Axes3D # silence pyflakes
498 fig = plt.figure()
499 ax = fig.add_subplot(projection=projection)
500 else:
501 if dims == 3 and not hasattr(ax, 'set_zlim'):
502 raise ValueError('Cannot make 3d plot unless axes projection '
503 'is 3d')
505 if dims == 2:
506 if N == 2:
507 self.plot2d2(ax, **plotkwargs)
508 elif N == 3:
509 self.plot2d3(ax)
510 else:
511 raise ValueError('Can only make 2-d plots for 2 and 3 '
512 'component systems!')
513 else:
514 if N == 3:
515 self.plot3d3(ax)
516 elif N == 4:
517 self.plot3d4(ax)
518 else:
519 raise ValueError('Can only make 3-d plots for 3 and 4 '
520 'component systems!')
521 if show:
522 plt.show()
523 return ax
525 def plot2d2(self, ax=None,
526 only_label_simplices=False, only_plot_simplices=False):
527 x, e = self.points[:, 1:].T
528 names = [re.sub(r'(\d+)', r'$_{\1}$', ref[2])
529 for ref in self.references]
530 hull = self.hull
531 simplices = self.simplices
532 xlabel = self.symbols[1]
533 ylabel = 'energy [eV/atom]'
535 if ax:
536 for i, j in simplices:
537 ax.plot(x[[i, j]], e[[i, j]], '-b')
538 ax.plot(x[hull], e[hull], 'sg')
539 if not only_plot_simplices:
540 ax.plot(x[~hull], e[~hull], 'or')
542 if only_plot_simplices or only_label_simplices:
543 x = x[self.hull]
544 e = e[self.hull]
545 names = [name for name, h in zip(names, self.hull) if h]
546 for a, b, name in zip(x, e, names):
547 ax.text(a, b, name, ha='center', va='top')
549 ax.set_xlabel(xlabel)
550 ax.set_ylabel(ylabel)
552 return (x, e, names, hull, simplices, xlabel, ylabel)
554 def plot2d3(self, ax=None):
555 x, y = self.points[:, 1:-1].T.copy()
556 x += y / 2
557 y *= 3**0.5 / 2
558 names = [re.sub(r'(\d+)', r'$_{\1}$', ref[2])
559 for ref in self.references]
560 hull = self.hull
561 simplices = self.simplices
563 if ax:
564 for i, j, k in simplices:
565 ax.plot(x[[i, j, k, i]], y[[i, j, k, i]], '-b')
566 ax.plot(x[hull], y[hull], 'og')
567 ax.plot(x[~hull], y[~hull], 'sr')
568 for a, b, name in zip(x, y, names):
569 ax.text(a, b, name, ha='center', va='top')
571 return (x, y, names, hull, simplices)
573 def plot3d3(self, ax):
574 x, y, e = self.points[:, 1:].T
576 ax.scatter(x[self.hull], y[self.hull], e[self.hull],
577 c='g', marker='o')
578 ax.scatter(x[~self.hull], y[~self.hull], e[~self.hull],
579 c='r', marker='s')
581 for a, b, c, ref in zip(x, y, e, self.references):
582 name = re.sub(r'(\d+)', r'$_{\1}$', ref[2])
583 ax.text(a, b, c, name, ha='center', va='bottom')
585 for i, j, k in self.simplices:
586 ax.plot(x[[i, j, k, i]],
587 y[[i, j, k, i]],
588 zs=e[[i, j, k, i]], c='b')
590 ax.set_xlim3d(0, 1)
591 ax.set_ylim3d(0, 1)
592 ax.view_init(azim=115, elev=30)
593 ax.set_xlabel(self.symbols[1])
594 ax.set_ylabel(self.symbols[2])
595 ax.set_zlabel('energy [eV/atom]')
597 def plot3d4(self, ax):
598 x, y, z = self.points[:, 1:-1].T
599 a = x / 2 + y + z / 2
600 b = 3**0.5 * (x / 2 + y / 6)
601 c = (2 / 3)**0.5 * z
603 ax.scatter(a[self.hull], b[self.hull], c[self.hull],
604 c='g', marker='o')
605 ax.scatter(a[~self.hull], b[~self.hull], c[~self.hull],
606 c='r', marker='s')
608 for x, y, z, ref in zip(a, b, c, self.references):
609 name = re.sub(r'(\d+)', r'$_{\1}$', ref[2])
610 ax.text(x, y, z, name, ha='center', va='bottom')
612 for i, j, k, w in self.simplices:
613 ax.plot(a[[i, j, k, i, w, k, j, w]],
614 b[[i, j, k, i, w, k, j, w]],
615 zs=c[[i, j, k, i, w, k, j, w]], c='b')
617 ax.set_xlim3d(0, 1)
618 ax.set_ylim3d(0, 1)
619 ax.set_zlim3d(0, 1)
620 ax.view_init(azim=115, elev=30)
623_aqueous = """\
624-525700,SiF6--
625-514100,Rh(SO4)3----
626-504800,Ru(SO4)3----
627-499900,Pd(SO4)3----
628-495200,Ru(SO4)3---
629-485700,H4P2O7
630-483700,Rh(SO4)3---
631-483600,H3P2O7-
632-480400,H2P2O7--
633-480380,Pt(SO4)3----
634-471400,HP2O7---
635-458700,P2O7----
636-447500,LaF4-
637-437600,LaH2PO4++
638-377900,LaF3
639-376299,Ca(HSiO3)+
640-370691,BeF4--
641-355400,BF4-
642-353025,Mg(HSiO3)+
643-346900,LaSO4+
644-334100,Rh(SO4)2--
645-325400,Ru(SO4)2--
646-319640,Pd(SO4)2--
647-317900,Ru(SO4)2-
648-312970,Cr2O7--
649-312930,CaSO4
650-307890,NaHSiO3
651-307800,LaF2+
652-307000,LaHCO3++
653-306100,Rh(SO4)2-
654-302532,BeF3-
655-300670,Pt(SO4)2--
656-299900,LaCO3+
657-289477,MgSO4
658-288400,LaCl4-
659-281500,HZrO3-
660-279200,HHfO3-
661-276720,Sr(HCO3)+
662-275700,Ba(HCO3)+
663-273830,Ca(HCO3)+
664-273100,H3PO4
665-270140,H2PO4-
666-266500,S2O8--
667-264860,Sr(CO3)
668-264860,SrCO3
669-263830,Ba(CO3)
670-263830,BaCO3
671-262850,Ca(CO3)
672-262850,CaCO3
673-260310,HPO4--
674-257600,LaCl3
675-250200,Mg(HCO3)+
676-249200,H3VO4
677-248700,S4O6--
678-246640,KSO4-
679-243990,H2VO4-
680-243500,PO4---
681-243400,KHSO4
682-242801,HSiO3-
683-241700,HYO2
684-241476,NaSO4-
685-239700,HZrO2+
686-239300,LaO2H
687-238760,Mg(CO3)
688-238760,MgCO3
689-237800,HHfO2+
690-236890,Ag(CO3)2---
691-236800,HNbO3
692-236600,LaF++
693-235640,MnSO4
694-233400,ZrO2
695-233000,HVO4--
696-231600,HScO2
697-231540,B(OH)3
698-231400,HfO2
699-231386,BeF2
700-231000,S2O6--
701-229000,S3O6--
702-229000,S5O6--
703-228460,HTiO3-
704-227400,YO2-
705-227100,NbO3-
706-226700,LaCl2+
707-223400,HWO4-
708-221700,LaO2-
709-218500,WO4--
710-218100,ScO2-
711-214900,VO4---
712-210000,YOH++
713-208900,LaOH++
714-207700,HAlO2
715-206400,HMoO4-
716-204800,H3PO3
717-202350,H2PO3-
718-202290,SrF+
719-201807,BaF+
720-201120,BaF+
721-200400,MoO4--
722-200390,CaF+
723-199190,SiO2
724-198693,AlO2-
725-198100,YO+
726-195900,LaO+
727-195800,LaCl++
728-194000,CaCl2
729-194000,HPO3--
730-191300,LaNO3++
731-190400,ZrOH+++
732-189000,HfOH+++
733-189000,S2O5--
734-187600,ZrO++
735-186000,HfO++
736-183700,HCrO4-
737-183600,ScO+
738-183100,H3AsO4
739-180630,HSO4-
740-180010,H2AsO4-
741-177930,SO4--
742-177690,MgF+
743-174800,CrO4--
744-173300,SrOH+
745-172300,BaOH+
746-172200,HBeO2-
747-171300,CaOH+
748-170790,HAsO4--
749-166000,ReO4-
750-165800,SrCl+
751-165475,Al(OH)++
752-165475,AlOH++
753-164730,BaCl+
754-164000,La+++
755-163800,Y+++
756-163100,CaCl+
757-162240,BO2-
758-158493,BeF+
759-158188,AlO+
760-155700,VOOH+
761-155164,CdF2
762-154970,AsO4---
763-153500,Rh(SO4)
764-152900,BeO2--
765-152370,HSO5-
766-151540,RuCl6---
767-149255,MgOH+
768-147400,H2S2O4
769-146900,HS2O4-
770-146081,CdCl4--
771-145521,BeCl2
772-145200,Ru(SO4)
773-145056,PbF2
774-143500,S2O4--
775-140330,H2AsO3-
776-140300,VO2+
777-140282,HCO3-
778-140200,Sc+++
779-139900,BeOH+
780-139700,MgCl+
781-139200,Ru(SO4)+
782-139000,Pd(SO4)
783-138160,HF2-
784-138100,HCrO2
785-138000,TiO++
786-137300,HGaO2
787-136450,RbF
788-134760,Sr++
789-134030,Ba++
790-133270,Zr++++
791-133177,PbCl4--
792-132600,Hf++++
793-132120,Ca++
794-129310,ZnCl3-
795-128700,GaO2-
796-128600,BeO
797-128570,NaF
798-128000,H2S2O3
799-127500,Rh(SO4)+
800-127200,HS2O3-
801-126191,CO3--
802-126130,HSO3-
803-125300,CrO2-
804-125100,H3PO2
805-124900,S2O3--
806-123641,MnF+
807-122400,H2PO2-
808-121000,HMnO2-
809-120700,RuCl5--
810-120400,MnO4--
811-120300,Pt(SO4)
812-119800,HInO2
813-116300,SO3--
814-115971,CdCl3-
815-115609,Al+++
816-115316,BeCl+
817-112280,AgCl4---
818-111670,TiO2++
819-111500,VOH++
820-111430,Ag(CO3)-
821-110720,HZnO2-
822-108505,Mg++
823-108100,HSeO4-
824-108000,LiOH
825-107600,MnO4-
826-106988,HgCl4--
827-106700,InO2-
828-106700,VO++
829-106100,VO+
830-105500,SeO4--
831-105100,RbOH
832-105000,CsOH
833-104500,KOH
834-104109,ZnF+
835-103900,PdCl4--
836-103579,CuCl4--
837-102600,MnO2--
838-102150,PbCl3-
839-101850,H2SeO3
840-101100,HFeO2
841-100900,CsCl
842-100500,CrOH++
843-99900,NaOH
844-99800,VOH+
845-99250,LiCl
846-98340,HSeO3-
847-98300,ZnCl2
848-97870,RbCl
849-97400,HSbO2
850-97300,HSnO2-
851-97300,MnOH+
852-97016,InF++
853-96240,HAsO2
854-95430,KCl
855-95400,HFeO2-
856-94610,CsBr
857-93290,ZnO2--
858-93250,RhCl4--
859-92910,NaCl
860-92800,CrO+
861-92250,CO2
862-91210,PtCl4--
863-91157,FeF+
864-91100,GaOH++
865-91010,RbBr
866-90550,Be++
867-90010,KBr
868-89963,CuCl3--
869-89730,RuCl4-
870-88400,SeO3--
871-88000,FeO2-
872-87373,CdF+
873-86600,GaO+
874-86500,HCdO2-
875-86290,MnCl+
876-85610,NaBr
877-84851,CdCl2
878-83900,RuCl4--
879-83650,AsO2-
880-83600,Ti+++
881-83460,CsI
882-83400,HCoO2-
883-82710,AgCl3--
884-82400,SbO2-
885-81980,HNiO2-
886-81732,CoF+
887-81500,MnO
888-81190,ZnOH+
889-81000,HPbO2-
890-79768,NiF+
891-79645,FeF++
892-79300,HBiO2
893-78900,RbI
894-77740,KI
895-77700,La++
896-77500,RhCl4-
897-75860,PbF+
898-75338,CuCl3-
899-75216,TlF
900-75100,Ti++
901-74600,InOH++
902-74504,HgCl3-
903-73480,FeCl2
904-72900,NaI
905-71980,SO2
906-71662,HF
907-71600,RuO4--
908-71200,PbCl2
909-69933,Li+
910-69810,PdCl3-
911-69710,Cs+
912-69400,InO+
913-67811,AuCl3--
914-67800,Rb+
915-67510,K+
916-67420,ZnO
917-67340,F-
918-67300,CdO2--
919-66850,ZnCl+
920-65850,FeOH+
921-65550,TlOH
922-64200,NiO2--
923-63530,RhCl3-
924-63200,CoO2--
925-62591,Na+
926-61700,BiO2-
927-61500,CdOH+
928-60100,HCuO2-
929-59226,InCl++
930-58600,SnOH+
931-58560,RuCl3
932-58038,CuCl2-
933-57900,V+++
934-57800,FeOH++
935-57760,PtCl3-
936-57600,HTlO2
937-56690,H2O
938-56025,CoOH+
939-55100,Mn++
940-54380,RuCl3-
941-53950,PbOH+
942-53739,CuF+
943-53600,SnO
944-53100,FeO+
945-53030,FeCl+
946-52850,NiOH+
947-52627,CdCl+
948-52000,V++
949-51560,AgCl2-
950-50720,FeO
951-49459,AgF
952-49300,Cr+++
953-47500,CdO
954-46190,RhCl3
955-46142,CuCl2
956-45200,HHgO2-
957-45157,CoCl+
958-44000,CoO
959-42838,HgCl2
960-41600,TlO2-
961-41200,CuO2--
962-40920,NiCl+
963-39815,TlCl
964-39400,Cr++
965-39350,PbO
966-39340,NiO
967-39050,PbCl+
968-38000,Ga+++
969-37518,FeCl++
970-36781,AuCl2-
971-35332,AuCl4-
972-35200,Zn++
973-35160,PdCl2
974-33970,RhCl2
975-32300,BiOH++
976-31700,HIO3
977-31379,Cl-
978-30600,IO3-
979-30410,HCl
980-30204,HgF+
981-30200,CuOH+
982-29300,BiO+
983-28682,CO
984-26507,NO3-
985-26440,RuCl2+
986-25590,Br3-
987-25060,RuCl2
988-24870,Br-
989-24730,HNO3
990-23700,HIO
991-23400,In+++
992-23280,OCN-
993-23000,CoOH++
994-22608,CuCl
995-22290,PtCl2
996-21900,AgOH
997-21870,Fe++
998-20800,CuO
999-20300,Mn+++
1000-20058,Pb(HS)2
1001-19700,HBrO
1002-19100,HClO
1003-19100,ScOH++
1004-18990,NH4+
1005-18971,Pb(HS)3-
1006-18560,Cd++
1007-18290,Rh(OH)+
1008-17450,AgCl
1009-16250,CuCl+
1010-14780,RhCl2+
1011-14000,IO4-
1012-13130,Pd(OH)+
1013-13000,Co++
1014-12700,HgOH+
1015-12410,I-
1016-12300,I3-
1017-12190,Ru(OH)2++
1018-12100,HNO2
1019-11500,PdO
1020-10900,Ni++
1021-10470,Ru(OH)+
1022-10450,RuO+
1023-9200,IO-
1024-8900,HgO
1025-8800,ClO-
1026-8000,BrO-
1027-7740,Tl+
1028-7738,AgNO3
1029-7700,NO2-
1030-7220,RhO
1031-6673,H2S
1032-6570,Sn++
1033-6383,NH3
1034-5710,Pb++
1035-5500,AgO-
1036-4500,TlOH++
1037-4120,Fe+++
1038-3380,RhCl+
1039-3200,TlO+
1040-3184,AuCl
1041-2155,HgCl+
1042-2040,ClO4-
1043-1900,ClO3-
1044-1130,PtO
1045-820,Rh(OH)++
10460,Ag(HS)2-
10470,H+
1048230,RuO
10491400,HClO2
10501560,Pt(OH)+
10512429,Au(HS)2-
10522500,PdCl+
10532860,HS-
10543140,RhO+
10553215,Xe
10563554,Kr
10573890,Ar
10584100,ClO2-
10594347,N2
10604450,BrO3-
10614565,Ne
10624658,He
10635210,RuCl+
10647100,RuCl++
10658600,H2N2O2
10669375,TlCl++
106710500,HSe-
106811950,Cu+
106915675,Cu++
107015700,S5--
107116500,S4--
107217600,S3--
107318200,HN2O2-
107418330,RhCl++
107518380,PtCl+
107618427,Ag+
107719000,S2--
107819500,SeCN-
107919700,N2H5+
108021100,N2H6++
108122160,SCN-
108222880,Bi+++
108327700,Rh++
108428200,BrO4-
108528600,HCN
108632000,Co+++
108733200,N2O2--
108835900,Ru++
108936710,Hg2++
109039360,Hg++
109141200,CN-
109241440,Ru+++
109342200,Pd++
109451300,Tl+++
109552450,Rh+++
109661600,Pt++
109764300,Ag++
1098103600,Au+++"""