CMS 3D CMS Logo

Public Member Functions | Private Attributes

HcalPedestalMCWidths Class Reference

#include <HcalPedestalMCWidths.h>

Inheritance diagram for HcalPedestalMCWidths:
edm::EDAnalyzer edm::EDConsumerBase

List of all members.

Public Member Functions

void analyze (const edm::Event &event, const edm::EventSetup &eventSetup)
 HcalPedestalMCWidths (const edm::ParameterSet &ps)
virtual ~HcalPedestalMCWidths ()

Private Attributes

std::vector< MCWidthsBunchBunches
bool firsttime
edm::InputTag hbheDigiCollectionTag_
TProfile * HBMeans [10]
TProfile * HEMeans [10]
edm::InputTag hfDigiCollectionTag_
TProfile * HFMeans [10]
bool histflag
edm::InputTag hoDigiCollectionTag_
TProfile * HOMeans [10]
std::string pedsADCfilename
int runnum
TFile * theFile
const HcalTopologytheTopology
std::string widthsfilename

Detailed Description

Definition at line 69 of file HcalPedestalMCWidths.h.


Constructor & Destructor Documentation

HcalPedestalMCWidths::HcalPedestalMCWidths ( const edm::ParameterSet ps)

Definition at line 12 of file HcalPedestalMCWidths.cc.

References firsttime, edm::ParameterSet::getUntrackedParameter(), and histflag.

                                                                    :
   hbheDigiCollectionTag_(ps.getParameter<edm::InputTag>("hbheDigiCollectionTag")),
   hoDigiCollectionTag_(ps.getParameter<edm::InputTag>("hoDigiCollectionTag")),
   hfDigiCollectionTag_(ps.getParameter<edm::InputTag>("hfDigiCollectionTag")) {

   firsttime = true;
   histflag = ps.getUntrackedParameter<bool>("saveHists",true);
}
HcalPedestalMCWidths::~HcalPedestalMCWidths ( ) [virtual]

Definition at line 22 of file HcalPedestalMCWidths.cc.

References HcalCovarianceMatrices::addValues(), Bunches, gather_cfg::cout, CastorDbASCIIIO::dumpObject(), HcalObjRepresent::Fill(), HBMeans, HEMeans, HFMeans, histflag, HOMeans, i, edm::isNotFinite(), j, gen::k, EdgesToViz::outfile, HcalCovarianceMatrix::setValue(), theFile, theTopology, and widthsfilename.

{
   HcalCovarianceMatrices * rawCovItem = new HcalCovarianceMatrices(theTopology);
   std::ofstream outfile(widthsfilename.c_str());
   std::vector<MCWidthsBunch>::iterator bunch_it;
   for(bunch_it=Bunches.begin(); bunch_it != Bunches.end(); bunch_it++)
   {
      if(!bunch_it->usedflag) continue;

      HcalCovarianceMatrix item(bunch_it->detid.rawId());
      for(int i = 0; i != 4; i++){
         for(int j = 0; j != 10; j++){
            for(int k = 0; k != 10; k++){            
               bunch_it->sig[i][j][k] = (bunch_it->prod[i][j][k]/bunch_it->num[i][j][k]);//-(bunch_it->cap[i]*bunch_it->cap[(i+j)%4]);
               if(!edm::isNotFinite(bunch_it->sig[i][j][k]))
               {
                  item.setValue(i,j,k,bunch_it->sig[i][j][k]);
               }else{
                  item.setValue(i,j,k,0.0);
               }
            }
         }
      }

      rawCovItem->addValues(item);

   if(histflag)
   {
     if(bunch_it->detid.subdet() == 1){
         for(int i = 0; i != 4; i++){
            for(int j = 0; j != 10; j++){
               for(int k = j; k != 10; k++)
               if(!edm::isNotFinite(bunch_it->sig[i][j][k])) HBMeans[j]->Fill(k-j,bunch_it->sig[i][j][k]);
            }
         }
      }
      if(bunch_it->detid.subdet() == 2){
         for(int i = 0; i != 4; i++){
            for(int j = 0; j != 10; j++){
               for(int k = j; k != 10; k++)
               if(!edm::isNotFinite(bunch_it->sig[i][j][k])) HEMeans[j]->Fill(k-j,bunch_it->sig[i][j][k]);
            }
         }
      }
      if(bunch_it->detid.subdet() == 3){
         for(int i = 0; i != 4; i++){
            for(int j = 0; j != 10; j++){
               for(int k = j; k != 10; k++)
               if(!edm::isNotFinite(bunch_it->sig[i][j][k])) HFMeans[j]->Fill(k-j,bunch_it->sig[i][j][k]);
            }
         }
      }
      if(bunch_it->detid.subdet() == 4){
         for(int i = 0; i != 4; i++){
            for(int j = 0; j != 10; j++){
               for(int k = j; k != 10; k++)
               if(!edm::isNotFinite(bunch_it->sig[i][j][k])) HOMeans[j]->Fill(k-j,bunch_it->sig[i][j][k]);
            }
         }
      }

         
    }
    }
    HcalDbASCIIIO::dumpObject (outfile, (*rawCovItem) );    

    if(histflag)
    {
    theFile->cd();
    std::cout << "Writing histograms..." << std::endl;
    for(int i = 0; i != 10; i++){
       HBMeans[i]->Write();
       HFMeans[i]->Write();
       HEMeans[i]->Write();
       HOMeans[i]->Write();
    }
    }
    theFile->Close();
    std::cout << "ROOT file closed.\n";
}

