CMS 3D CMS Logo

Public Types | Public Member Functions | Private Attributes

L1ExtraDQM::L1ExtraMonElement< CollectionType > Class Template Reference

#include <L1ExtraDQM.h>

List of all members.

Public Types

typedef
CollectionType::const_iterator 
CIterColl

Public Member Functions

void bookHistograms (const edm::EventSetup &evSetup, DQMStore *dbe, const std::string &l1ExtraObject, const std::vector< L1GtObject > &l1GtObj, const bool bookPhi=true, const bool bookEta=true)
void fillCharge (const CollectionType *collType, const bool validColl, const bool isL1Coll, const int bxInEvent)
 fill charge
void fillEtPhiEta (const CollectionType *collType, const bool validColl, const bool bookPhi, const bool bookEta, const bool isL1Coll, const int bxInEvent)
 ET, eta, phi.
void fillEtTotal (const CollectionType *collType, const bool validColl, const bool isL1Coll, const int bxInEvent)
 fill ET total in energy sums
void fillHfBitCounts (const CollectionType *collType, const bool validColl, const int countIndex, const bool isL1Coll, const int bxInEvent)
 fill bit counts in HFRings collections
void fillHfRingEtSums (const CollectionType *collType, const bool validColl, const int countIndex, const bool isL1Coll, const int bxInEvent)
 fill energy sums in HFRings collections
void fillNrObjects (const CollectionType *collType, const bool validColl, const bool isL1Coll, const int bxInEvent)
 number of objects
void fillPtPhiEta (const CollectionType *collType, const bool validColl, const bool bookPhi, const bool bookEta, const bool isL1Coll, const int bxInEvent)
 PT, eta, phi.
 L1ExtraMonElement (const edm::EventSetup &, const int)
virtual ~L1ExtraMonElement ()

Private Attributes

int m_indexCharge
int m_indexEt
int m_indexEta
int m_indexEtTotal
int m_indexHfBitCounts
int m_indexHfRingEtSums
int m_indexNrObjects
 histogram index for each quantity, set during histogram booking
int m_indexPhi
int m_indexPt
std::vector< MonitorElement * > m_monElement

Detailed Description

template<class CollectionType>
class L1ExtraDQM::L1ExtraMonElement< CollectionType >

Definition at line 87 of file L1ExtraDQM.h.


Member Typedef Documentation

template<class CollectionType >
typedef CollectionType::const_iterator L1ExtraDQM::L1ExtraMonElement< CollectionType >::CIterColl

Definition at line 97 of file L1ExtraDQM.h.


Constructor & Destructor Documentation

template<class CollectionType >
L1ExtraDQM::L1ExtraMonElement< CollectionType >::L1ExtraMonElement ( const edm::EventSetup evSetup,
const int  nrElements 
)
template<class CollectionType >
L1ExtraDQM::L1ExtraMonElement< CollectionType >::~L1ExtraMonElement ( ) [virtual]

Definition at line 260 of file L1ExtraDQM.h.

                                                                {

    //empty

}

Member Function Documentation

template<class CollectionType >
void L1ExtraDQM::L1ExtraMonElement< CollectionType >::bookHistograms ( const edm::EventSetup evSetup,
DQMStore dbe,
const std::string &  l1ExtraObject,
const std::vector< L1GtObject > &  l1GtObj,
const bool  bookPhi = true,
const bool  bookEta = true 
)

Definition at line 268 of file L1ExtraDQM.h.

