예제 #1
0
    def read(self,arffile):
        """Read the effective area from an OGIP ARF file."""

        (data, header) = fits.getdata(arffile,'SPECRESP',header=True)

        self.LowEnergy = data['ENERG_LO']
        self.HighEnergy = data['ENERG_HI']
        self.EffArea = data['SPECRESP']

        self.EnergyUnits = header['TUNIT1']

        if header['TUNIT3'] == 'cm**2':
            self.ARFUnits = 'cm2'
        elif header['TUNIT3'] == 'cm2':
            self.ARFUnits = 'cm2'
        else:
            message.warning("ARF units are not recognized.")

        try:
            self.Order = header['TG_M']
            self.Grating = header['TG_PART']
        except:
            self.Order = 0
            self.Grating = 0

        # Check for NULL values
        nans = np.isnan(self.EffArea)
        if np.any(nans):
            for i in np.arange(self.EffArea.size):
                if nans[i]:
                   self.EffArea[i] = 0.0

        return 0
예제 #2
0
    def read_header(self, header):
        """Utility function to read the header from a SPECTRUM extension for both PHA and PHAII files.

        :param header: Header of the SPECTRUM extension.
        :type header: astropy.io.fits.Header
        """

        # Read Exposure information
        self.Exposure = header['EXPOSURE']

        # Read how the spectrum is stored (COUNTS or RATE)
        try:
            self.Spectrumtype = header['HDUCLAS2']
        except KeyError:
            self.Spectrumtype = 'TOTAL'
            message.warning(
                "HDUCLAS2 keyword not found. Assuming spectrumtype is TOTAL.")

        try:
            self.PhaType = header['HDUCLAS3']
        except KeyError:
            self.PhaType = 'COUNTS'
            message.warning(
                "HDUCLAS3 keyword not found. Assuming PHA type is COUNTS.")

        # Read the POISERR keyword
        try:
            self.Poisserr = header['POISSERR']
        except KeyError:
            self.Poisserr = False

        # Read Correction scaling factor
        self.CorrScal = header['CORRSCAL']

        # Read a background file, if available
        try:
            self.bkgfile = header['BACKFILE']
        except KeyError:
            self.bkgfile = None

        # Read an respoonse file, if available
        try:
            self.rmffile = header['RESPFILE']
        except KeyError:
            self.rmffile = None

        # Read an effective area file, if available
        try:
            self.arffile = header['ANCRFILE']
        except KeyError:
            self.arffile = None

        # Read a correction file, if available
        try:
            self.corfile = header['CORRFILE']
        except KeyError:
            self.corfile = None
예제 #3
0
    def check_filename(self,filename):
        """Check if the output filename has the correct .res extension. The method returns a correct file name."""
        resname, res_extension = os.path.splitext(filename)
        if res_extension != ".res":
            message.warning("Output filename does not have the correct .res extension.")
            print("Renaming file to end with '.res'.")
            print("")
            resfile = resname + '.res'
        else:
            resfile = filename

        return resfile
예제 #4
0
    def check_filename(self, filename):
        """Check if the output filename has the correct .spo extension. The method returns a correct file name."""
        sponame, spo_extension = os.path.splitext(filename)
        if spo_extension != ".spo":
            message.warning(
                "Output filename does not have the correct .spo extension.")
            print("Renaming file to end with '.spo'.")
            print("")
            spofile = sponame + '.spo'
        else:
            spofile = filename

        return spofile
