CMS 3D CMS Logo

Public Member Functions | Private Types | Private Attributes | Static Private Attributes

L1RCTProducer Class Reference

#include <L1RCTProducer.h>

Inheritance diagram for L1RCTProducer:
edm::EDProducer edm::ProducerBase edm::ProductRegistryHelper

List of all members.

Public Member Functions

virtual void beginLuminosityBlock (edm::LuminosityBlock &lumiSeg, const edm::EventSetup &context)
virtual void beginRun (edm::Run &r, const edm::EventSetup &c)
const std::vector< int > getFedVectorFromOmds (const edm::EventSetup &)
const std::vector< int > getFedVectorFromRunInfo (const edm::EventSetup &)
 L1RCTProducer (const edm::ParameterSet &ps)
void printFedVector (const std::vector< int >)
void printUpdatedFedMask ()
void printUpdatedFedMaskVerbose ()
virtual void produce (edm::Event &e, const edm::EventSetup &c)
void updateConfiguration (const edm::EventSetup &)
void updateFedVector (const edm::EventSetup &, bool getFromOmds, int)
virtual ~L1RCTProducer ()

Private Types

enum  crateSection {
  c_min, ebOddFed = c_min, ebEvenFed, eeFed,
  hbheFed, hfFed, c_max = hfFed
}

Private Attributes

std::vector< int > bunchCrossings
std::vector< edm::InputTagecalDigis
L1RCTChannelMaskfedUpdatedMask
bool getFedsFromOmds
std::vector< edm::InputTaghcalDigis
unsigned int queryDelayInLS
unsigned int queryIntervalInLS
L1RCTrct
L1RCTLookupTablesrctLookupTables
bool useEcal
bool useHcal

Static Private Attributes

static const int crateFED [18][5]
static const int maxBarrel = 17
static const int maxEndcap = 28
static const int maxHF = 32
static const int minBarrel = 1
static const int minEndcap = 17
static const int minHF = 29

Detailed Description

Definition at line 47 of file L1RCTProducer.h.


Member Enumeration Documentation

Enumerator:
c_min 
ebOddFed 
ebEvenFed 
eeFed 
hbheFed 
hfFed 
c_max 

Definition at line 84 of file L1RCTProducer.h.


Constructor & Destructor Documentation

L1RCTProducer::L1RCTProducer ( const edm::ParameterSet ps) [explicit]

Definition at line 38 of file L1RCTProducer.cc.

                                                        : 
  rctLookupTables(new L1RCTLookupTables),
  rct(new L1RCT(rctLookupTables)),
  useEcal(conf.getParameter<bool>("useEcal")),
  useHcal(conf.getParameter<bool>("useHcal")),
  ecalDigis(conf.getParameter<std::vector<edm::InputTag> >("ecalDigis")),
  hcalDigis(conf.getParameter<std::vector<edm::InputTag> >("hcalDigis")),
  bunchCrossings(conf.getParameter<std::vector<int> >("BunchCrossings")),
  getFedsFromOmds(conf.getParameter<bool>("getFedsFromOmds")),
  queryDelayInLS(conf.getParameter<unsigned int>("queryDelayInLS")),
  queryIntervalInLS(conf.getParameter<unsigned int>("queryIntervalInLS")),
  fedUpdatedMask(0)
{
  produces<L1CaloEmCollection>();
  produces<L1CaloRegionCollection>();

}
L1RCTProducer::~L1RCTProducer ( ) [virtual]

Definition at line 56 of file L1RCTProducer.cc.

References fedUpdatedMask, rct, and rctLookupTables.

{
  if(rct != 0) delete rct;
  if(rctLookupTables != 0) delete rctLookupTables;
  if(fedUpdatedMask != 0) delete fedUpdatedMask;
}

Member Function Documentation

void L1RCTProducer::beginLuminosityBlock ( edm::LuminosityBlock lumiSeg,
const edm::EventSetup context 
) [virtual]

Reimplemented from edm::EDProducer.

Definition at line 76 of file L1RCTProducer.cc.