References L1GetHistLimits::L1HistLimits::binThresholds, DQMStore::book1D(), CenJet, filterCSVwithJSON::copy, ForJet, HfBitCounts, HfRingEtSums, IsoEG, l1extra::L1HFRings::kNumRings, L1GetHistLimits::l1HistLimits(), LogDebug, LogTrace, L1GetHistLimits::L1HistLimits::lowerBinValue, Mu, NoIsoEG, L1GetHistLimits::L1HistLimits::nrBins, AlCaHLTBitMon_QueryRunRegistry::string, TauJet, and L1GetHistLimits::L1HistLimits::upperBinValue.

                            {

    // FIXME
    L1GtObject gtObj = l1GtObj.at(0);

    //
    std::string histName;
    std::string histTitle;
    std::string xAxisTitle;
    std::string yAxisTitle;

    std::string quantity = "";

    int indexHistogram = -1;

    if (gtObj == HfBitCounts) {

        L1GetHistLimits l1GetHistLimits(evSetup);
        const L1GetHistLimits::L1HistLimits& histLimits =
                l1GetHistLimits.l1HistLimits(gtObj, quantity);

        const int histNrBins = histLimits.nrBins;
        const double histMinValue = histLimits.lowerBinValue;
        const double histMaxValue = histLimits.upperBinValue;

        indexHistogram++;
        m_indexHfBitCounts = indexHistogram;

        for (int iCount = 0; iCount < l1extra::L1HFRings::kNumRings; ++iCount) {

            histName = l1ExtraObject + "_Count_" + boost::lexical_cast<
                    std::string>(iCount);
            histTitle = l1ExtraObject + ": count " + boost::lexical_cast<
                    std::string>(iCount);
            xAxisTitle = l1ExtraObject;
            yAxisTitle = "Entries";

            m_monElement.push_back(dbe->book1D(histName, histTitle, histNrBins,
                    histMinValue, histMaxValue));
            m_monElement[m_indexHfBitCounts + iCount]->setAxisTitle(xAxisTitle,
                    1);
            m_monElement[m_indexHfBitCounts + iCount]->setAxisTitle(yAxisTitle,
                    2);

        }

        return;

    }

    // number of objects per event
    if ((gtObj == Mu) || (gtObj == IsoEG) || (gtObj == NoIsoEG) || (gtObj
            == CenJet) || (gtObj == ForJet) || (gtObj == TauJet)) {

        quantity = "NrObjects";

        L1GetHistLimits l1GetHistLimits(evSetup);
        const L1GetHistLimits::L1HistLimits& histLimits =
                l1GetHistLimits.l1HistLimits(gtObj, quantity);

        const int histNrBins = histLimits.nrBins;
        const double histMinValue = histLimits.lowerBinValue;
        const double histMaxValue = histLimits.upperBinValue;

        histName = l1ExtraObject + "_NrObjectsPerEvent";
        histTitle = l1ExtraObject + ": number of objects per event";
        xAxisTitle = "Nr_" + l1ExtraObject;
        yAxisTitle = "Entries";

        m_monElement.push_back(dbe->book1D(histName, histTitle, histNrBins,
                histMinValue, histMaxValue));
        indexHistogram++;

        m_monElement[indexHistogram]->setAxisTitle(xAxisTitle, 1);
        m_monElement[indexHistogram]->setAxisTitle(yAxisTitle, 2);
        m_indexNrObjects = indexHistogram;

    }

    // transverse momentum (energy)  PT (ET) [GeV]


    quantity = "ET";
    std::string quantityLongName = " transverse energy ";

    if (gtObj == Mu) {
        quantity = "PT";
        quantityLongName = " transverse momentum ";
    }

    L1GetHistLimits l1GetHistLimits(evSetup);
    const L1GetHistLimits::L1HistLimits& histLimits =
            l1GetHistLimits.l1HistLimits(gtObj, quantity);

    const int histNrBinsET = histLimits.nrBins;
    const double histMinValueET = histLimits.lowerBinValue;
    const double histMaxValueET = histLimits.upperBinValue;
    const std::vector<float>& binThresholdsET = histLimits.binThresholds;

    float* binThresholdsETf;
    size_t sizeBinThresholdsET = binThresholdsET.size();
    binThresholdsETf = new float[sizeBinThresholdsET];
    copy(binThresholdsET.begin(), binThresholdsET.end(), binThresholdsETf);

    LogDebug("L1ExtraDQM") << "\n PT/ET histogram for " << l1ExtraObject
            << "\n histNrBinsET = " << histNrBinsET << "\n histMinValueET = "
            << histMinValueET << "\n histMaxValueET = " << histMaxValueET
            << "\n Last bin value represents the upper limit of the histogram"
            << std::endl;
    for (size_t iBin = 0; iBin < sizeBinThresholdsET; ++iBin) {
        LogTrace("L1ExtraDQM") << "Bin " << iBin << ": " << quantity << " = "
                << binThresholdsETf[iBin] << " GeV" << std::endl;

    }

    histName = l1ExtraObject + "_" + quantity;
    histTitle = l1ExtraObject + ":" + quantityLongName + quantity + " [GeV]";
    xAxisTitle = l1ExtraObject + "_" + quantity + " [GeV]";
    yAxisTitle = "Entries";

    if (gtObj == HfRingEtSums) {

        indexHistogram++;
        m_indexHfRingEtSums = indexHistogram;

        for (int iCount = 0; iCount < l1extra::L1HFRings::kNumRings; ++iCount) {

            histName = l1ExtraObject + "_Count_" + boost::lexical_cast<
                    std::string>(iCount);
            histTitle = l1ExtraObject + ": count " + boost::lexical_cast<
                    std::string>(iCount);
            xAxisTitle = l1ExtraObject;
            yAxisTitle = "Entries";

            m_monElement.push_back(dbe->book1D(histName, histTitle,
                    histNrBinsET, binThresholdsETf));

            m_monElement[m_indexHfRingEtSums + iCount]->setAxisTitle(xAxisTitle,
                    1);
            m_monElement[m_indexHfRingEtSums + iCount]->setAxisTitle(yAxisTitle,
                    2);

        }

    } else {

        m_monElement.push_back(dbe->book1D(histName, histTitle, histNrBinsET,
                binThresholdsETf));
        indexHistogram++;

        m_monElement[indexHistogram]->setAxisTitle(xAxisTitle, 1);
        m_monElement[indexHistogram]->setAxisTitle(yAxisTitle, 2);
        m_indexPt = indexHistogram;
        m_indexEt = indexHistogram;
        m_indexEtTotal = indexHistogram;
    }


    delete[] binThresholdsETf;

    //

    if (bookPhi) {

        quantity = "phi";

        // get limits and binning from L1Extra
        L1GetHistLimits l1GetHistLimits(evSetup);
        const L1GetHistLimits::L1HistLimits& histLimits =
                l1GetHistLimits.l1HistLimits(gtObj, quantity);

        const int histNrBinsPhi = histLimits.nrBins;
        const double histMinValuePhi = histLimits.lowerBinValue;
        const double histMaxValuePhi = histLimits.upperBinValue;
        const std::vector<float>& binThresholdsPhi = histLimits.binThresholds;

        float* binThresholdsPhif;
        size_t sizeBinThresholdsPhi = binThresholdsPhi.size();
        binThresholdsPhif = new float[sizeBinThresholdsPhi];
        copy(binThresholdsPhi.begin(), binThresholdsPhi.end(),
                binThresholdsPhif);

        LogDebug("L1ExtraDQM") << "\n phi histogram for " << l1ExtraObject
                << "\n histNrBinsPhi = " << histNrBinsPhi
                << "\n histMinValuePhi = " << histMinValuePhi
                << "\n histMaxValuePhi = " << histMaxValuePhi
                << "\n Last bin value represents the upper limit of the histogram"
                << std::endl;
        for (size_t iBin = 0; iBin < sizeBinThresholdsPhi; ++iBin) {
            LogTrace("L1ExtraDQM") << "Bin " << iBin << ": phi = "
                    << binThresholdsPhif[iBin] << " deg" << std::endl;

        }

        histName = l1ExtraObject + "_phi";
        histTitle = l1ExtraObject + ": phi distribution ";
        xAxisTitle = l1ExtraObject + "_phi [deg]";
        yAxisTitle = "Entries";

        m_monElement.push_back(dbe->book1D(histName, histTitle, histNrBinsPhi,
                histMinValuePhi, histMaxValuePhi));
        indexHistogram++;

        m_monElement[indexHistogram]->setAxisTitle(xAxisTitle, 1);
        m_monElement[indexHistogram]->setAxisTitle(yAxisTitle, 2);
        m_indexPhi = indexHistogram;

        delete[] binThresholdsPhif;
    }

    //


    if (bookEta) {

        quantity = "eta";

        // get limits and binning from L1Extra
        L1GetHistLimits l1GetHistLimits(evSetup);
        const L1GetHistLimits::L1HistLimits& histLimits =
                l1GetHistLimits.l1HistLimits(gtObj, quantity);

        const int histNrBinsEta = histLimits.nrBins;
        const double histMinValueEta = histLimits.lowerBinValue;
        const double histMaxValueEta = histLimits.upperBinValue;
        const std::vector<float>& binThresholdsEta = histLimits.binThresholds;

        //
        float* binThresholdsEtaf;
        size_t sizeBinThresholdsEta = binThresholdsEta.size();
        binThresholdsEtaf = new float[sizeBinThresholdsEta];
        copy(binThresholdsEta.begin(), binThresholdsEta.end(),
                binThresholdsEtaf);

        LogDebug("L1ExtraDQM") << "\n eta histogram for " << l1ExtraObject
                << "\n histNrBinsEta = " << histNrBinsEta
                << "\n histMinValueEta = " << histMinValueEta
                << "\n histMaxValueEta = " << histMaxValueEta
                << "\n Last bin value represents the upper limit of the histogram"
                << std::endl;
        for (size_t iBin = 0; iBin < sizeBinThresholdsEta; ++iBin) {
            LogTrace("L1ExtraDQM") << "Bin " << iBin << ": eta = "
                    << binThresholdsEtaf[iBin] << std::endl;

        }

        histName = l1ExtraObject + "_eta";
        histTitle = l1ExtraObject + ": eta distribution ";
        xAxisTitle = l1ExtraObject + "_eta";
        yAxisTitle = "Entries";

        m_monElement.push_back(dbe->book1D(histName, histTitle, histNrBinsEta,
                binThresholdsEtaf));
        indexHistogram++;

        m_monElement[indexHistogram]->setAxisTitle(xAxisTitle, 1);
        m_monElement[indexHistogram]->setAxisTitle(yAxisTitle, 2);
        m_indexEta = indexHistogram;

        delete[] binThresholdsEtaf;

    }

}
template<class CollectionType >
void L1ExtraDQM::L1ExtraMonElement< CollectionType >::fillCharge ( const CollectionType *  collType,
const bool  validColl,
const bool  isL1Coll,
const int  bxInEvent 
)