예제 #5
0
    def write(self, rmffile, telescop=None, instrume=None, filterkey=None, overwrite=False):
        """Method to write an OGIP format RMF file.

        :param rmffile: RMF file name to write.
        :type rmffile: str
        :param telescop: Name of the telescope to be put in the TELESCOP keyword.
        :type telescop: str
        :param instrume: Name of the instrument to be put in the INSTRUME keyword.
        :type instrume: str
        :param filterkey: Name of the filter to be put in the FILTER keyword.
        :type filterkey: str
        :param overwrite: Overwrite existing file names? (True/False)
        :type overwrite: bool
        """

        #
        # Generate warning if there are multiple groups per energy
        #
        if np.amax(self.matrix[0].NumberGroups) != 1:
            message.warning("This method has not been tested for responses with multiple response groups per energy.")

        #
        # Create Primary HDU
        #
        primary = fits.PrimaryHDU()

        #
        # Create the EBOUNDS extension
        #
        ecol1 = fits.Column(name='CHANNEL', format='J', array=self.ebounds.Channel)
        ecol2 = fits.Column(name='E_MIN', format='D', unit=self.ebounds.EnergyUnits, array=self.ebounds.ChannelLowEnergy)
        ecol3 = fits.Column(name='E_MAX', format='D', unit=self.ebounds.EnergyUnits, array=self.ebounds.ChannelHighEnergy)

        ebnds = fits.BinTableHDU.from_columns([ecol1, ecol2, ecol3])

        ehdr = ebnds.header
        ehdr.set('EXTNAME', 'EBOUNDS')
        ehdr.set('DETCHANS', self.ebounds.NumberChannels)

        # Set the TELESCOP keyword (optional)
        if telescop is None:
            ehdr.set('TELESCOP', 'None', 'Telescope name')
        else:
            ehdr.set('TELESCOP', telescop, 'Telescope name')

        # Set the INSTRUME keyword (optional)
        if instrume is None:
            ehdr.set('INSTRUME', 'None', 'Instrument name')
        else:
            ehdr.set('INSTRUME', instrume, 'Instrument name')

        # Set the FILTER keyword (optional)
        if filterkey is None:
            ehdr.set('FILTER', 'None', 'Filter setting')
        else:
            ehdr.set('FILTER', filterkey, 'Filter setting')

        ehdr.set('DETNAM ', 'None')
        ehdr.set('CHANTYPE', 'PI')
        ehdr.set('HDUCLASS', 'OGIP')
        ehdr.set('HDUCLAS1', 'RESPONSE')
        ehdr.set('HDUCLAS2', 'EBOUNDS ')
        ehdr.set('HDUVERS1', '1.2.0')
        ehdr.set('ORIGIN ', 'SRON')

        hdu = fits.HDUList([primary, ebnds])

        #
        # Create SPECRESP MATRIX extension
        #
        for e in range(self.NumberMatrixExt):
            print("Writing matrix for matrix number: {0}".format(e))

            mcol1 = fits.Column(name='ENERG_LO', format='D', unit=self.matrix[e].EnergyUnits, array=self.matrix[e].LowEnergy)
            mcol2 = fits.Column(name='ENERG_HI', format='D', unit=self.matrix[e].EnergyUnits, array=self.matrix[e].HighEnergy)
            mcol3 = fits.Column(name='N_GRP', format='J', array=self.matrix[e].NumberGroups)
            mcol4 = fits.Column(name='F_CHAN', format='J', array=self.matrix[e].FirstChannelGroup)
            mcol5 = fits.Column(name='N_CHAN', format='J', array=self.matrix[e].NumberChannelsGroup)

            # Determine the width of the matrix
            width = np.amax(self.matrix[e].NumberChannelsGroup)
            formatstr = str(width)+'D'

            #
            # THIS PART COULD BE UPDATED TO OPTIMIZE THE SIZE USING VARIABLE SIZE ARRAYS IN FITS.
            #

            # Building the MATRIX column
            newmatrix = np.zeros(self.matrix[e].NumberEnergyBins * width, dtype=float).reshape(self.matrix[e].NumberEnergyBins, width)

            re = 0
            for i in np.arange(self.matrix[e].NumberEnergyBins):
                for j in np.arange(self.matrix[e].NumberGroups[i]):
                    for k in np.arange(self.matrix[e].NumberChannelsGroup[i]):
                        newmatrix[i, k] = self.matrix[e].Matrix[re]
                        re = re + 1

            mcol6 = fits.Column(name='MATRIX', format=formatstr, array=newmatrix)

            matrix = fits.BinTableHDU.from_columns([mcol1, mcol2, mcol3, mcol4, mcol5, mcol6])

            mhdr = matrix.header

            if self.matrix[e].AreaIncluded:
                mhdr.set('EXTNAME', 'SPECRESP MATRIX')
            else:
                mhdr.set('EXTNAME', 'MATRIX')

            # Set the TELESCOP keyword (optional)
            if telescop is None:
                mhdr.set('TELESCOP', 'None', 'Telescope name')
            else:
                mhdr.set('TELESCOP', telescop, 'Telescope name')

            # Set the INSTRUME keyword (optional)
            if instrume is None:
                mhdr.set('INSTRUME', 'None', 'Instrument name')
            else:
                mhdr.set('INSTRUME', instrume, 'Instrument name')

            # Set the FILTER keyword (optional)
            if filterkey is None:
                mhdr.set('FILTER', 'None', 'Filter setting')
            else:
                mhdr.set('FILTER', filterkey, 'Filter setting')

            mhdr.set('DETCHANS', self.ebounds.NumberChannels)
            mhdr.set('LO_THRES', self.matrix[e].ResponseThreshold)

            mhdr.set('CHANTYPE', 'PI')
            mhdr.set('HDUCLASS', 'OGIP')
            mhdr.set('HDUCLAS1', 'RESPONSE')
            mhdr.set('HDUCLAS2', 'RSP_MATRIX')

            if self.matrix[e].AreaIncluded:
                mhdr.set('HDUCLAS3', 'FULL')
            else:
                mhdr.set('HDUCLAS3', 'REDIST')
            mhdr.set('HDUVERS1', '1.3.0')
            mhdr.set('ORIGIN  ', 'SRON')

            matrix.header['HISTORY'] = 'Created by pyspextools:'
            matrix.header['HISTORY'] = 'https://github.com/spex-xray/pyspextools'

            hdu.append(matrix)

        try:
            hdu.writeto(rmffile, overwrite=overwrite)
        except IOError:
            message.error("File {0} already exists. I will not overwrite it!".format(rmffile))
            return 1

        return 0