References getFedsFromOmds, edm::LuminosityBlockBase::luminosityBlock(), queryDelayInLS, queryIntervalInLS, edm::LuminosityBlockBase::run(), convertSQLiteXML::runNumber, and updateFedVector().

{
  // check LS number every LS, if the checkOMDS flag is set AND it's the right LS, update the FED vector from OMDS
  // can pass the flag as the bool??  but only check LS number if flag is true anyhow
  if (getFedsFromOmds)
    {
      unsigned int nLumi = lumiSeg.luminosityBlock(); // doesn't even need the (unsigned int) cast because LuminosityBlockNumber_t is already an unsigned int
      // LS count starts at 1, want to be able to delay 0 LS's intuitively
      if ( ( (nLumi - 1) == queryDelayInLS) 
           || (queryIntervalInLS > 0 && nLumi % queryIntervalInLS == 0 ) ) // to guard against problems if online DQM crashes; every 100 LS is ~20-30 minutes, not too big a load, hopefully not too long between
        {
          int runNumber = lumiSeg.run();
          //      std::cout << "Lumi section for this FED vector update is " << nLumi << std::endl;
          updateFedVector(context,true,runNumber); // OMDS
        }
      else if (queryIntervalInLS <= 0)
        {
          // don't do interval checking... cout message??
        }
    }
} 
void L1RCTProducer::beginRun ( edm::Run r,
const edm::EventSetup c 
) [virtual]

Reimplemented from edm::EDProducer.

Definition at line 64 of file L1RCTProducer.cc.

References edm::RunBase::run(), convertSQLiteXML::runNumber, updateConfiguration(), and updateFedVector().

{
  //  std::cout << "getFedsFromOmds is " << getFedsFromOmds << std::endl;

  updateConfiguration(eventSetup);
  
  int runNumber = run.run();
  updateFedVector(eventSetup,false,runNumber); // RUNINFO ONLY at beginning of run

}
const std::vector< int > L1RCTProducer::getFedVectorFromOmds ( const edm::EventSetup eventSetup)

Definition at line 300 of file L1RCTProducer.cc.

References edm::EventSetup::get(), getFedVectorFromRunInfo(), edm::ESHandleBase::isValid(), RunInfo::m_fed_in, edm::ESHandle< T >::product(), and summarizeEdmComparisonLogfiles::summary.

Referenced by updateFedVector().

{

  //  std::cout << "Getting FED vector from my specific ES RunInfo object" << std::endl;

  // get FULL FED vector from RunInfo object specifically created to have OMDS fed vector
  edm::ESHandle<RunInfo> sum;
  eventSetup.get<RunInfoRcd>().get("OmdsFedVector",sum); // using label to get my specific instance of RunInfo
  if (sum.isValid())
    {
      const RunInfo* summary=sum.product();
      const std::vector<int> fedvector = summary->m_fed_in;

      return fedvector;
    }
  else
    {
      return getFedVectorFromRunInfo(eventSetup);
    }
  
}
const std::vector< int > L1RCTProducer::getFedVectorFromRunInfo ( const edm::EventSetup eventSetup)

Definition at line 287 of file L1RCTProducer.cc.

References edm::EventSetup::get(), RunInfo::m_fed_in, edm::ESHandle< T >::product(), and summarizeEdmComparisonLogfiles::summary.

Referenced by getFedVectorFromOmds(), and updateFedVector().

{
  //  std::cout << "Getting FED vector from standard RunInfo object" << std::endl;
  // get FULL FED vector from RUNINFO
  edm::ESHandle<RunInfo> sum;
  eventSetup.get<RunInfoRcd>().get(sum);
  const RunInfo* summary=sum.product();
  const std::vector<int> fedvector = summary->m_fed_in;

  return fedvector;
}
void L1RCTProducer::printFedVector ( const std::vector< int >  fedVector)

Definition at line 397 of file L1RCTProducer.cc.

References filterCSVwithJSON::copy, and gather_cfg::cout.

{
  std::cout << "Contents of given fedVector: ";
  std::copy(fedVector.begin(), fedVector.end(), std::ostream_iterator<int>(std::cout, ", "));
  std::cout << std::endl;
}
void L1RCTProducer::printUpdatedFedMask ( )