fill charge

Definition at line 633 of file L1ExtraDQM.h.

                                                                                                        {

    if (validColl) {
        for (CIterColl iterColl = collType->begin(); iterColl
                != collType->end(); ++iterColl) {

            if (isL1Coll && (iterColl->bx() != bxInEvent)) {
                continue;
            }

            m_monElement[m_indexCharge]->Fill(iterColl->charge());
        }
    }

}
template<class CollectionType >
void L1ExtraDQM::L1ExtraMonElement< CollectionType >::fillEtPhiEta ( const CollectionType *  collType,
const bool  validColl,
const bool  bookPhi,
const bool  bookEta,
const bool  isL1Coll,
const int  bxInEvent 
)

ET, eta, phi.

Definition at line 587 of file L1ExtraDQM.h.

References alignCSCRings::e, and rad2deg().

                                                                                          {

    if (validColl) {
        for (CIterColl iterColl = collType->begin(); iterColl
                != collType->end(); ++iterColl) {

            if (isL1Coll && (iterColl->bx() != bxInEvent)) {
                continue;
            }

            m_monElement[m_indexEt]->Fill(iterColl->et());

            if (bookPhi) {
                // add a very small quantity to get off the bin edge
                m_monElement[m_indexPhi]->Fill(rad2deg(iterColl->phi()) + 1.e-6);
            }

            if (bookEta) {
                m_monElement[m_indexEta]->Fill(iterColl->eta());
            }

        }
    }
}
template<class CollectionType >
void L1ExtraDQM::L1ExtraMonElement< CollectionType >::fillEtTotal ( const CollectionType *  collType,
const bool  validColl,
const bool  isL1Coll,
const int  bxInEvent 
)