Member Function Documentation

void HcalPedestalMCWidths::analyze ( const edm::Event event,
const edm::EventSetup eventSetup 
) [virtual]

Implements edm::EDAnalyzer.

Definition at line 105 of file HcalPedestalMCWidths.cc.

References a, HcalQIESample::adc(), HcalElectronicsMap::allPrecisionId(), Bunches, HcalQIESample::capid(), MCWidthsBunch::detid, firsttime, edm::EventSetup::get(), edm::Event::getByLabel(), hbheDigiCollectionTag_, HBMeans, HEMeans, hfDigiCollectionTag_, HFMeans, hoDigiCollectionTag_, HOMeans, i, HODataFrame::id(), HBHEDataFrame::id(), edm::EventBase::id(), HFDataFrame::id(), j, gen::k, MCWidthsBunch::num, MCWidthsBunch::prod, edm::ESHandle< T >::product(), DetId::rawId(), edm::EventID::run(), runnum, HFDataFrame::sample(), HBHEDataFrame::sample(), HODataFrame::sample(), MCWidthsBunch::sig, AlCaHLTBitMon_QueryRunRegistry::string, theFile, theTopology, MCWidthsBunch::usedflag, and widthsfilename.

{
   using namespace edm;
   using namespace std;

   edm::Handle<HBHEDigiCollection> hbhe;              e.getByLabel(hbheDigiCollectionTag_, hbhe);
   edm::Handle<HODigiCollection> ho;                  e.getByLabel(hoDigiCollectionTag_, ho);
   edm::Handle<HFDigiCollection> hf;                  e.getByLabel(hfDigiCollectionTag_, hf);
   edm::ESHandle<HcalDbService> conditions;
   iSetup.get<HcalDbRecord>().get(conditions);

   edm::ESHandle<HcalPedestals> refPeds;
   iSetup.get<HcalPedestalsRcd>().get(refPeds);
   const HcalPedestals* myRefPeds = refPeds.product();

   if(firsttime)
   {
      theTopology=new HcalTopology(*(myRefPeds->topo()));
      runnum = e.id().run();
      std::string runnum_string;
      std::stringstream tempstringout;
      tempstringout << runnum;
      runnum_string = tempstringout.str();
      widthsfilename = runnum_string + "-MCwidths.txt";
      std::string rootfilename = runnum_string + "-MCwidths.root";
      theFile = new TFile(rootfilename.c_str(), "RECREATE");

      for(int i = 0; i!= 10; i++)
      {
         tempstringout.str("");
         tempstringout << i;
         std::string histname3 = tempstringout.str();
         std::string histname2 = "Mean covariance TS ";
         std::string histname1 = "HB ";
         std::string histname = histname1 + histname2 + histname3;
         HBMeans[i] = new TProfile(histname.c_str(),histname1.c_str(), 10, -0.50, 9.5);
         histname1 = "HE ";
         histname = histname1 + histname2 + histname3;
         HEMeans[i] = new TProfile(histname.c_str(),histname1.c_str(), 10, -0.50, 9.5);
         histname1 = "HF ";
         histname = histname1 + histname2 + histname3;
         HFMeans[i] = new TProfile(histname.c_str(),histname1.c_str(), 10, -0.50, 9.5);
         histname1 = "HO ";
         histname = histname1 + histname2 + histname3;
         HOMeans[i] = new TProfile(histname.c_str(),histname1.c_str(), 10, -0.50, 9.5);
      }

      edm::ESHandle<HcalElectronicsMap> refEMap;
      iSetup.get<HcalElectronicsMapRcd>().get(refEMap);
      const HcalElectronicsMap* myRefEMap = refEMap.product();
      std::vector<HcalGenericDetId> listEMap = myRefEMap->allPrecisionId();
      for (std::vector<HcalGenericDetId>::const_iterator it = listEMap.begin(); it != listEMap.end(); it++)
      {     
         HcalGenericDetId mygenid(it->rawId());
         if(mygenid.isHcalDetId())
         {
            MCWidthsBunch a;
            HcalDetId chanid(mygenid.rawId());
            a.detid = chanid;
            a.usedflag = false;
            for(int i = 0; i != 4; i++)
            {
               for(int j = 0; j != 10; j++){
               for(int k = 0; k != 10; k++)
               {
                  a.sig[i][j][k] = 0;
                  a.prod[i][j][k] = 0;
                  a.num[i][j][k] = 0;
               }}
            }
            Bunches.push_back(a);
         }
      }
      firsttime = false;
   }

   std::vector<MCWidthsBunch>::iterator bunch_it;

   for(HBHEDigiCollection::const_iterator j = hbhe->begin(); j != hbhe->end(); j++)
   {
      const HBHEDataFrame digi = (const HBHEDataFrame)(*j);
      for(bunch_it = Bunches.begin(); bunch_it != Bunches.end(); bunch_it++)
         if(bunch_it->detid.rawId() == digi.id().rawId()) break;
      DetId searchid(digi.id().rawId());
      const HcalPedestal* peds = (myRefPeds->getValues( searchid ));
      bunch_it->usedflag = true;
      int firstcapid = digi.sample(0).capid();
      for(int ts = 0; ts != 10; ts++)
      {
         for(int j = ts; j != 10; j++){
         bunch_it->num[firstcapid][ts][j] += 1;
         bunch_it->prod[firstcapid][ts][j] += 
                                           ((digi.sample(j).adc()-peds->getValue(firstcapid))
                                           *(digi.sample(ts).adc()-peds->getValue(digi.sample(ts).capid())));
         }
      }
   }

   for(HODigiCollection::const_iterator j = ho->begin(); j != ho->end(); j++)
   {
      const HODataFrame digi = (const HODataFrame)(*j);
      for(bunch_it = Bunches.begin(); bunch_it != Bunches.end(); bunch_it++)
         if(bunch_it->detid.rawId() == digi.id().rawId()) break;
      DetId searchid(digi.id().rawId());
      const HcalPedestal* peds = (myRefPeds->getValues( searchid ));
      bunch_it->usedflag = true;
      int firstcapid = digi.sample(0).capid();
      for(int ts = 0; ts != 10; ts++)
      {
         for(int j = ts; j != 10; j++){
         bunch_it->num[firstcapid][ts][j] += 1;
         bunch_it->prod[firstcapid][ts][j] += 
                                           ((digi.sample(j).adc()-peds->getValue(firstcapid))
                                           *(digi.sample(ts).adc()-peds->getValue(digi.sample(ts).capid())));
         }
      }
   }

   for(HFDigiCollection::const_iterator j = hf->begin(); j != hf->end(); j++)
   {
      const HFDataFrame digi = (const HFDataFrame)(*j);
      for(bunch_it = Bunches.begin(); bunch_it != Bunches.end(); bunch_it++)
         if(bunch_it->detid.rawId() == digi.id().rawId()) break;
      DetId searchid(digi.id().rawId());
      const HcalPedestal* peds = (myRefPeds->getValues( searchid ));
      bunch_it->usedflag = true;
      int firstcapid = digi.sample(0).capid();
      for(int ts = 0; ts != 10; ts++)
      {
         for(int j = ts; j != 10; j++){
         bunch_it->num[firstcapid][ts][j] += 1;
         bunch_it->prod[firstcapid][ts][j] += 
                                           ((digi.sample(j).adc()-peds->getValue(firstcapid))
                                           *(digi.sample(ts).adc()-peds->getValue(digi.sample(ts).capid())));
         }
      }

   }
//this is the last brace
}