Definition at line 405 of file L1RCTProducer.cc.

References gather_cfg::cout, fedUpdatedMask, and L1RCTChannelMask::print().

{
  if (fedUpdatedMask != 0)
    {
      fedUpdatedMask->print(std::cout);
    }
  else
    {
      std::cout << "Trying to print contents of fedUpdatedMask, but it doesn't exist!" << std::endl;
    }
}
void L1RCTProducer::printUpdatedFedMaskVerbose ( )

Definition at line 418 of file L1RCTProducer.cc.

References gather_cfg::cout, L1RCTChannelMask::ecalMask, fedUpdatedMask, L1RCTChannelMask::hcalMask, L1RCTChannelMask::hfMask, i, j, and gen::k.

{
  if (fedUpdatedMask != 0)
    {
      // print contents of fedvector
      std::cout << "Contents of fedUpdatedMask: ";
//       std::copy(fedUpdatedMask.begin(), fedUpdatedMask.end(), std::ostream_iterator<int>(std::cout, ", "));
      std::cout << "--> ECAL mask: " << std::endl;
      for (int i = 0; i < 18; i++)
        {
          for (int j = 0; j < 2; j++)
            {
              for (int k = 0; k < 28; k++)
                {
                  std::cout << fedUpdatedMask->ecalMask[i][j][k] << ", ";
                }
            }
        }
      std::cout << "--> HCAL mask: " << std::endl;
      for (int i = 0; i < 18; i++)
        {
          for (int j = 0; j < 2; j++)
            {
              for (int k = 0; k < 28; k++)
                {
                  std::cout << fedUpdatedMask->hcalMask[i][j][k] << ", ";
                }
            }
        }
      std::cout << "--> HF mask: " << std::endl;
      for (int i = 0; i < 18; i++)
        {
          for (int j = 0; j < 2; j++)
            {
              for (int k = 0; k < 4; k++)
                {
                  std::cout << fedUpdatedMask->hfMask[i][j][k] << ", ";
                }
            }
        }

      std::cout << std::endl;
    }
  else
    {
      //print error message
      std::cout << "Trying to print contents of fedUpdatedMask, but it doesn't exist!" << std::endl;
    }
}
void L1RCTProducer::produce ( edm::Event e,
const edm::EventSetup c 
) [virtual]

Implements edm::EDProducer.

Definition at line 324 of file L1RCTProducer.cc.

References bunchCrossings, L1RCT::digiInput(), patCandidatesForDimuonsSequences_cff::ecal, ecalDigis, Exception, edm::Event::getByLabel(), L1RCT::getIsolatedEGObjects(), L1RCT::getNonisolatedEGObjects(), L1RCT::getRegions(), patCandidatesForDimuonsSequences_cff::hcal, hcalDigis, i, j, L1RCT::processEvent(), rct, compare_using_db::sample, useEcal, and useHcal.

