pyerrors.input.hadrons
View Source
import os from collections import Counter import h5py import numpy as np from ..obs import Obs, CObs from ..correlators import Corr def _get_files(path, filestem, idl): ls = os.listdir(path) # Clean up file list files = list(filter(lambda x: x.startswith(filestem), ls)) if not files: raise Exception('No files starting with', filestem, 'in folder', path) def get_cnfg_number(n): return int(n[len(filestem) + 1:-3]) # Sort according to configuration number files.sort(key=get_cnfg_number) cnfg_numbers = [] filtered_files = [] for line in files: no = get_cnfg_number(line) if idl: if no in list(idl): filtered_files.append(line) cnfg_numbers.append(no) else: filtered_files.append(line) cnfg_numbers.append(no) if idl: if Counter(list(idl)) != Counter(cnfg_numbers): raise Exception("Not all configurations specified in idl found, configurations " + str(list(Counter(list(idl)) - Counter(cnfg_numbers))) + " are missing.") # Check that configurations are evenly spaced dc = np.unique(np.diff(cnfg_numbers)) if np.any(dc < 0): raise Exception("Unsorted files") if len(dc) == 1: idx = range(cnfg_numbers[0], cnfg_numbers[-1] + dc[0], dc[0]) else: raise Exception('Configurations are not evenly spaced.') return filtered_files, idx def read_meson_hd5(path, filestem, ens_id, meson='meson_0', tree='meson', idl=None): """Read hadrons meson hdf5 file and extract the meson labeled 'meson' Parameters ----------------- path : str path to the files to read filestem : str namestem of the files to read ens_id : str name of the ensemble, required for internal bookkeeping meson : str label of the meson to be extracted, standard value meson_0 which corresponds to the pseudoscalar pseudoscalar two-point function. idl : range If specified only configurations in the given range are read in. """ files, idx = _get_files(path, filestem, idl) tree = meson.rsplit('_')[0] corr_data = [] infos = [] for hd5_file in files: file = h5py.File(path + '/' + hd5_file, "r") raw_data = list(file[tree + '/' + meson + '/corr']) real_data = [o[0] for o in raw_data] corr_data.append(real_data) if not infos: for k, i in file[tree + '/' + meson].attrs.items(): infos.append(k + ': ' + i[0].decode()) file.close() corr_data = np.array(corr_data) l_obs = [] for c in corr_data.T: l_obs.append(Obs([c], [ens_id], idl=[idx])) corr = Corr(l_obs) corr.tag = r", ".join(infos) return corr class Npr_matrix(np.ndarray): def __new__(cls, input_array, mom_in=None, mom_out=None): obj = np.asarray(input_array).view(cls) obj.mom_in = mom_in obj.mom_out = mom_out return obj @property def g5H(self): """Gamma_5 hermitean conjugate Uses the fact that the propagator is gamma5 hermitean, so just the in and out momenta of the propagator are exchanged. """ return Npr_matrix(self, mom_in=self.mom_out, mom_out=self.mom_in) def _propagate_mom(self, other, name): s_mom = getattr(self, name, None) o_mom = getattr(other, name, None) if s_mom is not None and o_mom is not None: if not np.allclose(s_mom, o_mom): raise Exception(name + ' does not match.') return o_mom if o_mom is not None else s_mom def __matmul__(self, other): return self.__new__(Npr_matrix, super().__matmul__(other), self._propagate_mom(other, 'mom_in'), self._propagate_mom(other, 'mom_out')) def __array_finalize__(self, obj): if obj is None: return self.mom_in = getattr(obj, 'mom_in', None) self.mom_out = getattr(obj, 'mom_out', None) def read_ExternalLeg_hd5(path, filestem, ens_id, idl=None): """Read hadrons ExternalLeg hdf5 file and output an array of CObs Parameters ---------- path : str path to the files to read filestem : str namestem of the files to read ens_id : str name of the ensemble, required for internal bookkeeping idl : range If specified only configurations in the given range are read in. """ files, idx = _get_files(path, filestem, idl) mom = None corr_data = [] for hd5_file in files: file = h5py.File(path + '/' + hd5_file, "r") raw_data = file['ExternalLeg/corr'][0][0].view('complex') corr_data.append(raw_data) if mom is None: mom = np.array(str(file['ExternalLeg/info'].attrs['pIn'])[3:-2].strip().split(), dtype=float) file.close() corr_data = np.array(corr_data) rolled_array = np.rollaxis(corr_data, 0, 5) matrix = np.empty((rolled_array.shape[:-1]), dtype=object) for si, sj, ci, cj in np.ndindex(rolled_array.shape[:-1]): real = Obs([rolled_array[si, sj, ci, cj].real], [ens_id], idl=[idx]) imag = Obs([rolled_array[si, sj, ci, cj].imag], [ens_id], idl=[idx]) matrix[si, sj, ci, cj] = CObs(real, imag) return Npr_matrix(matrix, mom_in=mom) def read_Bilinear_hd5(path, filestem, ens_id, idl=None): """Read hadrons Bilinear hdf5 file and output an array of CObs Parameters ---------- path : str path to the files to read filestem : str namestem of the files to read ens_id : str name of the ensemble, required for internal bookkeeping idl : range If specified only configurations in the given range are read in. """ files, idx = _get_files(path, filestem, idl) mom_in = None mom_out = None corr_data = {} for hd5_file in files: file = h5py.File(path + '/' + hd5_file, "r") for i in range(16): name = file['Bilinear/Bilinear_' + str(i) + '/info'].attrs['gamma'][0].decode('UTF-8') if name not in corr_data: corr_data[name] = [] raw_data = file['Bilinear/Bilinear_' + str(i) + '/corr'][0][0].view('complex') corr_data[name].append(raw_data) if mom_in is None: mom_in = np.array(str(file['Bilinear/Bilinear_' + str(i) + '/info'].attrs['pIn'])[3:-2].strip().split(), dtype=float) if mom_out is None: mom_out = np.array(str(file['Bilinear/Bilinear_' + str(i) + '/info'].attrs['pOut'])[3:-2].strip().split(), dtype=float) file.close() result_dict = {} for key, data in corr_data.items(): local_data = np.array(data) rolled_array = np.rollaxis(local_data, 0, 5) matrix = np.empty((rolled_array.shape[:-1]), dtype=object) for si, sj, ci, cj in np.ndindex(rolled_array.shape[:-1]): real = Obs([rolled_array[si, sj, ci, cj].real], [ens_id], idl=[idx]) imag = Obs([rolled_array[si, sj, ci, cj].imag], [ens_id], idl=[idx]) matrix[si, sj, ci, cj] = CObs(real, imag) result_dict[key] = Npr_matrix(matrix, mom_in=mom_in, mom_out=mom_out) return result_dict def read_Fourquark_hd5(path, filestem, ens_id, idl=None, vertices=["VA", "AV"]): """Read hadrons FourquarkFullyConnected hdf5 file and output an array of CObs Parameters ---------- path : str path to the files to read filestem : str namestem of the files to read ens_id : str name of the ensemble, required for internal bookkeeping idl : range If specified only configurations in the given range are read in. vertices : list Vertex functions to be extracted. """ files, idx = _get_files(path, filestem, idl) mom_in = None mom_out = None vertex_names = [] for vertex in vertices: vertex_names += _get_lorentz_names(vertex) corr_data = {} tree = 'FourQuarkFullyConnected/FourQuarkFullyConnected_' for hd5_file in files: file = h5py.File(path + '/' + hd5_file, "r") for i in range(32): name = (file[tree + str(i) + '/info'].attrs['gammaA'][0].decode('UTF-8'), file[tree + str(i) + '/info'].attrs['gammaB'][0].decode('UTF-8')) if name in vertex_names: if name not in corr_data: corr_data[name] = [] raw_data = file[tree + str(i) + '/corr'][0][0].view('complex') corr_data[name].append(raw_data) if mom_in is None: mom_in = np.array(str(file[tree + str(i) + '/info'].attrs['pIn'])[3:-2].strip().split(), dtype=float) if mom_out is None: mom_out = np.array(str(file[tree + str(i) + '/info'].attrs['pOut'])[3:-2].strip().split(), dtype=float) file.close() intermediate_dict = {} for vertex in vertices: lorentz_names = _get_lorentz_names(vertex) for v_name in lorentz_names: if vertex not in intermediate_dict: intermediate_dict[vertex] = np.array(corr_data[v_name]) else: intermediate_dict[vertex] += np.array(corr_data[v_name]) result_dict = {} for key, data in intermediate_dict.items(): rolled_array = np.moveaxis(data, 0, 8) matrix = np.empty((rolled_array.shape[:-1]), dtype=object) for index in np.ndindex(rolled_array.shape[:-1]): real = Obs([rolled_array[index].real], [ens_id], idl=[idx]) imag = Obs([rolled_array[index].imag], [ens_id], idl=[idx]) matrix[index] = CObs(real, imag) result_dict[key] = Npr_matrix(matrix, mom_in=mom_in, mom_out=mom_out) return result_dict def _get_lorentz_names(name): assert len(name) == 2 res = [] if not set(name) <= set(['S', 'P', 'V', 'A', 'T']): raise Exception("Name can only contain 'S', 'P', 'V', 'A' or 'T'") if 'S' in name or 'P' in name: if not set(name) <= set(['S', 'P']): raise Exception("'" + name + "' is not a Lorentz scalar") g_names = {'S': 'Identity', 'P': 'Gamma5'} res.append((g_names[name[0]], g_names[name[1]])) elif 'T' in name: if not set(name) <= set(['T']): raise Exception("'" + name + "' is not a Lorentz scalar") raise Exception("Tensor operators not yet implemented.") else: if not set(name) <= set(['V', 'A']): raise Exception("'" + name + "' is not a Lorentz scalar") lorentz_index = ['X', 'Y', 'Z', 'T'] for ind in lorentz_index: res.append(('Gamma' + ind + (name[0] == 'A') * 'Gamma5', 'Gamma' + ind + (name[1] == 'A') * 'Gamma5')) return res
View Source
def read_meson_hd5(path, filestem, ens_id, meson='meson_0', tree='meson', idl=None): """Read hadrons meson hdf5 file and extract the meson labeled 'meson' Parameters ----------------- path : str path to the files to read filestem : str namestem of the files to read ens_id : str name of the ensemble, required for internal bookkeeping meson : str label of the meson to be extracted, standard value meson_0 which corresponds to the pseudoscalar pseudoscalar two-point function. idl : range If specified only configurations in the given range are read in. """ files, idx = _get_files(path, filestem, idl) tree = meson.rsplit('_')[0] corr_data = [] infos = [] for hd5_file in files: file = h5py.File(path + '/' + hd5_file, "r") raw_data = list(file[tree + '/' + meson + '/corr']) real_data = [o[0] for o in raw_data] corr_data.append(real_data) if not infos: for k, i in file[tree + '/' + meson].attrs.items(): infos.append(k + ': ' + i[0].decode()) file.close() corr_data = np.array(corr_data) l_obs = [] for c in corr_data.T: l_obs.append(Obs([c], [ens_id], idl=[idx])) corr = Corr(l_obs) corr.tag = r", ".join(infos) return corr
Read hadrons meson hdf5 file and extract the meson labeled 'meson'
Parameters
- path (str): path to the files to read
- filestem (str): namestem of the files to read
- ens_id (str): name of the ensemble, required for internal bookkeeping
- meson (str): label of the meson to be extracted, standard value meson_0 which corresponds to the pseudoscalar pseudoscalar two-point function.
- idl (range): If specified only configurations in the given range are read in.
View Source
class Npr_matrix(np.ndarray): def __new__(cls, input_array, mom_in=None, mom_out=None): obj = np.asarray(input_array).view(cls) obj.mom_in = mom_in obj.mom_out = mom_out return obj @property def g5H(self): """Gamma_5 hermitean conjugate Uses the fact that the propagator is gamma5 hermitean, so just the in and out momenta of the propagator are exchanged. """ return Npr_matrix(self, mom_in=self.mom_out, mom_out=self.mom_in) def _propagate_mom(self, other, name): s_mom = getattr(self, name, None) o_mom = getattr(other, name, None) if s_mom is not None and o_mom is not None: if not np.allclose(s_mom, o_mom): raise Exception(name + ' does not match.') return o_mom if o_mom is not None else s_mom def __matmul__(self, other): return self.__new__(Npr_matrix, super().__matmul__(other), self._propagate_mom(other, 'mom_in'), self._propagate_mom(other, 'mom_out')) def __array_finalize__(self, obj): if obj is None: return self.mom_in = getattr(obj, 'mom_in', None) self.mom_out = getattr(obj, 'mom_out', None)
ndarray(shape, dtype=float, buffer=None, offset=0, strides=None, order=None)
An array object represents a multidimensional, homogeneous array of fixed-size items. An associated data-type object describes the format of each element in the array (its byte-order, how many bytes it occupies in memory, whether it is an integer, a floating point number, or something else, etc.)
Arrays should be constructed using array
, zeros
or empty
(refer
to the See Also section below). The parameters given here refer to
a low-level method (ndarray(...)
) for instantiating an array.
For more information, refer to the numpy
module and examine the
methods and attributes of an array.
Parameters
- (for the __new__ method; see Notes below)
- shape (tuple of ints): Shape of created array.
- dtype (data-type, optional): Any object that can be interpreted as a numpy data type.
- buffer (object exposing buffer interface, optional): Used to fill the array with data.
- offset (int, optional): Offset of array data in buffer.
- strides (tuple of ints, optional): Strides of data in memory.
- order ({'C', 'F'}, optional): Row-major (C-style) or column-major (Fortran-style) order.
Attributes
- T (ndarray): Transpose of the array.
- data (buffer): The array's elements, in memory.
- dtype (dtype object): Describes the format of the elements in the array.
- flags (dict): Dictionary containing information related to memory use, e.g., 'C_CONTIGUOUS', 'OWNDATA', 'WRITEABLE', etc.
- flat (numpy.flatiter object):
Flattened version of the array as an iterator. The iterator
allows assignments, e.g.,
x.flat = 3
(Seendarray.flat
for assignment examples; TODO). - imag (ndarray): Imaginary part of the array.
- real (ndarray): Real part of the array.
- size (int): Number of elements in the array.
- itemsize (int): The memory use of each array element in bytes.
- nbytes (int):
The total number of bytes required to store the array data,
i.e.,
itemsize * size
. - ndim (int): The array's number of dimensions.
- shape (tuple of ints): Shape of the array.
- strides (tuple of ints):
The step-size required to move from one element to the next in
memory. For example, a contiguous
(3, 4)
array of typeint16
in C-order has strides(8, 2)
. This implies that to move from element to element in memory requires jumps of 2 bytes. To move from row-to-row, one needs to jump 8 bytes at a time (2 * 4
). - ctypes (ctypes object): Class containing properties of the array needed for interaction with ctypes.
- base (ndarray):
If the array is a view into another array, that array is its
base
(unless that array is also a view). Thebase
array is where the array data is actually stored.
See Also
array
: Construct an array.
zeros
: Create an array, each element of which is zero.
empty
: Create an array, but leave its allocated memory unchanged (i.e.,
it contains "garbage").
dtype
: Create a data-type.
numpy.typing.NDArray
: An ndarray alias :term:generic <generic type>
w.r.t. its dtype.type <numpy.dtype.type>
.
Notes
There are two modes of creating an array using __new__
:
- If
buffer
is None, then onlyshape
,dtype
, andorder
are used. - If
buffer
is an object exposing the buffer interface, then all keywords are interpreted.
No __init__
method is needed because the array is fully initialized
after the __new__
method.
Examples
These examples illustrate the low-level ndarray
constructor. Refer
to the See Also
section above for easier ways of constructing an
ndarray.
First mode, buffer
is None:
>>> np.ndarray(shape=(2,2), dtype=float, order='F')
array([[0.0e+000, 0.0e+000], # random
[ nan, 2.5e-323]])
Second mode:
>>> np.ndarray((2,), buffer=np.array([1,2,3]),
... offset=np.int_().itemsize,
... dtype=int) # offset = 1*itemsize, i.e. skip first element
array([2, 3])
Gamma_5 hermitean conjugate
Uses the fact that the propagator is gamma5 hermitean, so just the in and out momenta of the propagator are exchanged.
Inherited Members
- numpy.ndarray
- dumps
- dump
- all
- any
- argmax
- argmin
- argpartition
- argsort
- astype
- byteswap
- choose
- clip
- compress
- conj
- conjugate
- copy
- cumprod
- cumsum
- diagonal
- dot
- fill
- flatten
- getfield
- item
- itemset
- max
- mean
- min
- newbyteorder
- nonzero
- partition
- prod
- ptp
- put
- ravel
- repeat
- reshape
- resize
- round
- searchsorted
- setfield
- setflags
- sort
- squeeze
- std
- sum
- swapaxes
- take
- tobytes
- tofile
- tolist
- tostring
- trace
- transpose
- var
- view
- ndim
- flags
- shape
- strides
- data
- itemsize
- size
- nbytes
- base
- dtype
- real
- imag
- flat
- ctypes
- T
View Source
def read_ExternalLeg_hd5(path, filestem, ens_id, idl=None): """Read hadrons ExternalLeg hdf5 file and output an array of CObs Parameters ---------- path : str path to the files to read filestem : str namestem of the files to read ens_id : str name of the ensemble, required for internal bookkeeping idl : range If specified only configurations in the given range are read in. """ files, idx = _get_files(path, filestem, idl) mom = None corr_data = [] for hd5_file in files: file = h5py.File(path + '/' + hd5_file, "r") raw_data = file['ExternalLeg/corr'][0][0].view('complex') corr_data.append(raw_data) if mom is None: mom = np.array(str(file['ExternalLeg/info'].attrs['pIn'])[3:-2].strip().split(), dtype=float) file.close() corr_data = np.array(corr_data) rolled_array = np.rollaxis(corr_data, 0, 5) matrix = np.empty((rolled_array.shape[:-1]), dtype=object) for si, sj, ci, cj in np.ndindex(rolled_array.shape[:-1]): real = Obs([rolled_array[si, sj, ci, cj].real], [ens_id], idl=[idx]) imag = Obs([rolled_array[si, sj, ci, cj].imag], [ens_id], idl=[idx]) matrix[si, sj, ci, cj] = CObs(real, imag) return Npr_matrix(matrix, mom_in=mom)
Read hadrons ExternalLeg hdf5 file and output an array of CObs
Parameters
- path (str): path to the files to read
- filestem (str): namestem of the files to read
- ens_id (str): name of the ensemble, required for internal bookkeeping
- idl (range): If specified only configurations in the given range are read in.
View Source
def read_Bilinear_hd5(path, filestem, ens_id, idl=None): """Read hadrons Bilinear hdf5 file and output an array of CObs Parameters ---------- path : str path to the files to read filestem : str namestem of the files to read ens_id : str name of the ensemble, required for internal bookkeeping idl : range If specified only configurations in the given range are read in. """ files, idx = _get_files(path, filestem, idl) mom_in = None mom_out = None corr_data = {} for hd5_file in files: file = h5py.File(path + '/' + hd5_file, "r") for i in range(16): name = file['Bilinear/Bilinear_' + str(i) + '/info'].attrs['gamma'][0].decode('UTF-8') if name not in corr_data: corr_data[name] = [] raw_data = file['Bilinear/Bilinear_' + str(i) + '/corr'][0][0].view('complex') corr_data[name].append(raw_data) if mom_in is None: mom_in = np.array(str(file['Bilinear/Bilinear_' + str(i) + '/info'].attrs['pIn'])[3:-2].strip().split(), dtype=float) if mom_out is None: mom_out = np.array(str(file['Bilinear/Bilinear_' + str(i) + '/info'].attrs['pOut'])[3:-2].strip().split(), dtype=float) file.close() result_dict = {} for key, data in corr_data.items(): local_data = np.array(data) rolled_array = np.rollaxis(local_data, 0, 5) matrix = np.empty((rolled_array.shape[:-1]), dtype=object) for si, sj, ci, cj in np.ndindex(rolled_array.shape[:-1]): real = Obs([rolled_array[si, sj, ci, cj].real], [ens_id], idl=[idx]) imag = Obs([rolled_array[si, sj, ci, cj].imag], [ens_id], idl=[idx]) matrix[si, sj, ci, cj] = CObs(real, imag) result_dict[key] = Npr_matrix(matrix, mom_in=mom_in, mom_out=mom_out) return result_dict
Read hadrons Bilinear hdf5 file and output an array of CObs
Parameters
- path (str): path to the files to read
- filestem (str): namestem of the files to read
- ens_id (str): name of the ensemble, required for internal bookkeeping
- idl (range): If specified only configurations in the given range are read in.
View Source
def read_Fourquark_hd5(path, filestem, ens_id, idl=None, vertices=["VA", "AV"]): """Read hadrons FourquarkFullyConnected hdf5 file and output an array of CObs Parameters ---------- path : str path to the files to read filestem : str namestem of the files to read ens_id : str name of the ensemble, required for internal bookkeeping idl : range If specified only configurations in the given range are read in. vertices : list Vertex functions to be extracted. """ files, idx = _get_files(path, filestem, idl) mom_in = None mom_out = None vertex_names = [] for vertex in vertices: vertex_names += _get_lorentz_names(vertex) corr_data = {} tree = 'FourQuarkFullyConnected/FourQuarkFullyConnected_' for hd5_file in files: file = h5py.File(path + '/' + hd5_file, "r") for i in range(32): name = (file[tree + str(i) + '/info'].attrs['gammaA'][0].decode('UTF-8'), file[tree + str(i) + '/info'].attrs['gammaB'][0].decode('UTF-8')) if name in vertex_names: if name not in corr_data: corr_data[name] = [] raw_data = file[tree + str(i) + '/corr'][0][0].view('complex') corr_data[name].append(raw_data) if mom_in is None: mom_in = np.array(str(file[tree + str(i) + '/info'].attrs['pIn'])[3:-2].strip().split(), dtype=float) if mom_out is None: mom_out = np.array(str(file[tree + str(i) + '/info'].attrs['pOut'])[3:-2].strip().split(), dtype=float) file.close() intermediate_dict = {} for vertex in vertices: lorentz_names = _get_lorentz_names(vertex) for v_name in lorentz_names: if vertex not in intermediate_dict: intermediate_dict[vertex] = np.array(corr_data[v_name]) else: intermediate_dict[vertex] += np.array(corr_data[v_name]) result_dict = {} for key, data in intermediate_dict.items(): rolled_array = np.moveaxis(data, 0, 8) matrix = np.empty((rolled_array.shape[:-1]), dtype=object) for index in np.ndindex(rolled_array.shape[:-1]): real = Obs([rolled_array[index].real], [ens_id], idl=[idx]) imag = Obs([rolled_array[index].imag], [ens_id], idl=[idx]) matrix[index] = CObs(real, imag) result_dict[key] = Npr_matrix(matrix, mom_in=mom_in, mom_out=mom_out) return result_dict
Read hadrons FourquarkFullyConnected hdf5 file and output an array of CObs
Parameters
- path (str): path to the files to read
- filestem (str): namestem of the files to read
- ens_id (str): name of the ensemble, required for internal bookkeeping
- idl (range): If specified only configurations in the given range are read in.
- vertices (list): Vertex functions to be extracted.