Coverage for /builds/debichem-team/python-ase/ase/phasediagram.py: 68.38%

370 statements  

« prev     ^ index     » next       coverage.py v7.5.3, created at 2025-03-06 04:00 +0000

1import fractions 

2import functools 

3import re 

4from collections import OrderedDict 

5from typing import Dict, List, Tuple 

6 

7import numpy as np 

8from scipy.spatial import ConvexHull 

9 

10import ase.units as units 

11from ase.formula import Formula 

12from ase.utils import deprecated 

13 

14_solvated: List[Tuple[str, Dict[str, int], float, bool, float]] = [] 

15 

16 

17def parse_formula(formula): 

18 aq = formula.endswith('(aq)') 

19 if aq: 

20 formula = formula[:-4] 

21 charge = formula.count('+') - formula.count('-') 

22 if charge: 

23 formula = formula.rstrip('+-') 

24 count = Formula(formula).count() 

25 return count, charge, aq 

26 

27 

28def float2str(x): 

29 f = fractions.Fraction(x).limit_denominator(100) 

30 n = f.numerator 

31 d = f.denominator 

32 if abs(n / d - f) > 1e-6: 

33 return f'{f:.3f}' 

34 if d == 0: 

35 return '0' 

36 if f.denominator == 1: 

37 return str(n) 

38 return f'{f.numerator}/{f.denominator}' 

39 

40 

41def solvated(symbols): 

42 """Extract solvation energies from database. 

43 

44 symbols: str 

45 Extract only those molecules that contain the chemical elements 

46 given by the symbols string (plus water and H+). 

47 

48 Data from: 

49 

50 Johnson JW, Oelkers EH, Helgeson HC (1992) 

51 Comput Geosci 18(7):899. 

52 doi:10.1016/0098-3004(92)90029-Q 

53 

54 and: 

55 

56 Pourbaix M (1966) 

57 Atlas of electrochemical equilibria in aqueous solutions. 

58 No. v. 1 in Atlas of Electrochemical Equilibria in Aqueous Solutions. 

59 Pergamon Press, New York. 

60 

61 Returns list of (name, energy) tuples. 

62 """ 

63 

64 if isinstance(symbols, str): 

65 symbols = Formula(symbols).count().keys() 

66 if len(_solvated) == 0: 

67 for line in _aqueous.splitlines(): 

68 energy, formula = line.split(',') 

69 name = formula + '(aq)' 

70 count, charge, aq = parse_formula(name) 

71 energy = float(energy) * 0.001 * units.kcal / units.mol 

72 _solvated.append((name, count, charge, aq, energy)) 

73 references = [] 

74 for name, count, charge, aq, energy in _solvated: 

75 for symbol in count: 

76 if symbol not in 'HO' and symbol not in symbols: 

77 break 

78 else: 

79 references.append((name, energy)) 

80 return references 

81 

82 

83def bisect(A, X, Y, f): 

84 a = [] 

85 for i in [0, -1]: 

86 for j in [0, -1]: 

87 if A[i, j] == -1: 

88 A[i, j] = f(X[i], Y[j]) 

89 a.append(A[i, j]) 

90 

91 if np.ptp(a) == 0: 

92 A[:] = a[0] 

93 return 

94 if a[0] == a[1]: 

95 A[0] = a[0] 

96 if a[1] == a[3]: 

97 A[:, -1] = a[1] 

98 if a[3] == a[2]: 

99 A[-1] = a[3] 

100 if a[2] == a[0]: 

101 A[:, 0] = a[2] 

102 if not (A == -1).any(): 

103 return 

104 i = len(X) // 2 

105 j = len(Y) // 2 

106 bisect(A[:i + 1, :j + 1], X[:i + 1], Y[:j + 1], f) 

107 bisect(A[:i + 1, j:], X[:i + 1], Y[j:], f) 

108 bisect(A[i:, :j + 1], X[i:], Y[:j + 1], f) 

109 bisect(A[i:, j:], X[i:], Y[j:], f) 

110 

111 

112def print_results(results): 

113 total_energy = 0.0 

114 print('reference coefficient energy') 

115 print('------------------------------------') 

116 for name, coef, energy in results: 

117 total_energy += coef * energy 

118 if abs(coef) < 1e-7: 

119 continue 

120 print(f'{name:14}{float2str(coef):>10}{energy:12.3f}') 

121 print('------------------------------------') 

122 print(f'Total energy: {total_energy:22.3f}') 

123 print('------------------------------------') 

124 

125 

126class Pourbaix: 

127 @deprecated( 

128 'Use ase.pourbaix.Pourbaix. ' 

129 'This class will be removed in a future version of ASE.') 

130 def __init__(self, references, formula=None, T=300.0, **kwargs): 