예제 #6
0
    def read(self, rmfhdu):

        # Read the Matrix table
        data = rmfhdu.data
        header = rmfhdu.header

        if rmfhdu.name == 'MATRIX':
            pass
        elif rmfhdu.name == 'SPECRESP MATRIX':
            message.warning("This is an RSP file with the effective area included.")
            print("Do not read an ARF file, unless you know what you are doing.")
            self.AreaIncluded = True
        else:
            message.error("MATRIX extension not successfully found in RMF file.")
            return

        self.LowEnergy = data['ENERG_LO']
        self.HighEnergy = data['ENERG_HI']
        self.NumberEnergyBins = self.LowEnergy.size
        self.EnergyUnits = header['TUNIT1']

        self.NumberGroups = data['N_GRP']
        self.NumberTotalGroups = np.sum(self.NumberGroups)

        self.FirstGroup = np.zeros(self.NumberEnergyBins, dtype=int)

        self.FirstChannelGroup = np.zeros(self.NumberTotalGroups, dtype=int)
        self.NumberChannelsGroup = np.zeros(self.NumberTotalGroups, dtype=int)
        self.FirstElement = np.zeros(self.NumberTotalGroups, dtype=int)

        self.Matrix = np.array([], dtype=float)

        try:
            self.Order = header['ORDER']
        except KeyError:
            pass

        fgroup = 0  # Count total number of groups
        felem = 0   # Count total number of response elements
        nelem = np.zeros(self.NumberEnergyBins, dtype=int)   # Count number of response elements per energy bin
        k = 0

        fchan_local = data['F_CHAN']
        nchan_local = data['N_CHAN']
        matrix_local = data['MATRIX']

        for i in np.arange(self.NumberEnergyBins):
            self.FirstGroup[i] = fgroup
            fgroup = fgroup + self.NumberGroups[i]

            if self.NumberGroups[i] != 0:
                for j in np.arange(self.NumberGroups[i]):
                    try:
                        self.FirstChannelGroup[k] = fchan_local[i][j]
                        self.NumberChannelsGroup[k] = nchan_local[i][j]
                    except IndexError:
                        self.FirstChannelGroup[k] = fchan_local[i]
                        self.NumberChannelsGroup[k] = nchan_local[i]
                    self.FirstElement[k] = felem
                    felem = felem + self.NumberChannelsGroup[k]
                    nelem[i] = nelem[i] + self.NumberChannelsGroup[k]

                    k = k + 1

        self.Matrix = np.zeros(felem, dtype=float)

        r = 0
        for i in np.arange(self.NumberEnergyBins):
            if nelem[i] != 0:
                for j in np.arange(nelem[i]):
                    self.Matrix[r] = matrix_local[i][j]
                    r = r + 1

        self.NumberTotalElements = self.Matrix.size
        self.ResponseThreshold = np.amin(self.Matrix)