fill ET total in energy sums

Definition at line 615 of file L1ExtraDQM.h.

                                                                                                        {

    if (validColl) {
        for (CIterColl iterColl = collType->begin(); iterColl
                != collType->end(); ++iterColl) {

            if (isL1Coll && (iterColl->bx() != bxInEvent)) {
                continue;
            }

            m_monElement[m_indexEtTotal]->Fill(iterColl->etTotal());
        }
    }

}
template<class CollectionType >
void L1ExtraDQM::L1ExtraMonElement< CollectionType >::fillHfBitCounts ( const CollectionType *  collType,
const bool  validColl,
const int  countIndex,
const bool  isL1Coll,
const int  bxInEvent 
)

fill bit counts in HFRings collections

Definition at line 651 of file L1ExtraDQM.h.

                                                                        {

    if (validColl) {
        for (CIterColl iterColl = collType->begin(); iterColl
                != collType->end(); ++iterColl) {

            if (isL1Coll && (iterColl->bx() != bxInEvent)) {
                continue;
            }

            m_monElement[m_indexHfBitCounts + countIndex]->Fill(
                    iterColl->hfBitCount(
                            (l1extra::L1HFRings::HFRingLabels) countIndex));
        }
    }

}
template<class CollectionType >
void L1ExtraDQM::L1ExtraMonElement< CollectionType >::fillHfRingEtSums ( const CollectionType *  collType,
const bool  validColl,
const int  countIndex,
const bool  isL1Coll,
const int  bxInEvent 
)