131 """Pourbaix object. 

132 

133 references: list of (name, energy) tuples 

134 Examples of names: ZnO2, H+(aq), H2O(aq), Zn++(aq), ... 

135 formula: str 

136 Stoichiometry. Example: ``'ZnO'``. Can also be given as 

137 keyword arguments: ``Pourbaix(refs, Zn=1, O=1)``. 

138 T: float 

139 Temperature in Kelvin. 

140 

141 .. deprecated:: 3.24.0 

142 """ 

143 

144 if formula: 

145 assert not kwargs 

146 kwargs = parse_formula(formula)[0] 

147 

148 if 'O' not in kwargs: 

149 kwargs['O'] = 0 

150 if 'H' not in kwargs: 

151 kwargs['H'] = 0 

152 

153 self.kT = units.kB * T 

154 self.references = [] 

155 for name, energy in references: 

156 if name == 'O': 

157 continue 

158 count, charge, aq = parse_formula(name) 

159 if all(symbol in kwargs for symbol in count): 

160 self.references.append((count, charge, aq, energy, name)) 

161 

162 self.references.append(({}, -1, False, 0.0, 'e-')) # an electron 

163 

164 self.count = kwargs 

165 

166 self.N = {'e-': 0} 

167 for symbol in kwargs: 

168 if symbol not in self.N: 

169 self.N[symbol] = len(self.N) 

170 

171 def decompose(self, U, pH, verbose=True, concentration=1e-6): 

172 """Decompose material. 

173 

174 U: float 

175 Potential in V. 

176 pH: float 

177 pH value. 

178 verbose: bool 

179 Default is True. 

180 concentration: float 

181 Concentration of solvated references. 

182 

183 Returns optimal coefficients and energy: 

184 

185 >>> from ase.phasediagram import Pourbaix, solvated 

186 >>> refs = solvated('CoO') + [ 

187 ... ('Co', 0.0), 

188 ... ('CoO', -2.509), 

189 ... ('Co3O4', -9.402)] 

190 >>> pb = Pourbaix(refs, Co=3, O=4) 

191 >>> coefs, energy = pb.decompose(U=1.5, pH=0, 

192 ... concentration=1e-6, 

193 ... verbose=True) 

194 0 HCoO2-(aq) -3.974 

195 1 CoO2--(aq) -3.098 

196 2 H2O(aq) -2.458 

197 3 CoOH+(aq) -2.787 

198 4 CoO(aq) -2.265 

199 5 CoOH++(aq) -1.355 

200 6 Co++(aq) -0.921 

201 7 H+(aq) 0.000 

202 8 Co+++(aq) 1.030 

203 9 Co 0.000 

204 10 CoO -2.509 

205 11 Co3O4 -9.402 

206 12 e- -1.500 

207 reference coefficient energy 

208 ------------------------------------ 

209 H2O(aq) 4 -2.458 

210 Co++(aq) 3 -0.921 

211 H+(aq) -8 0.000 

212 e- -2 -1.500 

213 ------------------------------------ 

214 Total energy: -9.596 

215 ------------------------------------ 

216 """ 

217 

218 alpha = np.log(10) * self.kT 

219 entropy = -np.log(concentration) * self.kT 

220 

221 # We want to minimize np.dot(energies, x) under the constraints: 

222 # 

223 # np.dot(x, eq2) == eq1 

224 # 

225 # with bounds[i,0] <= x[i] <= bounds[i, 1]. 

226 # 

227 # First two equations are charge and number of hydrogens, and 

228 # the rest are the remaining species. 

229 

230 eq1 = [0] + list(self.count.values()) 

231 eq2 = [] 

232 energies = [] 

233 bounds = [] 

234 names = [] 

235 for count, charge, aq, energy, name in self.references: 

236 eq = np.zeros(len(self.N)) 

237 eq[0] = charge 

238 for symbol, n in count.items(): 

239 eq[self.N[symbol]] = n 

240 eq2.append(eq) 

241 if name in ['H2O(aq)', 'H+(aq)', 'e-']: 

242 bounds.append((-np.inf, np.inf)) 

243 if name == 'e-': 

244 energy = -U 

245 elif name == 'H+(aq)': 

246 energy = -pH * alpha 

247 else: 

248 bounds.append((0, np.inf)) 

249 if aq: 

250 energy -= entropy 

251 if verbose: 

252 print('{:<5}{:10}{:10.3f}'.format(len(energies), 

253 name, energy)) 

254 energies.append(energy) 

255 names.append(name) 

256 

257 from scipy.optimize import linprog 

258 

259 result = linprog(c=energies, 

260 A_eq=np.transpose(eq2), 

261 b_eq=eq1, 

262 bounds=bounds) 

263 

264 if verbose: 

265 print_results(zip(names, result.x, energies)) 

266 

267 return result.x, result.fun 

268 

269 def diagram(self, U, pH, plot=True, show=False, ax=None): 

