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