fill energy sums in HFRings collections

Definition at line 672 of file L1ExtraDQM.h.

                                                                        {

    if (validColl) {
        for (CIterColl iterColl = collType->begin(); iterColl
                != collType->end(); ++iterColl) {

            if (isL1Coll && (iterColl->bx() != bxInEvent)) {
                continue;
            }

            m_monElement[m_indexHfRingEtSums + countIndex]->Fill(
                    iterColl->hfEtSum(
                            (l1extra::L1HFRings::HFRingLabels) countIndex));
        }
    }

}
template<class CollectionType >
void L1ExtraDQM::L1ExtraMonElement< CollectionType >::fillNrObjects ( const CollectionType *  collType,
const bool  validColl,
const bool  isL1Coll,
const int  bxInEvent 
)

number of objects

Definition at line 538 of file L1ExtraDQM.h.

                                                  {

    if (validColl && isL1Coll) {
        size_t collSize = 0;
        for (CIterColl iterColl = collType->begin(); iterColl
                != collType->end(); ++iterColl) {

            if (iterColl->bx() == bxInEvent) {
                collSize++;
            }
        }
        m_monElement[m_indexNrObjects]->Fill(collSize);
    } else {
        size_t collSize = collType->size();
        m_monElement[m_indexNrObjects]->Fill(collSize);
    }
}
template<class CollectionType >
void L1ExtraDQM::L1ExtraMonElement< CollectionType >::fillPtPhiEta ( const CollectionType *  collType,
const bool  validColl,
const bool  bookPhi,
const bool  bookEta,
const bool  isL1Coll,
const int  bxInEvent 
)

PT, eta, phi.

Definition at line 559 of file L1ExtraDQM.h.

References alignCSCRings::e, and rad2deg().

                                                                                          {

    if (validColl) {
        for (CIterColl iterColl = collType->begin(); iterColl
                != collType->end(); ++iterColl) {

            if (isL1Coll && (iterColl->bx() != bxInEvent)) {
                continue;
            }

            m_monElement[m_indexPt]->Fill(iterColl->pt());

            if (bookPhi) {
                // add a very small quantity to get off the bin edge
                m_monElement[m_indexPhi]->Fill(rad2deg(iterColl->phi()) + 1.e-6);
            }

            if (bookEta) {
                m_monElement[m_indexEta]->Fill(iterColl->eta());
            }

        }
    }
}

Member Data Documentation

template<class CollectionType >
int L1ExtraDQM::L1ExtraMonElement< CollectionType >::m_indexCharge [private]

Definition at line 147 of file L1ExtraDQM.h.

template<class CollectionType >
int L1ExtraDQM::L1ExtraMonElement< CollectionType >::m_indexEt [private]

Definition at line 143 of file L1ExtraDQM.h.

template<class CollectionType >
int L1ExtraDQM::L1ExtraMonElement< CollectionType >::m_indexEta [private]

Definition at line 145 of file L1ExtraDQM.h.

template<class CollectionType >
int L1ExtraDQM::L1ExtraMonElement< CollectionType >::m_indexEtTotal [private]

Definition at line 146 of file L1ExtraDQM.h.

template<class CollectionType >
int L1ExtraDQM::L1ExtraMonElement< CollectionType >::m_indexHfBitCounts [private]

Definition at line 148 of file L1ExtraDQM.h.

template<class CollectionType >
int L1ExtraDQM::L1ExtraMonElement< CollectionType >::m_indexHfRingEtSums [private]

Definition at line 149 of file L1ExtraDQM.h.

template<class CollectionType >
int L1ExtraDQM::L1ExtraMonElement< CollectionType >::m_indexNrObjects [private]

histogram index for each quantity, set during histogram booking

Definition at line 141 of file L1ExtraDQM.h.

template<class CollectionType >
int L1ExtraDQM::L1ExtraMonElement< CollectionType >::m_indexPhi [private]

Definition at line 144 of file L1ExtraDQM.h.

template<class CollectionType >
int L1ExtraDQM::L1ExtraMonElement< CollectionType >::m_indexPt [private]

Definition at line 142 of file L1ExtraDQM.h.

template<class CollectionType >
std::vector<MonitorElement*> L1ExtraDQM::L1ExtraMonElement< CollectionType >::m_monElement [private]