예제 #7
0
    def read(self, filename):
        """Read a spectrum from a PHA file."""

        # Read the data and header from the SPECTRUM extension
        (data, header) = fits.getdata(filename, 'SPECTRUM', header=True)

        # Read the number of channels (outside the header call because of PHAII files).
        self.DetChans = header['DETCHANS']

        # Read the header
        self.read_header(header)

        # Read Channel information
        self.Channel = data['CHANNEL']
        self.FirstChannel = self.Channel[0]

        # Read the spectrum and convert to rate if necessary
        if self.PhaType == 'RATE':
            self.Rate = data['RATE']
        else:
            self.Rate = np.zeros(self.DetChans, dtype=float)
            for i in np.arange(self.DetChans):
                self.Rate[i] = float(data['COUNTS'][i]) / self.Exposure

        # See if there are Statistical Errors present
        if not self.Poiserr:
            try:
                self.StatError = data['STAT_ERR']
            except:
                self.StatError = None
                message.warning(
                    "No Poisson errors, but no STAT_ERR keyword found.")
        else:
            self.StatError = np.zeros(self.DetChans, dtype=float)
            for i in np.arange(self.DetChans):
                self.StatError[i] = math.sqrt(self.Rate[i] / self.Exposure)

        # Are there systematic errors?
        try:
            self.SysError = data['SYS_ERR']
        except:
            self.SysError = np.zeros(self.DetChans, dtype=float)

        if self.PhaType == 'RATE':
            self.SysError = self.SysError / self.Exposure

        # Are there quality flags?
        try:
            self.Quality = data['QUALITY']
        except:
            self.Quality = np.zeros(self.DetChans, dtype=int)

        # Are there grouping flags?
        try:
            self.Grouping = data['GROUPING']
        except:
            self.Grouping = np.zeros(self.DetChans, dtype=int)

        # Is there a backscale column?
        try:
            self.BackScaling = data['BACKSCAL']
        except:
            self.BackScaling = np.ones(self.DetChans,
                                       dtype=float) * header['BACKSCAL']

        # Is there an areascale column?
        try:
            self.AreaScaling = data['AREASCAL']
        except:
            self.AreaScaling = np.ones(self.DetChans,
                                       dtype=float) * header['AREASCAL']

        return 0
예제 #8
0
    def read(self, rmffile):
        """Method to read OGIP RMF files. The variable naming is made consistent with the HEASOFT HEASP module by
        Keith Arnaud."""

        # Read the Ebounds table
        (data, header) = fits.getdata(rmffile, 'EBOUNDS', header=True)

        self.Channel = data['CHANNEL']
        self.ChannelLowEnergy = data['E_MIN']
        self.ChannelHighEnergy = data['E_MAX']
        self.NumberChannels = self.Channel.size
        self.FirstChannel = self.Channel[0]

        # Read the Matrix table
        rmf = fits.open(rmffile)

        try:
            data = rmf['MATRIX'].data
            header = rmf['MATRIX'].header
        except:
            data = rmf['SPECRESP MATRIX'].data
            header = rmf['SPECRESP MATRIX'].header
            message.warning(
                "This is an RSP file with the effective area included.")
            print(
                "Do not read an ARF file, unless you know what you are doing.")
            self.AreaIncluded = True

        self.LowEnergy = data['ENERG_LO']
        self.HighEnergy = data['ENERG_HI']
        self.NumberEnergyBins = self.LowEnergy.size
        self.EnergyUnits = header['TUNIT1']

        self.NumberGroups = data['N_GRP']
        self.NumberTotalGroups = np.sum(self.NumberGroups)

        self.FirstGroup = np.zeros(self.NumberEnergyBins, dtype=int)

        self.FirstChannelGroup = np.zeros(self.NumberTotalGroups, dtype=int)
        self.NumberChannelsGroup = np.zeros(self.NumberTotalGroups, dtype=int)
        self.FirstElement = np.zeros(self.NumberTotalGroups, dtype=int)

        self.Matrix = np.array([], dtype=float)

        try:
            self.Order = header['ORDER']
        except:
            pass

        fgroup = 0  # Count total number of groups
        felem = 0  # Count total number of response elements
        nelem = np.zeros(
            self.NumberEnergyBins,
            dtype=int)  # Count number of response elements per energy bin
        k = 0

        fchan_local = data['F_CHAN']
        nchan_local = data['N_CHAN']
        matrix_local = data['MATRIX']

        for i in np.arange(self.NumberEnergyBins):
            self.FirstGroup[i] = fgroup
            fgroup = fgroup + self.NumberGroups[i]

            if self.NumberGroups[i] != 0:
                for j in np.arange(self.NumberGroups[i]):
                    try:
                        self.FirstChannelGroup[k] = fchan_local[i][j]
                        self.NumberChannelsGroup[k] = nchan_local[i][j]
                    except IndexError:
                        self.FirstChannelGroup[k] = fchan_local[i]
                        self.NumberChannelsGroup[k] = nchan_local[i]
                    self.FirstElement[k] = felem
                    felem = felem + self.NumberChannelsGroup[k]
                    nelem[i] = nelem[i] + self.NumberChannelsGroup[k]

                    k = k + 1

        self.Matrix = np.zeros(felem, dtype=float)

        r = 0
        for i in np.arange(self.NumberEnergyBins):
            if nelem[i] != 0:
                for j in np.arange(nelem[i]):
                    self.Matrix[r] = matrix_local[i][j]
                    r = r + 1

        self.NumberTotalElements = self.Matrix.size
        self.ResponseThreshold = np.amin(self.Matrix)

        rmf.close()

        return 0