270 """Calculate Pourbaix diagram. 

271 

272 U: list of float 

273 Potentials in V. 

274 pH: list of float 

275 pH values. 

276 plot: bool 

277 Create plot. 

278 show: bool 

279 Open graphical window and show plot. 

280 ax: matplotlib axes object 

281 When creating plot, plot onto the given axes object. 

282 If none given, plot onto the current one. 

283 """ 

284 a = np.empty((len(U), len(pH)), int) 

285 a[:] = -1 

286 colors = {} 

287 f = functools.partial(self.colorfunction, colors=colors) 

288 bisect(a, U, pH, f) 

289 compositions = [None] * len(colors) 

290 names = [ref[-1] for ref in self.references] 

291 for indices, color in colors.items(): 

292 compositions[color] = ' + '.join(names[i] for i in indices 

293 if names[i] not in 

294 ['H2O(aq)', 'H+(aq)', 'e-']) 

295 text = [] 

296 for i, name in enumerate(compositions): 

297 b = (a == i) 

298 x = np.dot(b.sum(1), U) / b.sum() 

299 y = np.dot(b.sum(0), pH) / b.sum() 

300 name = re.sub(r'(\S)([+-]+)', r'\1$^{\2}$', name) 

301 name = re.sub(r'(\d+)', r'$_{\1}$', name) 

302 text.append((x, y, name)) 

303 

304 if plot: 

305 import matplotlib.cm as cm 

306 import matplotlib.pyplot as plt 

307 if ax is None: 

308 ax = plt.gca() 

309 

310 # rasterized pcolormesh has a bug which leaves a tiny 

311 # white border. Unrasterized pcolormesh produces 

312 # unreasonably large files. Avoid this by using the more 

313 # general imshow. 

314 ax.imshow(a, cmap=cm.Accent, 

315 extent=[min(pH), max(pH), min(U), max(U)], 

316 origin='lower', 

317 aspect='auto') 

318 

319 for x, y, name in text: 

320 ax.text(y, x, name, horizontalalignment='center') 

321 ax.set_xlabel('pH') 

322 ax.set_ylabel('potential [V]') 

323 ax.set_xlim(min(pH), max(pH)) 

324 ax.set_ylim(min(U), max(U)) 

325 if show: 

326 plt.show() 

327 

328 return a, compositions, text 

329 

330 def colorfunction(self, U, pH, colors): 

331 coefs, _energy = self.decompose(U, pH, verbose=False) 

332 indices = tuple(sorted(np.where(abs(coefs) > 1e-3)[0])) 

333 color = colors.get(indices) 

334 if color is None: 

335 color = len(colors) 

336 colors[indices] = color 

337 return color 

338 

339 

340class PhaseDiagram: 

341 def __init__(self, references, filter='', verbose=True): 

342 """Phase-diagram. 

343 

344 references: list of (name, energy) tuples 

345 List of references. The energy must be the total energy and not 

346 energy per atom. The names can also be dicts like 

347 ``{'Zn': 1, 'O': 2}`` which would be equivalent to ``'ZnO2'``. 

348 filter: str or list of str 

349 Use only those references that match the given filter. 

350 Example: ``filter='ZnO'`` will select those that 

351 contain zinc or oxygen. 

352 verbose: bool 

353 Write information. 

354 """ 

355 

356 if not references: 

357 raise ValueError("You must provide a non-empty list of references" 

358 " for the phase diagram! " 

359 "You have provided '{}'".format(references)) 

360 filter = parse_formula(filter)[0] 

361 

362 self.verbose = verbose 

363 

364 self.species = OrderedDict() 

365 self.references = [] 

366 for name, energy in references: 

367 if isinstance(name, str): 

368 count = parse_formula(name)[0] 

369 else: 

370 count = name 

371 

372 if filter and any(symbol not in filter for symbol in count): 

373 continue 

374 

375 if not isinstance(name, str): 

376 name = Formula.from_dict(count).format('metal') 

377 

378 natoms = 0 

379 for symbol, n in count.items(): 

380 natoms += n 

381 if symbol not in self.species: 

382 self.species[symbol] = len(self.species) 

383 self.references.append((count, energy, name, natoms)) 

384 

385 ns = len(self.species) 

386 self.symbols = [None] * ns 

387 for symbol, id in self.species.items(): 

388 self.symbols[id] = symbol 

389 

390 if verbose: 

391 print('Species:', ', '.join(self.symbols)) 

392 print('References:', len(self.references)) 

393 for i, (count, energy, name, natoms) in enumerate(self.references): 

394 print(f'{i:<5}{name:10}{energy:10.3f}') 

395 

396 self.points = np.zeros((len(self.references), ns + 1)) 

397 for s, (count, energy, name, natoms) in enumerate(self.references): 

398 for symbol, n in count.items(): 

399 self.points[s, self.species[symbol]] = n / natoms 

