def unpack_4byte_integer(file, count, endian='>'): """ Unpacks 4 byte integers. """ # Read as 4 byte integer so bit shifting works. data = from_buffer(file.read(count * 4), dtype=np.int32) # Swap the byte order if necessary. if BYTEORDER != endian: data = data.byteswap() return data
def unpack_4byte_ibm(file, count, endian='>'): """ Unpacks 4 byte IBM floating points. """ # Read as 4 byte integer so bit shifting works. data = from_buffer(file.read(count * 4), dtype=np.float32) # Swap the byte order if necessary. if BYTEORDER != endian: data = data.byteswap() length = len(data) # Call the C code which transforms the data inplace. clibsegy.ibm2ieee(data, length) return data
def parse_next_trace(self): """ Parse the next trace in the trace pointer list and return a Trace object. """ trace_descriptor_block = self.file_pointer.read(32) # Check if the trace descriptor block id is valid. if unpack(self.endian + b'H', trace_descriptor_block[0:2])[0] != \ 0x4422: msg = 'Invalid trace descriptor block id.' raise SEG2InvalidFileError(msg) size_of_this_block, = unpack_from(self.endian + b'H', trace_descriptor_block, 2) number_of_samples_in_data_block, = \ unpack_from(self.endian + b'L', trace_descriptor_block, 8) data_format_code, = unpack_from(b'B', trace_descriptor_block, 12) # Parse the data format code. if data_format_code == 4: dtype = self.endian + b'f4' sample_size = 4 elif data_format_code == 5: dtype = self.endian + b'f8' sample_size = 8 elif data_format_code == 1: dtype = self.endian + b'i2' sample_size = 2 elif data_format_code == 2: dtype = self.endian + b'i4' sample_size = 4 elif data_format_code == 3: dtype = self.endian + b'i2' sample_size = 2.5 if number_of_samples_in_data_block % 4 != 0: raise SEG2InvalidFileError( 'Data format code 3 requires that the number of samples ' 'is divisible by 4, but sample count is %d' % (number_of_samples_in_data_block, )) else: msg = 'Unrecognized data format code' raise SEG2InvalidFileError(msg) # The rest of the trace block is free form. header = {} header['seg2'] = AttribDict() self.parse_free_form(self.file_pointer.read(size_of_this_block - 32), header['seg2']) header['delta'] = float(header['seg2']['SAMPLE_INTERVAL']) # Set to the file's start time. header['starttime'] = deepcopy(self.starttime) if 'DELAY' in header['seg2']: if float(header['seg2']['DELAY']) != 0: msg = "Non-zero value found in Trace's 'DELAY' field. " + \ "This is not supported/tested yet and might lead " + \ "to a wrong starttime of the Trace. Please contact " + \ "the ObsPy developers with a sample file." warnings.warn(msg) if "DESCALING_FACTOR" in header["seg2"]: header['calib'] = float(header['seg2']['DESCALING_FACTOR']) # Unpack the data. data = from_buffer(self.file_pointer.read( int(number_of_samples_in_data_block * sample_size)), dtype=dtype) if data_format_code == 3: # Convert one's complement to two's complement by adding one to # negative numbers. one_to_two = (data < 0) # The first two bytes (1 word) of every 10 bytes (5 words) contains # a 4-bit exponent for each of the 4 remaining 2-byte (int16) # samples. exponents = data[0::5].view(self.endian + b'u2') result = np.empty(number_of_samples_in_data_block, dtype=np.int32) # Apply the negative correction, then multiply by correct exponent. result[0::4] = ((data[1::5] + one_to_two[1::5]) * 2**((exponents & 0x000f) >> 0)) result[1::4] = ((data[2::5] + one_to_two[2::5]) * 2**((exponents & 0x00f0) >> 4)) result[2::4] = ((data[3::5] + one_to_two[3::5]) * 2**((exponents & 0x0f00) >> 8)) result[3::4] = ((data[4::5] + one_to_two[4::5]) * 2**((exponents & 0xf000) >> 12)) data = result # Integrate SEG2 file header into each trace header tmp = self.stream.stats.seg2.copy() tmp.update(header['seg2']) header['seg2'] = tmp return Trace(data=data, header=header)