CMS 3D CMS Logo

Public Member Functions | Static Public Attributes | Protected Member Functions | Private Member Functions | Private Attributes

L1TdeECAL Class Reference

#include <L1TdeECAL.h>

Inheritance diagram for L1TdeECAL:
edm::EDAnalyzer

List of all members.

Public Member Functions

int iEtaiPhiToSMid (int, int)
 L1TdeECAL (const edm::ParameterSet &)
int TCCidToSMid (int)
 ~L1TdeECAL ()

Static Public Attributes

static const int nSM = 36
static const int nTTEta = 17
static const int nTTPhi = 4

Protected Member Functions

virtual void analyze (const edm::Event &, const edm::EventSetup &)
virtual void beginJob (void)
virtual void endJob ()

Private Member Functions

int verbose ()

Private Attributes

DQMStoredbe
edm::InputTag DEsource_
MonitorElementEcalEtMapDiff
MonitorElementEcalFGMapDiff
std::vector< MonitorElement * > etmapData
std::vector< MonitorElement * > etmapDiff
std::vector< MonitorElement * > etmapEmul
bool hasRecord_
std::string histFile_
std::string histFolder_
bool monitorDaemon_
int verbose_

Detailed Description

Definition at line 29 of file L1TdeECAL.h.


Constructor & Destructor Documentation

L1TdeECAL::L1TdeECAL ( const edm::ParameterSet iConfig) [explicit]

Definition at line 12 of file L1TdeECAL.cc.

References gather_cfg::cout, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), NULL, cmsCodeRules::cppFunctionSkipper::operator, and validate_alignment_devdb10_cfg::verbose.

                                                   {

  verbose_ = iConfig.getUntrackedParameter<int>("VerboseFlag",0);

  if(verbose())
    std::cout << "L1TdeECAL::L1TdeECAL()...\n" << std::flush;
  
  DEsource_ = iConfig.getParameter<edm::InputTag>("DataEmulCompareSource");
  histFolder_ = iConfig.getUntrackedParameter<std::string>("HistFolder", "L1TEMU//ECALexpert/");
  
  dbe = NULL;
  if (iConfig.getUntrackedParameter<bool>("DQMStore", false)) { 
    dbe = edm::Service<DQMStore>().operator->();
    dbe->setVerbose(0);
  }
  
  histFile_ = iConfig.getUntrackedParameter<std::string>("HistFile", "");
  if(iConfig.getUntrackedParameter<bool> ("disableROOToutput", true))
    histFile_ = "";

  if (histFile_.size()!=0) {
    edm::LogInfo("OutputRootFile") 
      << "L1TEmulator ECAL specific histograms will be saved to " 
      << histFile_.c_str() 
      << std::endl;
  }

  if(dbe!=NULL)
    dbe->setCurrentFolder(histFolder_);
  
  hasRecord_=true;

  if(verbose())
    std::cout << "L1TdeECAL::L1TdeECAL()...done.\n" << std::flush;
}
L1TdeECAL::~L1TdeECAL ( )

Definition at line 48 of file L1TdeECAL.cc.

{}

Member Function Documentation

void L1TdeECAL::analyze ( const edm::Event iEvent,
const edm::EventSetup iSetup 
) [protected, virtual]

get the comparison results

get the de candidates

--- Fill histograms(me) ---

Alternatively:

Implements edm::EDAnalyzer.

Definition at line 125 of file L1TdeECAL.cc.