400 self.points[s, -1] = energy / natoms 

401 

402 if len(self.points) == ns: 

403 # Simple case that qhull would choke on: 

404 self.simplices = np.arange(ns).reshape((1, ns)) 

405 self.hull = np.ones(ns, bool) 

406 elif ns == 1: 

407 # qhull also doesn't like ns=1: 

408 i = self.points[:, 1].argmin() 

409 self.simplices = np.array([[i]]) 

410 self.hull = np.zeros(len(self.points), bool) 

411 self.hull[i] = True 

412 else: 

413 hull = ConvexHull(self.points[:, 1:]) 

414 

415 # Find relevant simplices: 

416 ok = hull.equations[:, -2] < 0 

417 self.simplices = hull.simplices[ok] 

418 

419 # Create a mask for those points that are on the convex hull: 

420 self.hull = np.zeros(len(self.points), bool) 

421 for simplex in self.simplices: 

422 self.hull[simplex] = True 

423 

424 if verbose: 

425 print('Simplices:', len(self.simplices)) 

426 

427 def decompose(self, formula=None, **kwargs): 

428 """Find the combination of the references with the lowest energy. 

429 

430 formula: str 

431 Stoichiometry. Example: ``'ZnO'``. Can also be given as 

432 keyword arguments: ``decompose(Zn=1, O=1)``. 

433 

434 Example:: 

435 

436 pd = PhaseDiagram(...) 

437 pd.decompose(Zn=1, O=3) 

438 

439 Returns energy, indices of references and coefficients.""" 

440 

441 if formula: 

442 assert not kwargs 

443 kwargs = parse_formula(formula)[0] 

444 

445 point = np.zeros(len(self.species)) 

446 N = 0 

447 for symbol, n in kwargs.items(): 

448 point[self.species[symbol]] = n 

449 N += n 

450 

451 # Find coordinates within each simplex: 

452 X = self.points[self.simplices, 1:-1] - point[1:] / N 

453 

454 # Find the simplex with positive coordinates that sum to 

455 # less than one: 

456 eps = 1e-14 

457 candidates = [] 

458 for i, Y in enumerate(X): 

459 try: 

460 x = np.linalg.solve((Y[1:] - Y[:1]).T, -Y[0]) 

461 except np.linalg.linalg.LinAlgError: 

462 continue 

463 if (x > -eps).all() and x.sum() < 1 + eps: 

464 indices = self.simplices[i] 

465 points = self.points[indices] 

466 

467 scaledcoefs = [1 - x.sum()] 

468 scaledcoefs.extend(x) 

469 

470 energy = N * np.dot(scaledcoefs, points[:, -1]) 

471 candidates.append((energy, indices, points, scaledcoefs)) 

472 

473 # Pick the one with lowest energy: 

474 energy, indices, points, scaledcoefs = min( 

475 candidates, key=lambda x: x[0]) 

476 

477 coefs = [] 

478 results = [] 

479 for coef, s in zip(scaledcoefs, indices): 

480 _count, e, name, natoms = self.references[s] 

481 coef *= N / natoms 

482 coefs.append(coef) 

483 results.append((name, coef, e)) 

484 

485 if self.verbose: 

486 print_results(results) 

487 

488 return energy, indices, np.array(coefs) 

489 

490 def plot(self, ax=None, dims=None, show=False, **plotkwargs): 

491 """Make 2-d or 3-d plot of datapoints and convex hull. 

492 

493 Default is 2-d for 2- and 3-component diagrams and 3-d for a 

494 4-component diagram. 

495 """ 

496 import matplotlib.pyplot as plt 

497 

498 N = len(self.species) 

499 

500 if dims is None: 

501 if N <= 3: 

502 dims = 2 

503 else: 

504 dims = 3 

505 

506 if ax is None: 

507 projection = None 

508 if dims == 3: 

509 projection = '3d' 

510 from mpl_toolkits.mplot3d import Axes3D 

511 Axes3D # silence pyflakes 

512 fig = plt.figure() 

513 ax = fig.add_subplot(projection=projection) 

514 else: 

515 if dims == 3 and not hasattr(ax, 'set_zlim'): 

516 raise ValueError('Cannot make 3d plot unless axes projection ' 

517 'is 3d') 

518 

519 if dims == 2: 

520 if N == 2: 

521 self.plot2d2(ax, **plotkwargs) 

522 elif N == 3: 

523 self.plot2d3(ax) 

524 else: 

525 raise ValueError('Can only make 2-d plots for 2 and 3 ' 

526 'component systems!') 

527 else: 

528 if N == 3: 

529 self.plot3d3(ax) 

530 elif N == 4: 

531 self.plot3d4(ax) 

532 else: 

533 raise ValueError('Can only make 3-d plots for 3 and 4 ' 

534 'component systems!') 

535 if show: 

536 plt.show() 

537 return ax 

538 