{


  std::auto_ptr<L1CaloEmCollection> rctEmCands (new L1CaloEmCollection);
  std::auto_ptr<L1CaloRegionCollection> rctRegions (new L1CaloRegionCollection);


  if(!(ecalDigis.size()==hcalDigis.size()&&hcalDigis.size()==bunchCrossings.size()))
      throw cms::Exception("BadInput")
        << "From what I see the number of your your ECAL input digi collections.\n"
        <<"is different from the size of your HCAL digi input collections\n"
        <<"or the size of your BX factor collection" 
        <<"They must be the same to correspond to the same Bxs\n"
        << "It does not matter if one of them is empty\n"; 




  // loop through and process each bx
    for (unsigned short sample = 0; sample < bunchCrossings.size(); sample++)
      {
        edm::Handle<EcalTrigPrimDigiCollection> ecal;
        edm::Handle<HcalTrigPrimDigiCollection> hcal;

        EcalTrigPrimDigiCollection ecalIn;
        HcalTrigPrimDigiCollection hcalIn;


        if(useHcal&&event.getByLabel(hcalDigis[sample], hcal))
          hcalIn = *hcal;

        if(useEcal&&event.getByLabel(ecalDigis[sample],ecal))
          ecalIn = *ecal;

        rct->digiInput(ecalIn,hcalIn);
        rct->processEvent();

      // Stuff to create
        for (int j = 0; j<18; j++)
          {
            L1CaloEmCollection isolatedEGObjects = rct->getIsolatedEGObjects(j);
            L1CaloEmCollection nonisolatedEGObjects = rct->getNonisolatedEGObjects(j);
            for (int i = 0; i<4; i++) 
              {
                isolatedEGObjects.at(i).setBx(bunchCrossings[sample]);
                nonisolatedEGObjects.at(i).setBx(bunchCrossings[sample]);
                rctEmCands->push_back(isolatedEGObjects.at(i));
                rctEmCands->push_back(nonisolatedEGObjects.at(i));
              }
          }
      
      
        for (int i = 0; i < 18; i++)
          {
            std::vector<L1CaloRegion> regions = rct->getRegions(i);
            for (int j = 0; j < 22; j++)
              {
                regions.at(j).setBx(bunchCrossings[sample]);
                rctRegions->push_back(regions.at(j));
              }
          }

      }

  
  //putting stuff back into event
  event.put(rctEmCands);
  event.put(rctRegions);
  
}
void L1RCTProducer::updateConfiguration ( const edm::EventSetup eventSetup)

Definition at line 100 of file L1RCTProducer.cc.

References alignCSCRings::e, edm::EventSetup::get(), h, edm::ESHandle< T >::product(), alignCSCRings::r, rctLookupTables, alignCSCRings::s, L1RCTLookupTables::setEcalScale(), L1RCTLookupTables::setHcalScale(), L1RCTLookupTables::setL1CaloEtScale(), and L1RCTLookupTables::setRCTParameters().

Referenced by beginRun().

{
  // Refresh configuration information every event
  // Hopefully, this does not take too much time
  // There should be a call back function in future to
  // handle changes in configuration
  // parameters to configure RCT (thresholds, etc)
  edm::ESHandle<L1RCTParameters> rctParameters;
  eventSetup.get<L1RCTParametersRcd>().get(rctParameters);
  const L1RCTParameters* r = rctParameters.product();

  //SCALES

  // energy scale to convert eGamma output
  edm::ESHandle<L1CaloEtScale> emScale;
  eventSetup.get<L1EmEtScaleRcd>().get(emScale);
  const L1CaloEtScale* s = emScale.product();

 // get energy scale to convert input from ECAL
  edm::ESHandle<L1CaloEcalScale> ecalScale;
  eventSetup.get<L1CaloEcalScaleRcd>().get(ecalScale);
  const L1CaloEcalScale* e = ecalScale.product();
  
  // get energy scale to convert input from HCAL
  edm::ESHandle<L1CaloHcalScale> hcalScale;
  eventSetup.get<L1CaloHcalScaleRcd>().get(hcalScale);
  const L1CaloHcalScale* h = hcalScale.product();

  // set scales
  rctLookupTables->setEcalScale(e);
  rctLookupTables->setHcalScale(h);

  rctLookupTables->setRCTParameters(r);
  rctLookupTables->setL1CaloEtScale(s);
}
void L1RCTProducer::updateFedVector ( const edm::EventSetup eventSetup,
bool  getFromOmds,
int  runNumber 
)

Definition at line 137 of file L1RCTProducer.cc.

References c_max, c_min, crateFED, ebEvenFed, ebOddFed, patCandidatesForDimuonsSequences_cff::ecal, L1RCTChannelMask::ecalMask, eeFed, fedUpdatedMask, spr::find(), edm::EventSetup::get(), getFedVectorFromOmds(), getFedVectorFromRunInfo(), hbheFed, L1RCTChannelMask::hcalMask, hfFed, L1RCTChannelMask::hfMask, i, j, gen::k, maxBarrel, maxEndcap, maxHF, minBarrel, minEndcap, minHF, edm::ESHandle< T >::product(), rctLookupTables, L1RCTLookupTables::setChannelMask(), and L1RCTLookupTables::setNoisyChannelMask().

