/
volumetric_data.py
292 lines (259 loc) · 11.8 KB
/
volumetric_data.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
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
# coding: utf-8
# Copyright (c) Max-Planck-Institut für Eisenforschung GmbH - Computational Materials Design (CM) Department
# Distributed under the terms of "New BSD License", see the LICENSE file.
import math
import numpy as np
import os
from pyiron.base.settings.generic import Settings
from pyiron.vasp.structure import atoms_from_string, get_species_list_from_potcar
from pyiron.atomistics.volumetric.generic import VolumetricData
__author__ = "Sudarsan Surendralal"
__copyright__ = "Copyright 2019, Max-Planck-Institut für Eisenforschung GmbH - " \
"Computational Materials Design (CM) Department"
__version__ = "1.0"
__maintainer__ = "Sudarsan Surendralal"
__email__ = "surendralal@mpie.de"
__status__ = "production"
__date__ = "Sep 1, 2017"
class VaspVolumetricData(VolumetricData):
"""
General class for parsing and manipulating volumetric static within VASP. The basic idea of the Base class is
adapted from the pymatgen vasp VolumtricData class
http://pymatgen.org/_modules/pymatgen/io/vasp/outputs.html#VolumetricData
"""
def __init__(self):
super(VaspVolumetricData, self).__init__()
self.atoms = None
self._diff_data = None
self._total_data = None
def from_file(self, filename, normalize=True):
"""
Parsing the contents of from a file
Args:
filename (str): Path of file to parse
normalize (boolean): Flag to normalize by the volume of the cell
"""
try:
self.atoms, vol_data_list = self._read_vol_data(filename=filename, normalize=normalize)
except (ValueError, IndexError, TypeError):
try:
self.atoms, vol_data_list = self._read_vol_data_old(filename=filename, normalize=normalize)
except (ValueError, IndexError, TypeError):
raise ValueError("Unable to parse file: {}".format(filename))
if self.atoms is not None:
self._total_data = vol_data_list[0]
if len(vol_data_list) > 1:
self._diff_data = vol_data_list[1]
@staticmethod
def _read_vol_data_old(filename, normalize=True):
"""
Convenience method to parse a generic volumetric static file in the vasp like format.
Used by subclasses for parsing the file. This routine is adapted from the pymatgen vasp VolumetricData
class with very minor modifications. The new parser is faster
http://pymatgen.org/_modules/pymatgen/io/vasp/outputs.html#VolumetricData.
Args:
filename (str): Path of file to parse
normalize (boolean): Flag to normalize by the volume of the cell
"""
if os.stat(filename).st_size == 0:
s = Settings()
s.logger.warning("File:" + filename + "seems to be corrupted/empty")
return None, None
poscar_read = False
poscar_string = list()
dataset = list()
all_dataset = list()
dim = None
dimline = None
read_dataset = False
ngrid_pts = 0
data_count = 0
atoms = None
volume = None
with open(filename) as f:
for line in f:
line = line.strip()
if read_dataset:
toks = line.split()
for tok in toks:
if data_count < ngrid_pts:
# This complicated procedure is necessary because
# vasp outputs x as the fastest index, followed by y
# then z.
x = data_count % dim[0]
y = int(math.floor(data_count / dim[0])) % dim[1]
z = int(math.floor(data_count / dim[0] / dim[1]))
dataset[x, y, z] = float(tok)
data_count += 1
if data_count >= ngrid_pts:
read_dataset = False
data_count = 0
all_dataset.append(dataset)
elif not poscar_read:
if line != "" or len(poscar_string) == 0:
poscar_string.append(line)
elif line == "":
try:
atoms = atoms_from_string(poscar_string)
except ValueError:
pot_str = filename.split("/")
pot_str[-1] = "POTCAR"
potcar_file = "/".join(pot_str)
species = get_species_list_from_potcar(potcar_file)
atoms = atoms_from_string(poscar_string, species_list=species)
volume = atoms.get_volume()
poscar_read = True
elif not dim:
dim = [int(i) for i in line.split()]
ngrid_pts = dim[0] * dim[1] * dim[2]
dimline = line
read_dataset = True
dataset = np.zeros(dim)
elif line == dimline:
read_dataset = True
dataset = np.zeros(dim)
if not normalize:
volume = 1.0
if len(all_dataset) == 0:
s = Settings()
s.logger.warning("File:" + filename + "seems to be corrupted/empty")
return None, None
if len(all_dataset) == 2:
data = {"total": all_dataset[0] / volume, "diff": all_dataset[1] / volume}
return atoms, [data["total"], data["diff"]]
else:
data = {"total": all_dataset[0] / volume}
return atoms, [data["total"]]
def _read_vol_data(self, filename, normalize=True):
"""
Parses the VASP volumetric type files (CHGCAR, LOCPOT, PARCHG etc). Rather than looping over individual values,
this function utilizes numpy indexing resulting in a parsing efficiency of at least 10%.
Args:
filename (str): File to be parsed
normalize (bool): Normalize the data with respect to the volume (Recommended for CHGCAR files)
Returns:
pyiron.atomistics.structure.atoms.Atoms: The structure of the volumetric snapshot
list: A list of the volumetric data (length >1 for CHGCAR files with spin)
"""
if os.stat(filename).st_size == 0:
s = Settings()
s.logger.warning("File:" + filename + "seems to be corrupted/empty")
return None, None
with open(filename, "r") as f:
struct_lines = list()
get_grid = False
n_x = 0
n_y = 0
n_z = 0
n_grid = 0
n_grid_str = None
total_data_list = list()
atoms = None
for line in f:
strip_line = line.strip()
if not get_grid:
if strip_line == "":
get_grid = True
struct_lines.append(strip_line)
elif n_grid_str is None:
n_x, n_y, n_z = [int(val) for val in strip_line.split()]
n_grid = n_x * n_y * n_z
n_grid_str = " ".join([str(val) for val in [n_x, n_y, n_z]])
load_txt = np.genfromtxt(f, max_rows=int(n_grid / 5))
load_txt = np.hstack(load_txt)
if n_grid % 5 != 0:
add_line = np.genfromtxt(f, max_rows=1)
load_txt = np.append(load_txt, np.hstack(add_line))
total_data = self._fastest_index_reshape(load_txt, [n_x, n_y, n_z])
try:
atoms = atoms_from_string(struct_lines)
except ValueError:
pot_str = filename.split("/")
pot_str[-1] = "POTCAR"
potcar_file = "/".join(pot_str)
species = get_species_list_from_potcar(potcar_file)
atoms = atoms_from_string(struct_lines, species_list=species)
if normalize:
total_data /= atoms.get_volume()
total_data_list.append(total_data)
elif atoms is not None:
grid_str = n_grid_str.replace(" ", "")
if grid_str == strip_line.replace(" ", ""):
load_txt = np.genfromtxt(f, max_rows=int(n_grid / 5))
load_txt = np.hstack(load_txt)
if n_grid % 5 != 0:
add_line = np.genfromtxt(f, max_rows=1)
load_txt = np.hstack(np.append(load_txt, np.hstack(add_line)))
total_data = self._fastest_index_reshape(load_txt, [n_x, n_y, n_z])
if normalize:
total_data /= atoms.get_volume()
total_data_list.append(total_data)
if len(total_data_list) == 0:
s = Settings()
s.logger.warning("File:" + filename + "seems to be corrupted/empty")
return None, None
return atoms, total_data_list
@staticmethod
def _fastest_index_reshape(raw_data, grid):
"""
Helper function to parse volumetric data with x-axis as the fastest index into a 3D numpy array
Args:
raw_data (numpy.ndarray): Raw unprocessed volumetric data which is flattened
grid (list/turple/numpy.ndarray): Sequence of the integer grid points [Nx, Ny, Nz]
Returns:
numpy.ndarray: A Nx $\times$ Ny $\times$ Nz numpy array
"""
n_x, n_y, n_z = grid
total_data = np.zeros((n_x, n_y, n_z))
all_data = raw_data[0:np.prod(grid)]
all_indices = np.arange(len(all_data), dtype=int)
x_indices = all_indices % n_x
y_indices = all_indices / n_x % n_y
y_indices = np.array(y_indices, dtype=int)
z_indices = all_indices / (n_x * n_y)
z_indices = np.array(z_indices, dtype=int)
total_data[x_indices, y_indices, z_indices] = all_data
return total_data
@property
def total_data(self):
"""
numpy.ndarray: Total volumtric data (3D)
"""
return self._total_data
@total_data.setter
def total_data(self, val):
self._total_data = val
@property
def diff_data(self):
"""
numpy.ndarray: Volumtric difference data (3D)
"""
return self._diff_data
@diff_data.setter
def diff_data(self, val):
self._diff_data = val
def to_hdf(self, hdf5, group_name="volumetric_data"):
"""
Writes the data as a group to a HDF5 file
Args:
hdf5 (pyiron.base.generic.hdfio.ProjectHDFio): The HDF file/path to write the data to
group_name (str): The name of the group under which the data must be stored as
"""
with hdf5.open(group_name) as hdf_vd:
hdf_vd["TYPE"] = str(type(self))
hdf_vd["total"] = self.total_data
if self.diff_data is not None:
hdf_vd["diff"] = self.diff_data
def from_hdf(self, hdf5, group_name="volumetric_data"):
"""
Recreating the VolumetricData instance by reading data from the HDF5 files
Args:
hdf5 (pyiron.base.generic.hdfio.ProjectHDFio): The HDF file/path to write the data to
group_name (str): The name of the group under which the data must be stored as
Returns:
pyiron.atomistics.volumetric.generic.VolumetricData: The VolumetricData instance
"""
with hdf5.open(group_name) as hdf_vd:
self._total_data = hdf_vd["total"]
if "diff" in hdf_vd.list_nodes():
self._diff_data = hdf_vd["diff"]