Member Data Documentation

Definition at line 81 of file HcalPedestalMCWidths.h.

Referenced by analyze(), and ~HcalPedestalMCWidths().

Definition at line 96 of file HcalPedestalMCWidths.h.

Referenced by analyze(), and HcalPedestalMCWidths().

Definition at line 99 of file HcalPedestalMCWidths.h.

Referenced by analyze().

TProfile* HcalPedestalMCWidths::HBMeans[10] [private]

Definition at line 88 of file HcalPedestalMCWidths.h.

Referenced by analyze(), and ~HcalPedestalMCWidths().

TProfile* HcalPedestalMCWidths::HEMeans[10] [private]

Definition at line 89 of file HcalPedestalMCWidths.h.

Referenced by analyze(), and ~HcalPedestalMCWidths().

Definition at line 101 of file HcalPedestalMCWidths.h.

Referenced by analyze().

TProfile* HcalPedestalMCWidths::HFMeans[10] [private]

Definition at line 90 of file HcalPedestalMCWidths.h.

Referenced by analyze(), and ~HcalPedestalMCWidths().

Definition at line 97 of file HcalPedestalMCWidths.h.

Referenced by HcalPedestalMCWidths(), and ~HcalPedestalMCWidths().

Definition at line 100 of file HcalPedestalMCWidths.h.

Referenced by analyze().

TProfile* HcalPedestalMCWidths::HOMeans[10] [private]

Definition at line 91 of file HcalPedestalMCWidths.h.

Referenced by analyze(), and ~HcalPedestalMCWidths().

Definition at line 83 of file HcalPedestalMCWidths.h.

Definition at line 86 of file HcalPedestalMCWidths.h.

Referenced by analyze().

Definition at line 95 of file HcalPedestalMCWidths.h.

Referenced by analyze(), and ~HcalPedestalMCWidths().

Definition at line 93 of file HcalPedestalMCWidths.h.

Referenced by analyze(), and ~HcalPedestalMCWidths().

std::string HcalPedestalMCWidths::widthsfilename [private]

Definition at line 84 of file HcalPedestalMCWidths.h.

Referenced by analyze(), and ~HcalPedestalMCWidths().