References abs, gather_cfg::cout, diffTreeTool::diff, EcalBarrel, dedefs::ECALtp, dedefs::ETP, edm::Event::getByLabel(), ecalpyutils::ism(), edm::HandleBase::isValid(), j, LogDebug, phi, runTheMatrix::raw, evf::utils::sid, EcalElectronicsMapping::TCCid(), and validate_alignment_devdb10_cfg::verbose.

                                                                        {
  
  if(!hasRecord_)
    return;

  if(verbose())
    std::cout << "L1TdeECAL::analyze()  start\n" << std::flush;

  edm::Handle<L1DataEmulRecord> deRecord;
  iEvent.getByLabel(DEsource_, deRecord);

  if (!deRecord.isValid()) {
    edm::LogInfo("DataNotFound") 
      << "Cannot find L1DataEmulRecord with label "
      << DEsource_.label() 
      << " Please verify that comparator was successfully executed."
      << " Emulator DQM for ECAL will be skipped!"
      << std::endl;
    hasRecord_=false;
    return;
  }

  bool isComp = deRecord->get_isComp(ETP);
  if(!isComp) {
    if(verbose()) 
      std::cout << "[L1TdeECAL] Ecal information not generated in de-record."
                << " Skiping event!\n";
    return;
  }

  int DEncand[2];
  for(int j=0; j<2; j++) 
    DEncand[j] = deRecord->getNCand(ETP,j);
  
  if(verbose()) 
    std::cout << "[L1TdeECAL] ncands" 
              << " data: " << DEncand[0]
              << " emul: " << DEncand[1]
              << std::endl;

  
  L1DEDigiCollection deColl;
  deColl = deRecord->getColl();

  if(verbose()) {
    std::cout << "[L1TdeECAL] digis: \n";
    for(L1DEDigiCollection::const_iterator it=deColl.begin(); it!=deColl.end(); it++)
      if(it->sid()==ETP)
        std::cout << "\t" << *it << std::endl;
  }
  


  EcalElectronicsMapping emap;

  // d|e candidate loop
  for(L1DEDigiCollection::const_iterator it=deColl.begin(); 
      it!=deColl.end(); it++) {
    
    int    sid = it->sid();
    int    cid = it->cid();
    
    if(sid!=ETP)
      continue;

    if(it->empty())
      continue;

    assert(cid==ECALtp);
    
    if(verbose()) 
      std::cout << "[L1TdeECAL] processing digi: \n\t"
                << *it << "\n" << std::flush;
    
    //get (global) tt coordinates
    int iphi = (int)it->x1();
    int ieta = (int)it->x2();
    
    //get sm id
    int ism = iEtaiPhiToSMid(ieta,iphi); //1..36

    //get local indices
    int zside = (ieta>0?1:-1);
    int iet = abs(ieta);
    int ipt = iphi;
    ipt = ipt + 2;
    if ( ipt > 72 ) ipt = ipt - 72;
    ipt = (ipt-1)%4 + 1;
    if ( zside > 0 ) ipt = 5 - ipt;


    EcalTrigTowerDetId idt(zside, EcalBarrel, abs(ieta), iphi);
    // ... from EcalElectronicsMapping (needs addt'l lib)
    //int itt = map.iTT  (idt); //1..68
    int itcc = emap.TCCid(idt); //37..54(eb+) 55..72(eb-)
    // need here converter tcc->sm id
    int smid = TCCidToSMid(itcc);
    // ... or as in EBTriggerTowerTask (needs addt'l lib)
    //int ismt = Numbers::iSM( idt );

    if(verbose())
      std::cout << "L1TdeECAL \t"
                << " smid:" << smid 
                << " ism:"  << ism 
                << " itcc:" << itcc 
                << " local phi:" << ipt << " eta:" << iet 
                << "\n" << std::flush
                << *it
                << "\n" << std::flush;
    if(ism!=smid)
      LogDebug("L1TdeECAL") << "consistency check failure\n\t"
                            << " smid:" << smid 
                            << " ism:"  << ism 
                            << " itcc:" << itcc 
                            << std::endl;
    
    float xiet = iet+0.5;
    float xipt = ipt+0.5;
    
    //get energy values
    float rankarr[2]; 
    it->rank(rankarr);
    // get FG values
    unsigned int raw[2] ;
    it->data(raw) ;
    int FG[2] = { (raw[0] & 0x1000000)!=0, (raw[1] & 0x1000000)!=0 } ;

    int type = it->type(); //3 data only, 4 emul only
    if(type!=4 && etmapData[ism-1])
      etmapData[ism-1]->Fill(xiet-1, xipt-1, rankarr[0]);
    if(type!=3 && etmapEmul[ism-1])
      etmapEmul[ism-1]->Fill(xiet-1, xipt-1, rankarr[1]);
    if(type<2 && etmapDiff[ism-1]) {
      float diff = fabs(rankarr[0]-rankarr[1]);
      etmapDiff[ism-1]->Fill(xiet-1, xipt-1, diff);
      float phi = iphi ;
      if (phi>70) phi -= 73 ;
      phi *= 5 ;
      if (phi>0) phi -= 5 ;
      EcalEtMapDiff->Fill(ieta, phi, diff) ;
      diff = fabs(FG[0]-FG[1]);
      EcalFGMapDiff->Fill(ieta, phi, diff) ;
    }
  }//close loop over dedigi-cands


  if(verbose())
    std::cout << "L1TdeECAL::analyze() end.\n" << std::flush;

}
void L1TdeECAL::beginJob ( void  ) [protected, virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 51 of file L1TdeECAL.cc.

References DQMStore::book3D(), DQMStore::bookProfile2D(), gather_cfg::cout, j, cmsCodeRules::cppFunctionSkipper::operator, DQMStore::rmdir(), DQMStore::setCurrentFolder(), tmp, and validate_alignment_devdb10_cfg::verbose.

                        {

  if(verbose())
    std::cout << "L1TdeECAL::beginJob()  start\n";

  DQMStore* dbe = 0;
  dbe = edm::Service<DQMStore>().operator->();
  if(dbe) {
    dbe->setCurrentFolder(histFolder_);
    dbe->rmdir(histFolder_);
  }

  if(dbe) {
    dbe->setCurrentFolder(histFolder_);

    etmapData.reserve(nSM);
    etmapEmul.reserve(nSM);
    etmapDiff.reserve(nSM);
    etmapData.resize( nSM, static_cast<MonitorElement*>(0) );
    etmapEmul.resize( nSM, static_cast<MonitorElement*>(0) );
    etmapDiff.resize( nSM, static_cast<MonitorElement*>(0) );
    
    
    std::string lbl("");
    char tmp[100];
    for(int j=0; j<nSM; j++) {
      lbl.clear();
      sprintf(tmp, "etmapDataSM%d", j+1);
      lbl+=tmp;
      etmapData[j] = dbe->book3D(lbl.c_str(),lbl.c_str(),
                                 nTTEta, 0, nTTEta,
                                 nTTPhi, 0, nTTPhi,
                                 256, 0, 256.);
      sprintf(tmp, "etmapEmulSM%d", j+1);
      lbl.clear(); lbl+=tmp;
      etmapEmul[j] = dbe->book3D(lbl.c_str(),lbl.c_str(),
                                 nTTEta, 0, nTTEta,
                                 nTTPhi, 0, nTTPhi,
                                 256, 0, 256.);
      sprintf(tmp, "etmapDiffSM%d", j+1);
      lbl.clear(); lbl+=tmp;
      etmapDiff[j] = dbe->book3D(lbl.c_str(),lbl.c_str(),
                                 nTTEta, 0, nTTEta,
                                 nTTPhi, 0, nTTPhi,
                                 256, 0, 256.);
    }   
    lbl= "EcalEtMapDiff" ;
    EcalEtMapDiff = dbe->bookProfile2D(lbl.c_str(),lbl.c_str(),
                                       35, -17.5, 17.5,
                                       72, -10., 350.,
                                       256, 0, 256.);
    lbl= "EcalFGMapDiff" ;
    EcalFGMapDiff = dbe->bookProfile2D(lbl.c_str(),lbl.c_str(),
                                       35, -17.5, 17.5,
                                       72, -10., 350.,
                                       2, 0, 2.);
  }
  
  if(verbose())
    std::cout << "L1TdeECAL::beginJob()  end.\n" << std::flush;
}
void L1TdeECAL::endJob ( void  ) [protected, virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 114 of file L1TdeECAL.cc.

References gather_cfg::cout, and validate_alignment_devdb10_cfg::verbose.

                  {
  if(verbose())
    std::cout << "L1TdeECAL::endJob()...\n" << std::flush;
  if(histFile_.size()!=0  && dbe) 
    dbe->save(histFile_);
  //tbd delete emap;
}
int L1TdeECAL::iEtaiPhiToSMid ( int  ieta,
int  iphi 
)

Definition at line 283 of file L1TdeECAL.cc.

                                                {
  // barrel
  int iz = (ieta<0)?-1:1;
  iphi += 2;
  if (iphi > 72) iphi -= 72;
  const int kEBTowersInPhi = 4; //EBDetId::kTowersInPhi
  int sm = ( iphi - 1 ) / kEBTowersInPhi;
  if ( iz < 0 ) 
    sm += 19;
  else
    sm += 1;
  return sm;
}
int L1TdeECAL::TCCidToSMid ( int  tccid)

Definition at line 299 of file L1TdeECAL.cc.

                                    {
  // barrel
  if      ( tccid>37-1 && tccid<54+1) return tccid-37+19;
  else if ( tccid>55-1 && tccid<72+1) return tccid-55+ 1;
  else return 999;
}
int L1TdeECAL::verbose ( ) [inline, private]

Definition at line 56 of file L1TdeECAL.h.

References verbose_.

{return verbose_;}

Member Data Documentation

Definition at line 65 of file L1TdeECAL.h.

Definition at line 51 of file L1TdeECAL.h.

Definition at line 72 of file L1TdeECAL.h.

Definition at line 73 of file L1TdeECAL.h.

std::vector<MonitorElement*> L1TdeECAL::etmapData [private]

Definition at line 69 of file L1TdeECAL.h.

std::vector<MonitorElement*> L1TdeECAL::etmapDiff [private]

Definition at line 71 of file L1TdeECAL.h.

std::vector<MonitorElement*> L1TdeECAL::etmapEmul [private]

Definition at line 70 of file L1TdeECAL.h.

bool L1TdeECAL::hasRecord_ [private]

Definition at line 52 of file L1TdeECAL.h.

std::string L1TdeECAL::histFile_ [private]

Definition at line 59 of file L1TdeECAL.h.

std::string L1TdeECAL::histFolder_ [private]

Definition at line 62 of file L1TdeECAL.h.

bool L1TdeECAL::monitorDaemon_ [private]

Definition at line 66 of file L1TdeECAL.h.

const int L1TdeECAL::nSM = 36 [static]

Definition at line 44 of file L1TdeECAL.h.

const int L1TdeECAL::nTTEta = 17 [static]

Definition at line 45 of file L1TdeECAL.h.

const int L1TdeECAL::nTTPhi = 4 [static]

Definition at line 46 of file L1TdeECAL.h.

int L1TdeECAL::verbose_ [private]

Definition at line 55 of file L1TdeECAL.h.

Referenced by verbose().