539 def plot2d2(self, ax=None, 

540 only_label_simplices=False, only_plot_simplices=False): 

541 x, e = self.points[:, 1:].T 

542 names = [re.sub(r'(\d+)', r'$_{\1}$', ref[2]) 

543 for ref in self.references] 

544 hull = self.hull 

545 simplices = self.simplices 

546 xlabel = self.symbols[1] 

547 ylabel = 'energy [eV/atom]' 

548 

549 if ax: 

550 for i, j in simplices: 

551 ax.plot(x[[i, j]], e[[i, j]], '-b') 

552 ax.plot(x[hull], e[hull], 'sg') 

553 if not only_plot_simplices: 

554 ax.plot(x[~hull], e[~hull], 'or') 

555 

556 if only_plot_simplices or only_label_simplices: 

557 x = x[self.hull] 

558 e = e[self.hull] 

559 names = [name for name, h in zip(names, self.hull) if h] 

560 for a, b, name in zip(x, e, names): 

561 ax.text(a, b, name, ha='center', va='top') 

562 

563 ax.set_xlabel(xlabel) 

564 ax.set_ylabel(ylabel) 

565 

566 return (x, e, names, hull, simplices, xlabel, ylabel) 

567 

568 def plot2d3(self, ax=None): 

569 x, y = self.points[:, 1:-1].T.copy() 

570 x += y / 2 

571 y *= 3**0.5 / 2 

572 names = [re.sub(r'(\d+)', r'$_{\1}$', ref[2]) 

573 for ref in self.references] 

574 hull = self.hull 

575 simplices = self.simplices 

576 

577 if ax: 

578 for i, j, k in simplices: 

579 ax.plot(x[[i, j, k, i]], y[[i, j, k, i]], '-b') 

580 ax.plot(x[hull], y[hull], 'og') 

581 ax.plot(x[~hull], y[~hull], 'sr') 

582 for a, b, name in zip(x, y, names): 

583 ax.text(a, b, name, ha='center', va='top') 

584 

585 return (x, y, names, hull, simplices) 

586 

587 def plot3d3(self, ax): 

588 x, y, e = self.points[:, 1:].T 

589 

590 ax.scatter(x[self.hull], y[self.hull], e[self.hull], 

591 c='g', marker='o') 

592 ax.scatter(x[~self.hull], y[~self.hull], e[~self.hull], 

593 c='r', marker='s') 

594 

595 for a, b, c, ref in zip(x, y, e, self.references): 

596 name = re.sub(r'(\d+)', r'$_{\1}$', ref[2]) 

597 ax.text(a, b, c, name, ha='center', va='bottom') 

598 

599 for i, j, k in self.simplices: 

600 ax.plot(x[[i, j, k, i]], 

601 y[[i, j, k, i]], 

602 zs=e[[i, j, k, i]], c='b') 

603 

604 ax.set_xlim3d(0, 1) 

605 ax.set_ylim3d(0, 1) 

606 ax.view_init(azim=115, elev=30) 

607 ax.set_xlabel(self.symbols[1]) 

608 ax.set_ylabel(self.symbols[2]) 

609 ax.set_zlabel('energy [eV/atom]') 

610 

611 def plot3d4(self, ax): 

612 x, y, z = self.points[:, 1:-1].T 

613 a = x / 2 + y + z / 2 

614 b = 3**0.5 * (x / 2 + y / 6) 

615 c = (2 / 3)**0.5 * z 

616 

617 ax.scatter(a[self.hull], b[self.hull], c[self.hull], 

618 c='g', marker='o') 

619 ax.scatter(a[~self.hull], b[~self.hull], c[~self.hull], 

620 c='r', marker='s') 

621 

622 for x, y, z, ref in zip(a, b, c, self.references): 

623 name = re.sub(r'(\d+)', r'$_{\1}$', ref[2]) 

624 ax.text(x, y, z, name, ha='center', va='bottom') 

625 

626 for i, j, k, w in self.simplices: 

627 ax.plot(a[[i, j, k, i, w, k, j, w]], 

628 b[[i, j, k, i, w, k, j, w]], 

629 zs=c[[i, j, k, i, w, k, j, w]], c='b') 

630 

631 ax.set_xlim3d(0, 1) 

632 ax.set_ylim3d(0, 1) 

633 ax.set_zlim3d(0, 1) 

634 ax.view_init(azim=115, elev=30) 

635 

636 

