/
thermo.py
205 lines (164 loc) · 6.72 KB
/
thermo.py
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
""" Core definition of a Thermo Document """
from collections import defaultdict
from typing import Dict, List, Union
from datetime import datetime
from pydantic import BaseModel, Field
from pymatgen.analysis.phase_diagram import PhaseDiagram
from pymatgen.entries.computed_entries import ComputedEntry, ComputedStructureEntry
from emmet.core.material_property import PropertyDoc
from emmet.core.mpid import MPID
class DecompositionProduct(BaseModel):
"""
Entry metadata for a decomposition process
"""
material_id: MPID = Field(
None,
description="The Materials Project ID for the material this decomposition points to.",
)
formula: str = Field(
None,
description="The formula of the decomposed material this material decomposes to.",
)
amount: float = Field(
None,
description="The amount of the decomposed material by formula units this this material decomposes to.",
)
class ThermoDoc(PropertyDoc):
"""
A thermo entry document
"""
property_name = "thermo"
uncorrected_energy_per_atom: float = Field(
..., description="The total DFT energy of this material per atom in eV/atom."
)
energy_per_atom: float = Field(
...,
description="The total corrected DFT energy of this material per atom in eV/atom.",
)
energy_uncertainy_per_atom: float = Field(None, description="")
formation_energy_per_atom: float = Field(
None, description="The formation energy per atom in eV/atom."
)
energy_above_hull: float = Field(
..., description="The energy above the hull in eV/Atom."
)
is_stable: bool = Field(
False,
description="Flag for whether this material is on the hull and therefore stable.",
)
equilibrium_reaction_energy_per_atom: float = Field(
None,
description="The reaction energy of a stable entry from the neighboring equilibrium stable materials in eV."
" Also known as the inverse distance to hull.",
)
decomposes_to: List[DecompositionProduct] = Field(
None,
description="List of decomposition data for this material. Only valid for metastable or unstable material.",
)
decomposition_enthalpy: float = Field(
None,
description="Decomposition enthalpy as defined by `get_decomp_and_phase_separation_energy` in pymatgen."
)
decomposition_enthalpy_decomposes_to: List[DecompositionProduct] = Field(
None,
description="List of decomposition data associated with the decomposition_enthalpy quantity."
)
energy_type: str = Field(
...,
description="The type of calculation this energy evaluation comes from.",
)
entry_types: List[str] = Field(
description="List of available energy types computed for this material."
)
entries: Dict[str, Union[ComputedEntry, ComputedStructureEntry]] = Field(
...,
description="List of all entries that are valid for this material."
" The keys for this dictionary are names of various calculation types.",
)
@classmethod
def from_entries(
cls, entries: List[Union[ComputedEntry, ComputedStructureEntry]], **kwargs
):
entries_by_comp = defaultdict(list)
for e in entries:
entries_by_comp[e.composition.reduced_formula].append(e)
# Only use lowest entry per composition to speed up QHull in Phase Diagram
reduced_entries = [
sorted(comp_entries, key=lambda e: e.energy_per_atom)[0]
for comp_entries in entries_by_comp.values()
]
pd = PhaseDiagram(reduced_entries)
docs = []
for e in entries:
(decomp, ehull) = pd.get_decomp_and_e_above_hull(e)
d = {
"material_id": e.entry_id,
"uncorrected_energy_per_atom": e.uncorrected_energy
/ e.composition.num_atoms,
"energy_per_atom": e.energy / e.composition.num_atoms,
"formation_energy_per_atom": pd.get_form_energy_per_atom(e),
"energy_above_hull": ehull,
"is_stable": e in pd.stable_entries,
}
if "last_updated" in e.data:
d["last_updated"] = e.data["last_updated"]
# Store different info if stable vs decomposes
if d["is_stable"]:
d[
"equilibrium_reaction_energy_per_atom"
] = pd.get_equilibrium_reaction_energy(e)
else:
d["decomposes_to"] = [
{
"material_id": de.entry_id,
"formula": de.composition.formula,
"amount": amt,
}
for de, amt in decomp.items()
]
try:
decomp, energy = pd.get_decomp_and_phase_separation_energy(e)
d["decomposition_enthalpy"] = energy
d["decomposition_enthalpy_decomposes_to"] = [
{
"material_id": de.entry_id,
"formula": de.composition.formula,
"amount": amt,
}
for de, amt in decomp.items()
]
except ValueError:
# try/except so this quantity does not take down the builder if it fails:
# it includes an optimization step that can be fragile in some instances,
# most likely failure is ValueError, "invalid value encountered in true_divide"
d["warnings"] = ["Could not calculate decomposition enthalpy for this entry."]
d["energy_type"] = e.parameters.get("run_type", "Unknown")
d["entry_types"] = [e.parameters.get("run_type", "Unknown")]
d["entries"] = {e.parameters.get("run_type", ""): e}
for k in ["last_updated"]:
if k in e.parameters:
d[k] = e.parameters[k]
elif k in e.data:
d[k] = e.data[k]
docs.append(
ThermoDoc.from_structure(meta_structure=e.structure, **d, **kwargs)
)
return docs, pd
class PhaseDiagramDoc(BaseModel):
"""
A phase diagram document
"""
property_name = "phase_diagram"
chemsys: str = Field(
...,
title="Chemical System",
description="Dash-delimited string of elements in the material",
)
phase_diagram: PhaseDiagram = Field(
...,
description="Phase diagram for the chemical system.",
)
last_updated: datetime = Field(
description="Timestamp for the most recent calculation update for this property",
default_factory=datetime.utcnow,
)