예제 #9
0
    def write(self,
              rmffile,
              telescop=None,
              instrume=None,
              filter=None,
              overwrite=False):
        '''Method to write an OGIP format RMF file.'''

        #
        # Generate warning if there are multiple groups per energy
        #
        if (np.amax(self.NumberGroups) != 1):
            message.warning(
                "This method has not been tested for responses with multiple response groups per energy."
            )

        #
        # Create Primary HDU
        #
        primary = fits.PrimaryHDU()

        #
        # Create the EBOUNDS extension
        #
        ecol1 = fits.Column(name='CHANNEL', format='J', array=self.Channel)
        ecol2 = fits.Column(name='E_MIN',
                            format='D',
                            unit=self.EnergyUnits,
                            array=self.ChannelLowEnergy)
        ecol3 = fits.Column(name='E_MAX',
                            format='D',
                            unit=self.EnergyUnits,
                            array=self.ChannelHighEnergy)

        ebounds = fits.BinTableHDU.from_columns([ecol1, ecol2, ecol3])

        ehdr = ebounds.header
        ehdr.set('EXTNAME', 'EBOUNDS')
        ehdr.set('DETCHANS', self.NumberChannels)

        # Set the TELESCOP keyword (optional)
        if telescop == None:
            ehdr.set('TELESCOP', 'None', 'Telescope name')
        else:
            ehdr.set('TELESCOP', telescop, 'Telescope name')

        # Set the INSTRUME keyword (optional)
        if instrume == None:
            ehdr.set('INSTRUME', 'None', 'Instrument name')
        else:
            ehdr.set('INSTRUME', instrume, 'Instrument name')

        # Set the FILTER keyword (optional)
        if filter == None:
            ehdr.set('FILTER', 'None', 'Filter setting')
        else:
            ehdr.set('FILTER', filter, 'Filter setting')

        ehdr.set('DETNAM ', 'None')
        ehdr.set('CHANTYPE', 'PI')
        ehdr.set('HDUCLASS', 'OGIP')
        ehdr.set('HDUCLAS1', 'RESPONSE')
        ehdr.set('HDUCLAS2', 'EBOUNDS ')
        ehdr.set('HDUVERS1', '1.2.0')
        ehdr.set('ORIGIN ', 'SRON')

        #
        # Create SPECRESP MATRIX extension
        #
        mcol1 = fits.Column(name='ENERG_LO',
                            format='D',
                            unit=self.EnergyUnits,
                            array=self.LowEnergy)
        mcol2 = fits.Column(name='ENERG_HI',
                            format='D',
                            unit=self.EnergyUnits,
                            array=self.HighEnergy)
        mcol3 = fits.Column(name='N_GRP', format='J', array=self.NumberGroups)
        mcol4 = fits.Column(name='F_CHAN',
                            format='J',
                            array=self.FirstChannelGroup)
        mcol5 = fits.Column(name='N_CHAN',
                            format='J',
                            array=self.NumberChannelsGroup)

        # Determine the width of the matrix
        width = np.amax(self.NumberChannelsGroup)
        format = str(width) + 'D'

        # Building the MATRIX column
        newmatrix = np.zeros(self.Matrix.size,
                             dtype=float).reshape(self.NumberEnergyBins, width)

        re = 0
        for i in np.arange(self.NumberEnergyBins):
            for j in np.arange(self.NumberGroups[i]):
                for k in np.arange(self.NumberChannelsGroup[i]):
                    newmatrix[i, k] = self.Matrix[re]
                    re = re + 1

        mcol6 = fits.Column(name='MATRIX', format=format, array=newmatrix)

        matrix = fits.BinTableHDU.from_columns(
            [mcol1, mcol2, mcol3, mcol4, mcol5, mcol6])

        mhdr = matrix.header

        if self.AreaIncluded:
            mhdr.set('EXTNAME', 'SPECRESP MATRIX')
        else:
            mhdr.set('EXTNAME', 'MATRIX')

        # Set the TELESCOP keyword (optional)
        if telescop == None:
            mhdr.set('TELESCOP', 'None', 'Telescope name')
        else:
            mhdr.set('TELESCOP', telescop, 'Telescope name')

        # Set the INSTRUME keyword (optional)
        if instrume == None:
            mhdr.set('INSTRUME', 'None', 'Instrument name')
        else:
            mhdr.set('INSTRUME', instrume, 'Instrument name')

        # Set the FILTER keyword (optional)
        if filter == None:
            mhdr.set('FILTER', 'None', 'Filter setting')
        else:
            mhdr.set('FILTER', filter, 'Filter setting')

        mhdr.set('DETCHANS', self.NumberChannels)
        mhdr.set('LO_THRES', self.ResponseThreshold)

        mhdr.set('CHANTYPE', 'PI')
        mhdr.set('HDUCLASS', 'OGIP')
        mhdr.set('HDUCLAS1', 'RESPONSE')
        mhdr.set('HDUCLAS2', 'RSP_MATRIX')

        if self.AreaIncluded:
            mhdr.set('HDUCLAS3', 'FULL')
        else:
            mhdr.set('HDUCLAS3', 'REDIST')
        mhdr.set('HDUVERS1', '1.3.0')
        mhdr.set('ORIGIN  ', 'SRON')

        matrix.header['HISTORY'] = 'Created by pyspextools:'
        matrix.header['HISTORY'] = 'https://github.com/spex-xray/pyspextools'

        #
        # Final HDU List
        #
        hdu = fits.HDUList([primary, ebounds, matrix])

        try:
            hdu.writeto(rmffile, overwrite=overwrite)
        except IOError:
            message.error(
                "File {0} already exists. I will not overwrite it!".format(
                    rmffile))
            return 1

        return 0