Referenced by beginLuminosityBlock(), and beginRun().

{
  // list of RCT channels to mask
  edm::ESHandle<L1RCTChannelMask> channelMask;
  eventSetup.get<L1RCTChannelMaskRcd>().get(channelMask);
  const L1RCTChannelMask* cEs = channelMask.product();


  // list of Noisy RCT channels to mask
  edm::ESHandle<L1RCTNoisyChannelMask> hotChannelMask;
  eventSetup.get<L1RCTNoisyChannelMaskRcd>().get(hotChannelMask);
  const L1RCTNoisyChannelMask* cEsNoise = hotChannelMask.product();
  rctLookupTables->setNoisyChannelMask(cEsNoise);


  
  //Update the channel mask according to the FED VECTOR
  //This is the beginning of run. We delete the old
  //create the new and set it in the LUTs

  if(fedUpdatedMask!=0) delete fedUpdatedMask;

  fedUpdatedMask = new L1RCTChannelMask();
  // copy a constant object
  for (int i = 0; i < 18; i++)
    {
      for (int j = 0; j < 2; j++)
        {
          for (int k = 0; k < 28; k++)
            {
              fedUpdatedMask->ecalMask[i][j][k] = cEs->ecalMask[i][j][k];
              fedUpdatedMask->hcalMask[i][j][k] = cEs->hcalMask[i][j][k] ;
            }
          for (int k = 0; k < 4; k++)
            {
              fedUpdatedMask->hfMask[i][j][k] = cEs->hfMask[i][j][k];
            }
        }
    }


//   // adding fed mask into channel mask
  
  const std::vector<int> Feds = getFromOmds ? getFedVectorFromOmds(eventSetup) : getFedVectorFromRunInfo(eventSetup); // so can create/initialize/assign const quantity in one line accounting for if statement
  // wikipedia says this is exactly what it's for: http://en.wikipedia.org/wiki/%3F:#C.2B.2B

//   std::cout << "Contents of ";
//   std::cout << (getFromOmds ? "OMDS RunInfo" : "standard RunInfo");
//   std::cout << " FED vector" << std::endl;
//   printFedVector(Feds);

  std::vector<int> caloFeds;  // pare down the feds to the interesting ones
  // is this unneccesary?
  // Mike B : This will decrease the find speed so better do it
  for(std::vector<int>::const_iterator cf = Feds.begin(); cf != Feds.end(); ++cf)
    {
      int fedNum = *cf;
      if(fedNum > 600 && fedNum <724) 
        caloFeds.push_back(fedNum);
    }

  for(int  cr = 0; cr < 18; ++cr)
    {
      
      for(crateSection cs = c_min; cs <= c_max; cs = crateSection(cs +1)) 
        {
          bool fedFound = false;
          
          
          //Try to find the FED
          std::vector<int>::iterator fv = std::find(caloFeds.begin(),caloFeds.end(),crateFED[cr][cs]);
          if(fv!=caloFeds.end())
            fedFound = true;
          
          if(!fedFound) {
            int eta_min=0;
            int eta_max=0;
            bool phi_even[2] = {false};//, phi_odd = false;
            bool ecal=false;
            
            switch (cs) {
            case ebEvenFed :
              eta_min = minBarrel;
              eta_max = maxBarrel;
              phi_even[0] = true;
              ecal = true;      
              break;
              
            case ebOddFed:
              eta_min = minBarrel;
              eta_max = maxBarrel;
              phi_even[1] = true;
              ecal = true;      
              break;
              
            case eeFed:
              eta_min = minEndcap;
              eta_max = maxEndcap;
              phi_even[0] = true;
              phi_even[1] = true;
              ecal = true;      
              break;    
              
            case hbheFed:
              eta_min = minBarrel;
              eta_max = maxEndcap;
              phi_even[0] = true;
              phi_even[1] = true;
              ecal = false;
              break;
              
            case hfFed: 
              eta_min = minHF;
              eta_max = maxHF;
              
              phi_even[0] = true;
              phi_even[1] = true;
              ecal = false;
              break;
            default:
              break;
              
            }
            for(int ieta = eta_min; ieta <= eta_max; ++ieta)
              {
                if(ieta<=28) // barrel and endcap
                  for(int even = 0; even<=1 ; even++)
                    {    
                      if(phi_even[even])
                        {
                          if(ecal)
                            fedUpdatedMask->ecalMask[cr][even][ieta-1] = true;
                          else
                            fedUpdatedMask->hcalMask[cr][even][ieta-1] = true;
                        }
                    }
                else
                  for(int even = 0; even<=1 ; even++)
                    if(phi_even[even])
                      fedUpdatedMask->hfMask[cr][even][ieta-29] = true;
                
              }
          }
        }
    }
  
  rctLookupTables->setChannelMask(fedUpdatedMask); 

}

