wreader.py
Go to the documentation of this file.
1 # Author: M. Reichert
2 # Date: 04.07.2024
3 import numpy as np
4 from tqdm import tqdm
5 from nucleus_multiple_class import nucleus_multiple
6 import h5py
7 import os
8 
9 
10 
11 class wreader(object):
12  """
13  Minimalistic class to lazily read WinNet data.
14  """
15 
16  def __init__(self, path):
17  """
18  Initialize the class
19  - path: Path to the WinNet data
20  """
21  # The path to the WinNet run
22  self.path = path
23 
24  # The path to the hdf5 file
25  self.filename = os.path.join(path, "WinNet_data.h5")
26 
27  # Check whether hdf5 mode was activated or not
28  if not os.path.exists(self.filename):
29  self.__mode = 'ascii'
30  self.__snapshot_path = os.path.join(path, "snaps")
31  else:
32  self.__mode = 'hdf5'
33 
34 
35  @property
36  def A(self):
37  """
38  Mass number
39  """
40  if not hasattr(self,"_wreader__A"):
41  self.__read_A_Z_N()
42  return self.__A
43 
44  def __read_A_Z_N(self):
45  """
46  Read the mass number A
47  """
48  # Check if Ascii (mode = 2), hdf5 (mode = 1), or not present (mode = 0)
49  mode = self.check_existence('snapshot')
50 
51  if mode == 2:
52  self.__N, self.__Z = np.loadtxt(os.path.join(self.__snapshot_path, "snapsh_0001.dat"),unpack=True, usecols=(0,1),
53  skiprows=3,dtype=int)
54  self.__A = self.__N + self.__Z
55  elif mode == 1:
56  with h5py.File(self.filename, 'r') as hf:
57  self.__A = hf['snapshots/A'][:]
58  self.__Z = hf['snapshots/Z'][:]
59  self.__N = hf['snapshots/N'][:]
60  elif mode == 0:
61  error_msg = 'Failed to read A, Z, and N in "'+str(self.path)+'". '+\
62  'Snapshots not present as Ascii nor as Hdf5!'
63  raise ValueError(error_msg)
64 
65  @property
66  def Z(self):
67  """
68  Atomic number
69  """
70  if not hasattr(self,"_wreader__Z"):
71  self.__read_A_Z_N()
72  return self.__Z
73 
74  @property
75  def N(self):
76  """
77  Neutron number
78  """
79  if not hasattr(self,"_wreader__N"):
80  self.__read_A_Z_N()
81  return self.__N
82 
83 
84  def __read_snapshots(self):
85  """
86  Read the snapshots
87  """
88  # Check if Ascii (mode = 2), hdf5 (mode = 1), or not present (mode = 0)
89  mode = self.check_existence('snapshot')
90 
91  if mode == 2:
92  # Get list of files (no directories)
93  snapshot_files = [f for f in os.listdir(self.__snapshot_path) if os.path.isfile(os.path.join(self.__snapshot_path, f))]
94  # Now only take files that start with a 'snapsh_'
95  snapshot_files = [f for f in snapshot_files if f.startswith('snapsh_')]
96  self.__snapshots_time = np.zeros(len(snapshot_files))
97  self.__snapshots_Y = np.zeros((len(snapshot_files), len(self.A)))
98  self.__snapshots_X = np.zeros((len(snapshot_files), len(self.A)))
99  for i, f in enumerate(tqdm(snapshot_files, desc='Reading snapshots')):
100  fname = 'snapsh_' + str(i+1).zfill(4) + '.dat'
101  with open(os.path.join(self.__snapshot_path, fname), 'r') as file:
102  lines = file.readlines()
103  self.__snapshots_time[i] = float(lines[1].split()[0])
104  for j, line in enumerate(lines[3:]):
105  self.__snapshots_Y[i, j] = float(line.split()[2])
106  self.__snapshots_X[i, j] = float(line.split()[3])
107  elif mode == 1:
108  with h5py.File(self.filename, 'r') as hf:
109  self.__snapshots_time = hf['snapshots/time'][:]
110  self.__snapshots_Y = hf['snapshots/Y'][:]
111  self.__snapshots_X = self.__snapshots_Y*self.A
112  elif mode == 0:
113  error_msg = 'Failed to read timescales in "'+str(self.path)+'". '+\
114  'Not present as Ascii nor as Hdf5!'
115  raise ValueError(error_msg)
116 
117 
118 
119  def check_existence(self, entry):
120  """
121  Check whether an entry is in the hdf5 format (return 1),
122  or in the ascii format (return 2), or does not exist (return 0)
123  """
124  if (entry == 'mainout'):
125  value = self.__check_files('mainout.dat','mainout')
126  elif (entry == 'timescales'):
127  value = self.__check_files('timescales.dat','timescales')
128  elif (entry == 'energy'):
129  value = self.__check_files('generated_energy.dat','energy')
130  elif (entry == 'tracked_nuclei'):
131  value = self.__check_files('tracked_nuclei.dat','tracked_nuclei')
132  elif (entry == 'snapshot'):
133  value = self.__check_files('snaps/snapsh_0001.dat','snapshots')
134  elif (entry == 'flows'):
135  value = self.__check_files('flow/flow_0001.dat','flows')
136  elif (entry == 'finab'):
137  value = self.__check_files('finab.dat','finab')
138  elif (entry == 'finabsum'):
139  value = self.__check_files('finabsum.dat','finab')
140  elif (entry == 'finabelem'):
141  value = self.__check_files('finabelem.dat','finab')
142  else:
143  error_msg = 'Checked for unknown entry "'+str(entry)+'". '
144  raise ValueError(error_msg)
145 
146  return value
147 
148 
149  def __check_files(self, ascii_file_path, hdf5_key):
150  '''
151  Check if something exists in hdf5 or ascii
152  '''
153  value = 0
154  if not os.path.exists(self.filename):
155  if not os.path.exists(os.path.join(self.path, ascii_file_path)):
156  value = 0
157  else:
158  value = 2
159  else:
160  with h5py.File(self.filename, 'r') as hf:
161  if hdf5_key in hf.keys():
162  value = 1
163  elif os.path.exists(os.path.join(self.path, ascii_file_path)):
164  value = 2
165  else:
166  value = 0
167 
168  return value
169 
170  @property
171  def nr_of_snaps(self):
172  """
173  Number of snapshots
174  """
175  if not hasattr(self,"_wreader__nr_of_snaps"):
176  self.__read_nr_of_snaps()
177  return self.__nr_of_snaps
178 
180  """
181  Read the number of snapshots
182  """
183  # Check if Ascii (mode = 2), hdf5 (mode = 1), or not present (mode = 0)
184  mode = self.check_existence('snapshot')
185 
186  if mode == 2:
187  snapshot_files = [f for f in os.listdir(self.__snapshot_path) if os.path.isfile(os.path.join(self.__snapshot_path, f))]
188  snapshot_files = [f for f in snapshot_files if f.startswith('snapsh_')]
189  self.__nr_of_snaps = len(snapshot_files)
190  elif mode == 1:
191  with h5py.File(self.filename, 'r') as hf:
192  self.__nr_of_snaps = len(hf['snapshots/time'][:])
193  elif mode == 0:
194  error_msg = 'Failed to read nr of snapshots in "'+str(self.path)+'". '+\
195  'Not present as Ascii nor as Hdf5!'
196  raise ValueError(error_msg)
197 
198  @property
199  def tracked_nuclei(self):
200  """
201  Get the tracked nuclei
202  """
203  if not hasattr(self,"_wreader__tracked_nuclei"):
204  self.__read_tracked_nuclei()
205  return self.__tracked_nuclei
206 
208  """
209  Read the tracked nuclei
210  """
211  # Check if Ascii (mode = 2), hdf5 (mode = 1), or not present (mode = 0)
212  mode = self.check_existence('tracked_nuclei')
214 
215  if mode == 2:
216  # self.__tracked_nuclei = np.loadtxt(os.path.join(self.path, "tracked_nuclei.dat"), unpack=True)
217  path = os.path.join(self.path, "tracked_nuclei.dat")
218  with open(path, 'r') as f:
219  first_line = f.readline()
220 
221  first_line = first_line.replace('Y(','')
222  first_line = first_line.replace(')','')
223  first_line = first_line.replace('#','')
224  track_nuclei_names = first_line.split()[1:] # Since the first entry is "time"
225  track_nuclei_data = np.loadtxt(path, skiprows=1)
226  nm = nucleus_multiple(names=track_nuclei_names)
227  Z = nm.Z
228  N = nm.N
229  A = nm.A
230  Xnuclei = track_nuclei_data[:,1:]*A
231  self.__tracked_nuclei['time'] = track_nuclei_data[:,0]
232  self.__tracked_nuclei['Z'] = Z
233  self.__tracked_nuclei['N'] = N
234  self.__tracked_nuclei['A'] = A
235  self.__tracked_nuclei['names'] = track_nuclei_names
236  elnames = nm.elnames
237  self.__tracked_nuclei['latex_names'] = [r"$^{"+str(A[i])+r"}$"+elnames[i].title() for i in range(len(A))]
238  for i, name in enumerate(track_nuclei_names):
239  self.__tracked_nuclei[name] = Xnuclei[:,i]
240  elif mode == 1:
241  with h5py.File(self.filename, 'r') as hf:
242  self.__tracked_nuclei['Z'] = hf['tracked_nuclei/Z'][:]
243  self.__tracked_nuclei['N'] = hf['tracked_nuclei/N'][:]
244  self.__tracked_nuclei['A'] = hf['tracked_nuclei/A'][:]
245  self.__tracked_nuclei['time'] = hf['tracked_nuclei/time'][:]
246  Xnuclei = hf['tracked_nuclei/Y'][:,:] * self.__tracked_nuclei['A']
247  nm = nucleus_multiple(Z=self.__tracked_nuclei['Z'], N=self.__tracked_nuclei['N'])
248  self.__tracked_nuclei['names'] = nm.names
249  for i, name in enumerate(nm.names):
250  self.__tracked_nuclei[name] = Xnuclei[:,i]
251 
252  elnames = nm.elnames
253  A = self.__tracked_nuclei['A']
254  self.__tracked_nuclei['latex_names'] = [r"$^{"+str(A[i])+r"}$"+elnames[i].title() for i in range(len(A))]
255  elif mode == 0:
256  error_msg = 'Failed to read tracked nuclei in "'+str(self.path)+'". '+\
257  'Not present as Ascii nor as Hdf5!'
258  raise ValueError(error_msg)
259 
260  @property
261  def Y(self):
262  """
263  Get abundance at snapshot idx
264  """
265  if not hasattr(self,"_wreader__snapshots_time"):
266  self.__read_snapshots()
267  return self.__snapshots_Y
268 
269  @property
270  def X(self):
271  """
272  Get mass fraction at snapshot idx
273  """
274  if not hasattr(self,"_wreader__snapshots_time"):
275  self.__read_snapshots()
276  return self.__snapshots_X
277 
278  @property
279  def snapshot_time(self):
280  """
281  Get time at snapshot idx
282  """
283  if not hasattr(self,"_wreader__snapshots_time"):
284  self.__read_snapshots()
285  return self.__snapshots_time
286 
287  @property
288  def tau(self):
289  """
290  Get the timescale of "tau", e.g., "tau_ag"
291  """
292  if not hasattr(self,"_wreader__timescales"):
293  self.__read_timescales()
294  return self.__timescales
295 
296 
297  def __read_timescales(self):
298  """
299  Read the timescales
300  """
301  # Check if Ascii (mode = 2), hdf5 (mode = 1), or not present (mode = 0)
302  mode = self.check_existence('timescales')
303 
304  if (mode == 2):
305  self.__timescales = {}
306  with open(os.path.join(self.path, "timescales.dat"), 'r') as f:
307  lines = f.readlines()
308  header = lines[0]
309  header = header.replace('#', '').replace('[s]', '').replace('[GK]', '')
310  key = header.split()
311  data = np.loadtxt(os.path.join(self.path, "timescales.dat"), unpack=True)
312  for i, k in enumerate(key):
313  self.__timescales[k] = data[i]
314 
315  elif (mode == 1):
316  with h5py.File(self.filename, 'r') as hf:
317  self.__timescales = {}
318  for key in hf['timescales'].keys():
319  self.__timescales[key] = hf['timescales'][key][:]
320  elif (mode == 0):
321  error_msg = 'Failed to read timescales in "'+str(self.path)+'". '+\
322  'Not present as Ascii nor as Hdf5!'
323  raise ValueError(error_msg)
324 
325 
326  @property
327  def mainout(self):
328  """
329  Get an entry from the mainout
330  """
331  if not hasattr(self,"_wreader__mainout"):
332  self.__read_mainout()
333  return self.__mainout
334 
335  def __read_mainout(self):
336  """
337  Read the mainout
338  """
339  # Check if Ascii (mode = 2), hdf5 (mode = 1), or not present (mode = 0)
340  mode = self.check_existence('mainout')
341 
342  if mode == 2:
343  columns = ['iteration', 'time', 'temp', 'dens', 'ye',
344  'rad', 'yn', 'yp', 'ya', 'ylight', 'yheavy',
345  'zbar', 'abar', 'entr']
346  data = np.loadtxt(os.path.join(self.path, "mainout.dat"), unpack=True)
347  self.__mainout = {}
348  for i, col in enumerate(columns):
349  self.__mainout[col] = data[i]
350  elif mode == 1:
351  with h5py.File(self.filename, 'r') as hf:
352  self.__mainout = {}
353  for key in hf['mainout'].keys():
354  self.__mainout[key] = hf['mainout'][key][:]
355  elif mode == 0:
356  error_msg = 'Failed to read mainout in "'+str(self.path)+'". '+\
357  'Not present as Ascii nor as Hdf5!'
358  raise ValueError(error_msg)
359 
360 
361  @property
362  def energy(self):
363  """
364  Get the energy
365  """
366  if not hasattr(self,"_wreader__energy"):
367  self.__read_energy()
368  return self.__energy
369 
370  def __read_energy(self):
371  """
372  Read the energy
373  """
374  # Check if Ascii (mode = 2), hdf5 (mode = 1), or not present (mode = 0)
375  mode = self.check_existence('energy')
376 
377  if mode == 2:
378  columns = ['time', 'engen_tot', 'S_src', 'engen_beta', 'engen_ng_gn', 'engen_pg_gp',
379  'engen_ag_ga', 'engen_np_pn', 'engen_na_an', 'engen_ap_pa', 'engen_fiss']
380  data = np.loadtxt(os.path.join(self.path, "generated_energy.dat"), unpack=True)
381  self.__energy = {}
382  for i, col in enumerate(columns):
383  self.__energy[col] = data[i]
384  elif mode == 1:
385  with h5py.File(self.filename, 'r') as hf:
386  self.__energy = {}
387  for key in hf['energy'].keys():
388  self.__energy[key] = hf['energy'][key][:]
389  elif mode == 0:
390  error_msg = 'Failed to read energy in "'+str(self.path)+'". '+\
391  'Not present as Ascii nor as Hdf5!'
392  raise ValueError(error_msg)
393 
394  @property
395  def finab(self):
396  """
397  Get the final abundances from the finab.dat file
398  """
399  if not hasattr(self,"_wreader__finab"):
400  self.__read_finab()
401  return self.__finab
402 
403  def __read_finab(self):
404  """
405  Reader of the finab
406  """
407  # Check if Ascii (mode = 2), hdf5 (mode = 1), or not present (mode = 0)
408  mode = self.check_existence('finab')
409  if mode == 2:
410  A,Z,N,Y,X = np.loadtxt(os.path.join(self.path, "finab.dat"), unpack=True)
411  elif mode == 1:
412  with h5py.File(self.filename, 'r') as hf:
413  A = hf['finab/finab/A'][:]
414  Z = hf['finab/finab/Z'][:]
415  N = A-Z
416  Y = hf['finab/finab/Y'][:]
417  X = hf['finab/finab/X'][:]
418  elif mode == 0:
419  error_msg = 'Failed to read finab in "'+str(self.path)+'". '+\
420  'Not present as Ascii nor as Hdf5!'
421  raise ValueError(error_msg)
422  self.__finab = {"A": A, "Z": Z, "N": N, "Y": Y, "X": X}
423 
424 
425  @property
426  def finabsum(self):
427  """
428  Get the final abundances from the finabsum.dat file
429  """
430  if not hasattr(self,"_wreader__finabsum"):
431  self.__read_finabsum()
432  return self.__finabsum
433 
434  def __read_finabsum(self):
435  """
436  Reader of the finabsum
437  """
438  # Check if Ascii (mode = 2), hdf5 (mode = 1), or not present (mode = 0)
439  mode = self.check_existence('finabsum')
440  if mode == 2:
441  A,Y,X = np.loadtxt(os.path.join(self.path, "finabsum.dat"), unpack=True)
442  elif mode == 1:
443  with h5py.File(self.filename, 'r') as hf:
444  A = hf['finab/finabsum/A'][:]
445  Y = hf['finab/finabsum/Y'][:]
446  X = hf['finab/finabsum/X'][:]
447  elif mode == 0:
448  error_msg = 'Failed to read finabsum in "'+str(self.path)+'". '+\
449  'Not present as Ascii nor as Hdf5!'
450  raise ValueError(error_msg)
451  self.__finabsum = {"A": A, "Y": Y, "X": X}
452 
453 
454  @property
455  def finabelem(self):
456  """
457  Get the final abundances from the finabelem.dat file
458  """
459  if not hasattr(self,"_wreader__finabelem"):
460  self.__read_finabelem()
461  return self.__finabelem
462 
463  def __read_finabelem(self):
464  """
465  Reader of the finabsum
466  """
467  # Check if Ascii (mode = 2), hdf5 (mode = 1), or not present (mode = 0)
468  mode = self.check_existence('finabelem')
469  if mode == 2:
470  Z, Y = np.loadtxt(os.path.join(self.path, "finabelem.dat"), unpack=True)
471  elif mode == 1:
472  with h5py.File(self.filename, 'r') as hf:
473  Z = hf['finab/finabelem/Z'][:]
474  Y = hf['finab/finabelem/Y'][:]
475  elif mode == 0:
476  error_msg = 'Failed to read finabelem in "'+str(self.path)+'". '+\
477  'Not present as Ascii nor as Hdf5!'
478  raise ValueError(error_msg)
479  self.__finabelem = {"Z": Z, "Y": Y}
480 
481  def flow_entry(self, iteration, flow_group='flows'):
482  """
483  Get the flow entry
484  """
485  mode = self.check_existence('flows')
486 
487  if mode == 2:
488  fname = os.path.join(self.path, "flow/flow_" + str(iteration).zfill(4) + ".dat")
489  with open(fname, 'r') as f:
490  lines = f.readlines()
491  header = lines[1]
492  time, temp, dens = np.array(header.split()).astype(float)
493 
494  # Read the whole flow
495  nin,zin,yin,nout,zout,yout,flows = np.loadtxt(fname, unpack=True, skiprows=3)
496  nin = nin.astype(int)
497  zin = zin.astype(int)
498  nout = nout.astype(int)
499  zout = zout.astype(int)
500  flow = {}
501  flow['time'] = time
502  flow['temp'] = temp
503  flow['dens'] = dens
504  flow['n_in'] = nin
505  flow['p_in'] = zin
506  flow['y_in'] = yin
507  flow['n_out'] = nout
508  flow['p_out'] = zout
509  flow['y_out'] = yout
510  flow['flow'] = flows
511  elif mode == 1:
512  with h5py.File(self.filename, 'r') as hf:
513  flow = {}
514  for key in hf[flow_group][str(iteration)].keys():
515  flow[key] = hf[flow_group][str(iteration)][key][:]
516  elif mode == 0:
517  error_msg = 'Failed to read flow in "'+str(self.path)+'". '+\
518  'Not present as Ascii nor as Hdf5!'
519  raise ValueError(error_msg)
520 
521  return flow
522 
523 if __name__ == "__main__":
524  w = wreader('/home/mreichert/data/Networks/comparison_winNet/WinNet-dev/runs/Example_MRSN_r_process_winteler')
525  w.tracked_nuclei
526  pass
src_files.wreader.wreader.__N
__N
Definition: wreader.py:59
src_files.wreader.wreader.__read_energy
def __read_energy(self)
Definition: wreader.py:370
src_files.wreader.wreader.__read_tracked_nuclei
def __read_tracked_nuclei(self)
Definition: wreader.py:207
src_files.wreader.wreader.finabelem
def finabelem(self)
Definition: wreader.py:455
src_files.wreader.wreader.__read_nr_of_snaps
def __read_nr_of_snaps(self)
Definition: wreader.py:179
src_files.wreader.wreader.path
path
Definition: wreader.py:22
src_files.wreader.wreader.flow_entry
def flow_entry(self, iteration, flow_group='flows')
Definition: wreader.py:481
src_files.wreader.wreader.tracked_nuclei
def tracked_nuclei(self)
Definition: wreader.py:199
src_files.wreader.wreader.__read_finabelem
def __read_finabelem(self)
Definition: wreader.py:463
src_files.wreader.wreader.A
def A(self)
Definition: wreader.py:36
src_files.wreader.wreader.__nr_of_snaps
__nr_of_snaps
Definition: wreader.py:189
src_files.wreader.wreader.__energy
__energy
Definition: wreader.py:381
src_files.wreader.wreader.__read_mainout
def __read_mainout(self)
Definition: wreader.py:335
src_files.wreader.wreader.tau
def tau(self)
Definition: wreader.py:288
src_files.wreader.wreader.__snapshots_time
__snapshots_time
Definition: wreader.py:96
src_files.wreader.wreader.__mode
__mode
Definition: wreader.py:29
src_files.wreader.wreader.finabsum
def finabsum(self)
Definition: wreader.py:426
src_files.wreader.wreader.__tracked_nuclei
__tracked_nuclei
Definition: wreader.py:213
src_files.wreader.wreader.filename
filename
Definition: wreader.py:25
src_files.wreader.wreader.Y
def Y(self)
Definition: wreader.py:261
src_files.wreader.wreader.__snapshots_Y
__snapshots_Y
Definition: wreader.py:97
src_files.wreader.wreader.__read_snapshots
def __read_snapshots(self)
Definition: wreader.py:84
src_files.wreader.wreader.__read_finab
def __read_finab(self)
Definition: wreader.py:403
src_files.wreader.wreader.__read_A_Z_N
def __read_A_Z_N(self)
Definition: wreader.py:44
src_files.wreader.wreader.__read_timescales
def __read_timescales(self)
Definition: wreader.py:297
src_files.wreader.wreader.check_existence
def check_existence(self, entry)
Definition: wreader.py:119
src_files.wreader.wreader.__A
__A
Definition: wreader.py:54
src_files.wreader.wreader.finab
def finab(self)
Definition: wreader.py:395
src_files.wreader.wreader.__mainout
__mainout
Definition: wreader.py:347
src_files.wreader.wreader.__finab
__finab
Definition: wreader.py:422
src_files.wreader.wreader.__init__
def __init__(self, path)
Definition: wreader.py:16
src_files.wreader.wreader.__snapshot_path
__snapshot_path
Definition: wreader.py:30
src_files.wreader.wreader.__timescales
__timescales
Definition: wreader.py:305
src_files.wreader.wreader.mainout
def mainout(self)
Definition: wreader.py:327
src_files.wreader.wreader.N
def N(self)
Definition: wreader.py:75
src_files.nucleus_multiple_class.nucleus_multiple
Definition: nucleus_multiple_class.py:10
src_files.wreader.wreader.__read_finabsum
def __read_finabsum(self)
Definition: wreader.py:434
src_files.wreader.wreader.energy
def energy(self)
Definition: wreader.py:362
src_files.wreader.wreader.__finabelem
__finabelem
Definition: wreader.py:479
src_files.wreader.wreader.__snapshots_X
__snapshots_X
Definition: wreader.py:98
src_files.wreader.wreader.__check_files
def __check_files(self, ascii_file_path, hdf5_key)
Definition: wreader.py:149
src_files.wreader.wreader.__finabsum
__finabsum
Definition: wreader.py:451
src_files.wreader.wreader
Definition: wreader.py:11
src_files.wreader.wreader.nr_of_snaps
def nr_of_snaps(self)
Definition: wreader.py:171
src_files.wreader.wreader.X
def X(self)
Definition: wreader.py:270
src_files.wreader.wreader.Z
def Z(self)
Definition: wreader.py:66
src_files.wreader.wreader.__Z
__Z
Definition: wreader.py:52
src_files.wreader.wreader.snapshot_time
def snapshot_time(self)
Definition: wreader.py:279