예제 #10
0
    def correct_possible_shift(self):
        """See if there is a shift in the response array. When the spectral channels start at 0 in OGIP responses,
        then there is a possibility that the channel indexing needs to be shifted by 1. The SPEX format first
        channel should always be 1. Run this function after the conversion of OGIP to SPEX had taken place
        and if the first channel in the OGIP spectrum is 0.
        The ogip_to_spex method calls this function by default."""

        if not isinstance(self.spo, Spo) or not isinstance(self.res, Res):
            message.error(
                "Could not find spo and res information in this region.")
            return 1

        # This code is not prepared for situations where the channel arrays are swapped, like for RGS.
        # Luckily, the combination of a 0 first channel and a swapped array are rare.
        # In that case, we stop with a warning:
        if self.spo.swap:
            message.warning(
                "Not auto-detecting shifts in the response array. ")
            return 1

        # Start with the OGIP response object
        # Check the channel indices for the first group with useful data
        # Find such a group first in OGIP response:
        i = 0
        while True:
            # Find an energy bin with at least one response group.
            if self.resp.NumberGroups[i] == 0:
                i = i + 1
            else:
                break

        # Save the energy boundaries and calculate the average model energy for the group
        elow = self.resp.LowEnergy[i]
        ehigh = self.resp.HighEnergy[i]
        target_energy = (elow + ehigh) / 2.0

        # For this group, save the first channel of the group (F_CHAN)
        fchan = self.resp.FirstChannelGroup[0]

        # Find the array index for this channel number
        j = 0
        while True:
            if self.resp.Channel[j] != fchan:
                j = j + 1
            else:
                break

        # For this array index, the corresponding channel energy boundaries should be:
        lchan = self.resp.ChannelLowEnergy[j]
        hchan = self.resp.ChannelHighEnergy[j]
        target_channel = (lchan + hchan) / 2.0

        # Now find the same group and channel in the SPEX format objects
        # Find the group in the res object for the same model energy bin (with target_energy):
        s = 0
        while True:
            if self.res.eg1[s] < target_energy and self.res.eg2[
                    s] > target_energy:
                break
            else:
                s = s + 1

        # Find the target channel number in spo file
        t = 0
        while True:
            if self.spo.echan1[t] < target_channel and self.spo.echan2[
                    t] > target_channel:
                break
            else:
                t = t + 1

        # Corresponding first channel of this group according to SPEX format
        ic1 = self.res.ic1[s]

        # Calculate the difference between the SPEX channel number and the OGIP one.
        shift = int(t + 1 - ic1)

        if shift != 0:
            message.warning("Shift in response array detected.")
            message.proc_start(
                "Trying to shift indices with {0} ".format(shift))
            stat = self.res.channel_shift(shift)
            message.proc_end(stat)
        else:
            print("No shift in response array detected. Continuing...")

        return 0