Member Data Documentation

std::vector<int> L1RCTProducer::bunchCrossings [private]

Definition at line 75 of file L1RCTProducer.h.

Referenced by produce().

const int L1RCTProducer::crateFED [static, private]
Initial value:
      {{613, 614, 603, 702, 718},
    {611, 612, 602, 700, 718},
    {627, 610, 601, 716, 722},
    {625, 626, 609, 714, 722},
    {623, 624, 608, 712, 722},
    {621, 622, 607, 710, 720},
    {619, 620, 606, 708, 720},
    {617, 618, 605, 706, 720},
    {615, 616, 604, 704, 718},
    {631, 632, 648, 703, 719},
    {629, 630, 647, 701, 719},
    {645, 628, 646, 717, 723},
    {643, 644, 654, 715, 723},
    {641, 642, 653, 713, 723},
    {639, 640, 652, 711, 721},
    {637, 638, 651, 709, 721},
    {635, 636, 650, 707, 721},
    {633, 634, 649, 705, 719}}

Definition at line 96 of file L1RCTProducer.h.

Referenced by updateFedVector().

std::vector<edm::InputTag> L1RCTProducer::ecalDigis [private]

Definition at line 73 of file L1RCTProducer.h.

Referenced by produce().

Definition at line 76 of file L1RCTProducer.h.

Referenced by beginLuminosityBlock().

std::vector<edm::InputTag> L1RCTProducer::hcalDigis [private]

Definition at line 74 of file L1RCTProducer.h.

Referenced by produce().

const int L1RCTProducer::maxBarrel = 17 [static, private]

Definition at line 98 of file L1RCTProducer.h.

Referenced by updateFedVector().

const int L1RCTProducer::maxEndcap = 28 [static, private]

Definition at line 100 of file L1RCTProducer.h.

Referenced by updateFedVector().

const int L1RCTProducer::maxHF = 32 [static, private]

Definition at line 102 of file L1RCTProducer.h.

Referenced by updateFedVector().

const int L1RCTProducer::minBarrel = 1 [static, private]

Definition at line 97 of file L1RCTProducer.h.

Referenced by updateFedVector().

const int L1RCTProducer::minEndcap = 17 [static, private]

Definition at line 99 of file L1RCTProducer.h.

Referenced by updateFedVector().

const int L1RCTProducer::minHF = 29 [static, private]

Definition at line 101 of file L1RCTProducer.h.

Referenced by updateFedVector().

unsigned int L1RCTProducer::queryDelayInLS [private]

Definition at line 77 of file L1RCTProducer.h.

Referenced by beginLuminosityBlock().

unsigned int L1RCTProducer::queryIntervalInLS [private]

Definition at line 78 of file L1RCTProducer.h.

Referenced by beginLuminosityBlock().

Definition at line 70 of file L1RCTProducer.h.

Referenced by produce(), and ~L1RCTProducer().

Definition at line 69 of file L1RCTProducer.h.

Referenced by updateConfiguration(), updateFedVector(), and ~L1RCTProducer().

bool L1RCTProducer::useEcal [private]

Definition at line 71 of file L1RCTProducer.h.

Referenced by produce().

bool L1RCTProducer::useHcal [private]

Definition at line 72 of file L1RCTProducer.h.

Referenced by produce().