637_aqueous = """\ 

638-525700,SiF6-- 

639-514100,Rh(SO4)3---- 

640-504800,Ru(SO4)3---- 

641-499900,Pd(SO4)3---- 

642-495200,Ru(SO4)3--- 

643-485700,H4P2O7 

644-483700,Rh(SO4)3--- 

645-483600,H3P2O7- 

646-480400,H2P2O7-- 

647-480380,Pt(SO4)3---- 

648-471400,HP2O7--- 

649-458700,P2O7---- 

650-447500,LaF4- 

651-437600,LaH2PO4++ 

652-377900,LaF3 

653-376299,Ca(HSiO3)+ 

654-370691,BeF4-- 

655-355400,BF4- 

656-353025,Mg(HSiO3)+ 

657-346900,LaSO4+ 

658-334100,Rh(SO4)2-- 

659-325400,Ru(SO4)2-- 

660-319640,Pd(SO4)2-- 

661-317900,Ru(SO4)2- 

662-312970,Cr2O7-- 

663-312930,CaSO4 

664-307890,NaHSiO3 

665-307800,LaF2+ 

666-307000,LaHCO3++ 

667-306100,Rh(SO4)2- 

668-302532,BeF3- 

669-300670,Pt(SO4)2-- 

670-299900,LaCO3+ 

671-289477,MgSO4 

672-288400,LaCl4- 

673-281500,HZrO3- 

674-279200,HHfO3- 

675-276720,Sr(HCO3)+ 

676-275700,Ba(HCO3)+ 

677-273830,Ca(HCO3)+ 

678-273100,H3PO4 

679-270140,H2PO4- 

680-266500,S2O8-- 

681-264860,Sr(CO3) 

682-264860,SrCO3 

683-263830,Ba(CO3) 

684-263830,BaCO3 

685-262850,Ca(CO3) 

686-262850,CaCO3 

687-260310,HPO4-- 

688-257600,LaCl3 

689-250200,Mg(HCO3)+ 

690-249200,H3VO4 

691-248700,S4O6-- 

692-246640,KSO4- 

693-243990,H2VO4- 

694-243500,PO4--- 

695-243400,KHSO4 

696-242801,HSiO3- 

697-241700,HYO2 

698-241476,NaSO4- 

699-239700,HZrO2+ 

700-239300,LaO2H 

701-238760,Mg(CO3) 

702-238760,MgCO3 

703-237800,HHfO2+ 

704-236890,Ag(CO3)2--- 

705-236800,HNbO3 

706-236600,LaF++ 

707-235640,MnSO4 

708-233400,ZrO2 

709-233000,HVO4-- 

710-231600,HScO2 

711-231540,B(OH)3 

712-231400,HfO2 

713-231386,BeF2 

714-231000,S2O6-- 

715-229000,S3O6-- 

716-229000,S5O6-- 

717-228460,HTiO3- 

718-227400,YO2- 

719-227100,NbO3- 

720-226700,LaCl2+ 

721-223400,HWO4- 

722-221700,LaO2- 

723-218500,WO4-- 

724-218100,ScO2- 

725-214900,VO4--- 

726-210000,YOH++ 

727-208900,LaOH++ 

728-207700,HAlO2 

729-206400,HMoO4- 

730-204800,H3PO3 

731-202350,H2PO3- 

732-202290,SrF+ 

733-201807,BaF+ 

734-201120,BaF+ 

735-200400,MoO4-- 

736-200390,CaF+ 

737-199190,SiO2 

738-198693,AlO2- 

739-198100,YO+ 

740-195900,LaO+ 

741-195800,LaCl++ 

742-194000,CaCl2 

743-194000,HPO3-- 

744-191300,LaNO3++ 

745-190400,ZrOH+++ 

746-189000,HfOH+++ 

747-189000,S2O5-- 

748-187600,ZrO++ 

749-186000,HfO++ 

750-183700,HCrO4- 

751-183600,ScO+ 

752-183100,H3AsO4 

753-180630,HSO4- 

754-180010,H2AsO4- 

755-177930,SO4-- 

756-177690,MgF+ 

757-174800,CrO4-- 

758-173300,SrOH+ 

759-172300,BaOH+ 

760-172200,HBeO2- 

761-171300,CaOH+ 

762-170790,HAsO4-- 

763-166000,ReO4- 

764-165800,SrCl+ 

765-165475,Al(OH)++ 

766-165475,AlOH++ 

767-164730,BaCl+ 

768-164000,La+++ 

769-163800,Y+++ 

770-163100,CaCl+ 

771-162240,BO2- 

772-158493,BeF+ 

773-158188,AlO+ 

774-155700,VOOH+ 

775-155164,CdF2 

776-154970,AsO4--- 

777-153500,Rh(SO4) 

778-152900,BeO2-- 

779-152370,HSO5- 

780-151540,RuCl6--- 

781-149255,MgOH+ 

782-147400,H2S2O4 

783-146900,HS2O4- 

784-146081,CdCl4-- 

785-145521,BeCl2 

786-145200,Ru(SO4) 

787-145056,PbF2 

788-143500,S2O4-- 

789-140330,H2AsO3- 

790-140300,VO2+ 

791-140282,HCO3- 

792-140200,Sc+++ 

793-139900,BeOH+ 

794-139700,MgCl+ 

795-139200,Ru(SO4)+ 

796-139000,Pd(SO4) 

797-138160,HF2- 

798-138100,HCrO2 

799-138000,TiO++ 

800-137300,HGaO2 

801-136450,RbF 

802-134760,Sr++ 

803-134030,Ba++ 

804-133270,Zr++++ 

805-133177,PbCl4-- 

806-132600,Hf++++ 

807-132120,Ca++ 

808-129310,ZnCl3- 

809-128700,GaO2- 

810-128600,BeO 

811-128570,NaF 

812-128000,H2S2O3 

813-127500,Rh(SO4)+ 

814-127200,HS2O3- 

815-126191,CO3-- 

816-126130,HSO3- 

817-125300,CrO2- 

818-125100,H3PO2 

819-124900,S2O3-- 

820-123641,MnF+ 

821-122400,H2PO2- 

822-121000,HMnO2- 

823-120700,RuCl5-- 

824-120400,MnO4-- 

825-120300,Pt(SO4) 

826-119800,HInO2 

827-116300,SO3-- 

828-115971,CdCl3- 

829-115609,Al+++ 

830-115316,BeCl+ 

831-112280,AgCl4--- 

832-111670,TiO2++ 

833-111500,VOH++ 

834-111430,Ag(CO3)- 

835-110720,HZnO2- 

836-108505,Mg++ 

837-108100,HSeO4- 

838-108000,LiOH 

839-107600,MnO4- 

840-106988,HgCl4-- 

841-106700,InO2- 

842-106700,VO++ 

843-106100,VO+ 

844-105500,SeO4-- 

845-105100,RbOH 

846-105000,CsOH 

847-104500,KOH 

848-104109,ZnF+ 

849-103900,PdCl4-- 

850-103579,CuCl4-- 

851-102600,MnO2-- 

852-102150,PbCl3- 

853-101850,H2SeO3 

854-101100,HFeO2 

855-100900,CsCl 

856-100500,CrOH++ 

857-99900,NaOH 

858-99800,VOH+ 

859-99250,LiCl 

860-98340,HSeO3- 

861-98300,ZnCl2 

862-97870,RbCl 

863-97400,HSbO2 

864-97300,HSnO2- 

865-97300,MnOH+ 

866-97016,InF++ 

867-96240,HAsO2 

868-95430,KCl 

869-95400,HFeO2- 

870-94610,CsBr 

871-93290,ZnO2-- 

872-93250,RhCl4-- 

873-92910,NaCl 

874-92800,CrO+ 

875-92250,CO2 

876-91210,PtCl4-- 

877-91157,FeF+ 

878-91100,GaOH++ 

879-91010,RbBr 

880-90550,Be++ 

881-90010,KBr 

882-89963,CuCl3-- 

883-89730,RuCl4- 

884-88400,SeO3-- 

885-88000,FeO2- 

886-87373,CdF+ 

887-86600,GaO+ 

888-86500,HCdO2- 

889-86290,MnCl+ 

890-85610,NaBr 

891-84851,CdCl2 

892-83900,RuCl4-- 

893-83650,AsO2- 

894-83600,Ti+++ 

895-83460,CsI 

896-83400,HCoO2- 

897-82710,AgCl3-- 

898-82400,SbO2- 

899-81980,HNiO2- 

900-81732,CoF+ 

901-81500,MnO 

902-81190,ZnOH+ 

903-81000,HPbO2- 

904-79768,NiF+ 

905-79645,FeF++ 

906-79300,HBiO2 

907-78900,RbI 

908-77740,KI 

909-77700,La++ 

910-77500,RhCl4- 

911-75860,PbF+ 

912-75338,CuCl3- 

913-75216,TlF 

914-75100,Ti++ 

915-74600,InOH++ 

916-74504,HgCl3- 

917-73480,FeCl2 

918-72900,NaI 

919-71980,SO2 

920-71662,HF 

921-71600,RuO4-- 

922-71200,PbCl2 

923-69933,Li+ 

924-69810,PdCl3- 

925-69710,Cs+ 

926-69400,InO+ 

927-67811,AuCl3-- 

928-67800,Rb+ 

929-67510,K+ 

930-67420,ZnO 

931-67340,F- 

932-67300,CdO2-- 

933-66850,ZnCl+ 

934-65850,FeOH+ 

935-65550,TlOH 

936-64200,NiO2-- 

937-63530,RhCl3- 

938-63200,CoO2-- 

939-62591,Na+ 

940-61700,BiO2- 

941-61500,CdOH+ 

942-60100,HCuO2- 

943-59226,InCl++ 

944-58600,SnOH+ 

945-58560,RuCl3 

946-58038,CuCl2- 

947-57900,V+++ 

948-57800,FeOH++ 

949-57760,PtCl3- 

950-57600,HTlO2 

951-56690,H2O 

952-56025,CoOH+ 

953-55100,Mn++ 

954-54380,RuCl3- 

955-53950,PbOH+ 

956-53739,CuF+ 

957-53600,SnO 

958-53100,FeO+ 

959-53030,FeCl+ 

960-52850,NiOH+ 

961-52627,CdCl+ 

962-52000,V++ 

963-51560,AgCl2- 

964-50720,FeO 

965-49459,AgF 

966-49300,Cr+++ 

967-47500,CdO 

968-46190,RhCl3 

969-46142,CuCl2 

970-45200,HHgO2- 

971-45157,CoCl+ 

972-44000,CoO 

973-42838,HgCl2 

974-41600,TlO2- 

975-41200,CuO2-- 

976-40920,NiCl+ 

977-39815,TlCl 

978-39400,Cr++ 

979-39350,PbO 

980-39340,NiO 

981-39050,PbCl+ 

982-38000,Ga+++ 

983-37518,FeCl++ 

984-36781,AuCl2- 

985-35332,AuCl4- 

986-35200,Zn++ 

987-35160,PdCl2 

988-33970,RhCl2 

989-32300,BiOH++ 

990-31700,HIO3 

991-31379,Cl- 

992-30600,IO3- 

993-30410,HCl 

994-30204,HgF+ 

995-30200,CuOH+ 

996-29300,BiO+ 

997-28682,CO 

998-26507,NO3- 

999-26440,RuCl2+ 

1000-25590,Br3- 

1001-25060,RuCl2 

1002-24870,Br- 

1003-24730,HNO3 

1004-23700,HIO 

1005-23400,In+++ 

1006-23280,OCN- 

1007-23000,CoOH++ 

1008-22608,CuCl 

1009-22290,PtCl2 

1010-21900,AgOH 

1011-21870,Fe++ 

1012-20800,CuO 

1013-20300,Mn+++ 

1014-20058,Pb(HS)2 

1015-19700,HBrO 

1016-19100,HClO 

1017-19100,ScOH++ 

1018-18990,NH4+ 

1019-18971,Pb(HS)3- 

1020-18560,Cd++ 

1021-18290,Rh(OH)+ 

1022-17450,AgCl 

1023-16250,CuCl+ 

1024-14780,RhCl2+ 

1025-14000,IO4- 

1026-13130,Pd(OH)+ 

1027-13000,Co++ 

1028-12700,HgOH+ 

1029-12410,I- 

1030-12300,I3- 

1031-12190,Ru(OH)2++ 

1032-12100,HNO2 

1033-11500,PdO 

1034-10900,Ni++ 

1035-10470,Ru(OH)+ 

1036-10450,RuO+ 

1037-9200,IO- 

1038-8900,HgO 

1039-8800,ClO- 

1040-8000,BrO- 

1041-7740,Tl+ 

1042-7738,AgNO3 

1043-7700,NO2- 

1044-7220,RhO 

1045-6673,H2S 

1046-6570,Sn++ 

1047-6383,NH3 

1048-5710,Pb++ 

1049-5500,AgO- 

1050-4500,TlOH++ 

1051-4120,Fe+++ 

1052-3380,RhCl+ 

1053-3200,TlO+ 

1054-3184,AuCl 

1055-2155,HgCl+ 

1056-2040,ClO4- 

1057-1900,ClO3- 

1058-1130,PtO 

1059-820,Rh(OH)++ 

10600,Ag(HS)2- 

10610,H+ 

1062230,RuO 

10631400,HClO2 

10641560,Pt(OH)+ 

10652429,Au(HS)2- 

10662500,PdCl+ 

10672860,HS- 

10683140,RhO+ 

10693215,Xe 

10703554,Kr 

10713890,Ar 

10724100,ClO2- 

10734347,N2 

10744450,BrO3- 

10754565,Ne 

10764658,He 

10775210,RuCl+ 

10787100,RuCl++ 

10798600,H2N2O2 

10809375,TlCl++ 

108110500,HSe- 

108211950,Cu+ 

108315675,Cu++ 

108415700,S5-- 

108516500,S4-- 

108617600,S3-- 

108718200,HN2O2- 

108818330,RhCl++ 

108918380,PtCl+ 

109018427,Ag+ 

109119000,S2-- 

109219500,SeCN- 

109319700,N2H5+ 

109421100,N2H6++ 

109522160,SCN- 

109622880,Bi+++ 

109727700,Rh++ 

109828200,BrO4- 

109928600,HCN 

110032000,Co+++ 

110133200,N2O2-- 

110235900,Ru++ 

110336710,Hg2++ 

110439360,Hg++ 

110541200,CN- 

110641440,Ru+++ 

110742200,Pd++ 

110851300,Tl+++ 

110952450,Rh+++ 

111061600,Pt++ 

111164300,Ag++ 

1112103600,Au+++"""