예제 #11
0
def pha_to_spo(src, rmf, back=None, corr=None, save_grouping=False):
    """Convert the source (src) and optional background (back) and correction spectra (corr) from OGIP to SPEX format. 
    Please also provide an OGIP rmf object from the Rmf class to this function to read the channel energies.
    When the save_grouping flag is true, the grouping information in the PHA file will be copied to the spo file.
    The default behaviour is to ignore the grouping.
    This method returns a pyspextools Spo object containing the source and background rates."""

    if not isinstance(src, Pha):
        message.error("Input source spectrum is not a PHA object.")
        return 1

    if not isinstance(rmf, Rmf):
        message.error("Input response matrix is not an RMF object.")
        return 1

    if back is not None:
        if not isinstance(back, Pha):
            message.error("Input background spectrum is not a PHA object.")
            return 1
        input_back = True
    else:
        input_back = False

    if corr is not None:
        if not isinstance(corr, Pha):
            message.error("Input correction spectrum is not a PHA object.")
            return 1
        input_corr = True
    else:
        input_corr = False

    spo = Spo()

    # Determine number of channels and add to spo
    spo.nchan = np.append(spo.nchan, src.DetChans)
    spo.sponame = None
    spo.nregion = 1

    # Create zero arrays of length nchan to fill in the loop later
    spo.zero_spo(src.DetChans)

    # Loop through all the bins for the relevant spectral arrays
    for i in np.arange(src.DetChans):
        spo.tints[i] = src.Exposure * src.AreaScaling[i]

        # Calculate the source rates and errors
        if spo.tints[i] > 0:
            spo.ochan[i] = src.Rate[i] / src.AreaScaling[i]
            spo.dochan[i] = (src.StatError[i])**2 / src.AreaScaling[
                i]  # Add the errors in square later
        else:
            spo.ochan[i] = 0.
            spo.dochan[i] = 0.

        # Subtract background if available
        if input_back:
            btints = back.Exposure * back.AreaScaling[i]
            # Calculate backscale ratio
            if back.BackScaling[i] > 0:
                fb = src.BackScaling[i] / back.BackScaling[i]
            else:
                fb = 0.

            # Subtract background and calculate errors
            spo.ochan[
                i] = spo.ochan[i] - back.Rate[i] * fb / back.AreaScaling[i]
            spo.dochan[i] = spo.dochan[i] + (back.StatError[i] * fb /
                                             back.AreaScaling[i])**2
            spo.mbchan[i] = back.Rate[i] * fb / back.AreaScaling[i]
            spo.dbchan[i] = (back.StatError[i] * fb / back.AreaScaling[i])**2

            # Calculate the Exp_Rate backscale ratio
            if fb > 0 and src.Exposure > 0:
                spo.brat[i] = btints / spo.tints[i] / fb
            else:
                spo.brat[i] = 0.

        # Subtract correction spectrum, if available
        if input_corr:
            ctints = corr.Exposure * corr.AreaScaling[i]
            # Note: The influence of brat on the corr spectrum is not taken into account!
            if corr.BackScaling[i] > 0:
                fc = src.BackScaling[i] / corr.BackScaling[i]
            else:
                fc = 0.

            # Subtract correction spectrum and calculate errors
            spo.ochan[i] = spo.ochan[i] - corr.Rate[i] * fc / ctints
            spo.dochan[i] = spo.dochan[i] + corr.Rate[i] * (fc / ctints)**2
            spo.mbchan[i] = spo.mbchan[i] + corr.Rate[i] * fc / ctints
            spo.dbchan[i] = spo.dbchan[i] + corr.Rate[i] * (fc / ctints)**2

        # Set background to zero for zero exposure bins
        if spo.tints[i] <= 0.:
            spo.mbchan[i] = 0.
            spo.dbchan[i] = 0.
            spo.brat[i] = 0.

        spo.dochan[i] = math.sqrt(spo.dochan[i])
        spo.dbchan[i] = math.sqrt(spo.dbchan[i])

        # Set first, last and used variables
        if src.Quality[i] != 0:
            spo.used[i] = False

        if input_back:
            if back.Quality[i] != 0:
                spo.used[i] = False

        if input_corr:
            if corr.Quality[i] != 0:
                spo.used[i] = False

        if save_grouping:
            if src.Grouping[i] == 1:
                spo.first[i] = True
                spo.last[i] = False
            if src.Grouping[i] == 0:
                spo.first[i] = False
                if i < src.DetChans - 1:
                    if src.Grouping[i + 1] == 1:
                        spo.last[i] = True
                elif i == src.DetChans - 1:
                    spo.last[i] = True
                else:
                    spo.last[i] = False

        # Get channel boundaries from response
        # Channel boundary cannot be 0.
        if rmf.EnergyUnits != "keV":
            message.warning(
                "Energy units of keV are expected in the response file!")

        if rmf.ChannelLowEnergy[i] <= 0.:
            spo.echan1[i] = 1e-5
            message.warning(
                "Lowest channel boundary energy is 0. Set to 1E-5 to avoid problems."
            )
        else:
            spo.echan1[i] = rmf.ChannelLowEnergy[i]
        spo.echan2[i] = rmf.ChannelHighEnergy[i]

    # Check if channel order needs to be swapped
    if src.DetChans > 1:
        if spo.echan1[0] > spo.echan1[1]:
            spo.swap = True
            spo.swap_order()

    spo.empty = False

    return spo
예제 #12
0
def rmf_to_res(rmf, arf=None):
    """Convert an response matrix object from OGIP to SPEX format. The response matrix is translated one-to-one
    without optimizations. Providing an ARF object is optional. All groups in the OGIP matrix are put into one
    SPEX response component. This method returns a pyspextools Res object containing the response matrix."""

    if not isinstance(rmf, Rmf):
        message.error("The input RMF object is not of type Rmf.")
        return 1

    if arf is not None:
        if not isinstance(arf, Arf):
            message.error("The input ARF object is not of type Arf.")
        input_area = True
    else:
        input_area = False

    try:
        rmf.NumberChannels
    except NameError:
        message.error("The OGIP response matrix has not been initialised yet.")
        return 0

    res = Res()

    # Read the number of energy bins and groups
    res.nchan = np.append(res.nchan, rmf.NumberChannels)
    res.nsector = 1
    res.nregion = 1
    res.sector = np.append(res.sector, 1)
    res.region = np.append(res.region, 1)
    res.resname = None

    # Read the total number of groups (which is neg in SPEX format)
    res.neg = np.append(res.neg, rmf.NumberTotalGroups)

    res.eg1 = np.zeros(res.neg, dtype=float)
    res.eg2 = np.zeros(res.neg, dtype=float)
    res.nc = np.zeros(res.neg, dtype=int)
    res.ic1 = np.zeros(res.neg, dtype=int)
    res.ic2 = np.zeros(res.neg, dtype=int)

    # Read the total number of matrix elements
    nm = rmf.NumberTotalElements
    res.resp = np.zeros(nm, dtype=float)

    # Set the number of components to 1 (no optimization or re-ordering)
    res.ncomp = 1

    # Read the energy bin boundaries and group information
    g = 0  # Index for the group number
    m = 0  # Index for the matrix element number
    for i in np.arange(rmf.NumberEnergyBins):
        # Number of response groups for this energy bin
        ngrp = rmf.NumberGroups[i]
        for j in np.arange(ngrp):
            # Energy bin boundaries
            if rmf.LowEnergy[i] <= 0.:
                res.eg1[g] = 1e-7
                message.warning(
                    "Lowest energy boundary is 0. Set to 1E-7 to avoid problems."
                )
            else:
                res.eg1[g] = rmf.LowEnergy[i]

            res.eg2[g] = rmf.HighEnergy[i]
            if res.eg2[g] <= res.eg1[g]:
                message.error(
                    "Discontinous bins in energy array in channel {0}. Please check the numbers."
                    .format(i + 1))
                return

            res.nc[g] = rmf.NumberChannelsGroup[g]
            # Add the start channel to the IC to correct for cases where we start at channel 0/1
            res.ic1[g] = rmf.FirstChannelGroup[g]
            ic2 = res.ic1[g] + res.nc[g] - 1
            res.ic2[g] = ic2

            if input_area:
                area = arf.EffArea[i]
            else:
                area = 1.0

            for k in np.arange(res.nc[g]):
                res.resp[m] = rmf.Matrix[m] * area
                if res.resp[m] < 0.0:
                    res.resp[m] = 0.0
                m = m + 1
            g = g + 1

    if g > res.neg:
        message.error("Mismatch between number of groups.")
        return 0

    if m > nm:
        message.error("Mismatch between number of matrix elements.")
        return 0

    # Convert matrix to m**2 units for SPEX
    if input_area:
        if arf.ARFUnits == "cm2":
            res.resp *= 1.E-4
    else:
        res.resp *= 1.E-4

    # Check if channel order needs to be swapped
    if res.nchan > 1:
        if rmf.ChannelLowEnergy[0] > rmf.ChannelLowEnergy[1]:
            res.swap = True
            res.swap_order()

    res.empty = False

    return res