CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_10_patch1/src/CalibCalorimetry/HcalStandardModules/src/HcalPedestalMCWidths.cc

Go to the documentation of this file.
00001 // Original Author:  Steven Won
00002 //         Created:  Fri May  2 15:34:43 CEST 2008
00003 // Written to replace the combination of HcalPedestalAnalyzer and HcalPedestalAnalysis 
00004 // This code runs 1000x faster and produces all outputs from a single run
00005 // (ADC, fC in .txt plus an .xml file)
00006 //
00007 #include <memory>
00008 #include "CalibCalorimetry/HcalStandardModules/interface/HcalPedestalMCWidths.h"
00009 
00010 HcalPedestalMCWidths::HcalPedestalMCWidths(const edm::ParameterSet& ps)
00011 {
00012    firsttime = true;
00013    histflag = ps.getUntrackedParameter<bool>("saveHists",true);
00014 }
00015 
00016 
00017 HcalPedestalMCWidths::~HcalPedestalMCWidths()
00018 {
00019    HcalCovarianceMatrices * rawCovItem = new HcalCovarianceMatrices();
00020    std::ofstream outfile(widthsfilename.c_str());
00021    std::vector<MCWidthsBunch>::iterator bunch_it;
00022    for(bunch_it=Bunches.begin(); bunch_it != Bunches.end(); bunch_it++)
00023    {
00024       if(!bunch_it->usedflag) continue;
00025 
00026       HcalCovarianceMatrix item(bunch_it->detid.rawId());
00027       for(int i = 0; i != 4; i++){
00028          for(int j = 0; j != 10; j++){
00029             for(int k = 0; k != 10; k++){            
00030                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]);
00031                if(!std::isnan(bunch_it->sig[i][j][k]))
00032                {
00033                   item.setValue(i,j,k,bunch_it->sig[i][j][k]);
00034                }else{
00035                   item.setValue(i,j,k,0.0);
00036                }
00037             }
00038          }
00039       }
00040 
00041       rawCovItem->addValues(item);
00042 
00043    if(histflag)
00044    {
00045      if(bunch_it->detid.subdet() == 1){
00046          for(int i = 0; i != 4; i++){
00047             for(int j = 0; j != 10; j++){
00048                for(int k = j; k != 10; k++)
00049                if(!std::isnan(bunch_it->sig[i][j][k])) HBMeans[j]->Fill(k-j,bunch_it->sig[i][j][k]);
00050             }
00051          }
00052       }
00053       if(bunch_it->detid.subdet() == 2){
00054          for(int i = 0; i != 4; i++){
00055             for(int j = 0; j != 10; j++){
00056                for(int k = j; k != 10; k++)
00057                if(!std::isnan(bunch_it->sig[i][j][k])) HEMeans[j]->Fill(k-j,bunch_it->sig[i][j][k]);
00058             }
00059          }
00060       }
00061       if(bunch_it->detid.subdet() == 3){
00062          for(int i = 0; i != 4; i++){
00063             for(int j = 0; j != 10; j++){
00064                for(int k = j; k != 10; k++)
00065                if(!std::isnan(bunch_it->sig[i][j][k])) HFMeans[j]->Fill(k-j,bunch_it->sig[i][j][k]);
00066             }
00067          }
00068       }
00069       if(bunch_it->detid.subdet() == 4){
00070          for(int i = 0; i != 4; i++){
00071             for(int j = 0; j != 10; j++){
00072                for(int k = j; k != 10; k++)
00073                if(!std::isnan(bunch_it->sig[i][j][k])) HOMeans[j]->Fill(k-j,bunch_it->sig[i][j][k]);
00074             }
00075          }
00076       }
00077 
00078          
00079     }
00080     }
00081     HcalDbASCIIIO::dumpObject (outfile, (*rawCovItem) );    
00082 
00083     if(histflag)
00084     {
00085     theFile->cd();
00086     std::cout << "Writing histograms..." << std::endl;
00087     for(int i = 0; i != 10; i++){
00088        HBMeans[i]->Write();
00089        HFMeans[i]->Write();
00090        HEMeans[i]->Write();
00091        HOMeans[i]->Write();
00092     }
00093     }
00094     theFile->Close();
00095     std::cout << "ROOT file closed.\n";
00096 }
00097 
00098 // ------------ method called to for each event  ------------
00099 void
00100 HcalPedestalMCWidths::analyze(const edm::Event& e, const edm::EventSetup& iSetup)
00101 {
00102    using namespace edm;
00103    using namespace std;
00104 
00105    edm::Handle<HBHEDigiCollection> hbhe;              e.getByType(hbhe);
00106    edm::Handle<HODigiCollection> ho;                  e.getByType(ho);
00107    edm::Handle<HFDigiCollection> hf;                  e.getByType(hf);
00108    edm::ESHandle<HcalDbService> conditions;
00109    iSetup.get<HcalDbRecord>().get(conditions);
00110 
00111    edm::ESHandle<HcalPedestals> refPeds;
00112    iSetup.get<HcalPedestalsRcd>().get(refPeds);
00113    const HcalPedestals* myRefPeds = refPeds.product();
00114 
00115    if(firsttime)
00116    {
00117       runnum = e.id().run();
00118       std::string runnum_string;
00119       std::stringstream tempstringout;
00120       tempstringout << runnum;
00121       runnum_string = tempstringout.str();
00122       widthsfilename = runnum_string + "-MCwidths.txt";
00123       std::string rootfilename = runnum_string + "-MCwidths.root";
00124       theFile = new TFile(rootfilename.c_str(), "RECREATE");
00125 
00126       for(int i = 0; i!= 10; i++)
00127       {
00128          tempstringout.str("");
00129          tempstringout << i;
00130          std::string histname3 = tempstringout.str();
00131          std::string histname2 = "Mean covariance TS ";
00132          std::string histname1 = "HB ";
00133          std::string histname = histname1 + histname2 + histname3;
00134          HBMeans[i] = new TProfile(histname.c_str(),histname1.c_str(), 10, -0.50, 9.5);
00135          histname1 = "HE ";
00136          histname = histname1 + histname2 + histname3;
00137          HEMeans[i] = new TProfile(histname.c_str(),histname1.c_str(), 10, -0.50, 9.5);
00138          histname1 = "HF ";
00139          histname = histname1 + histname2 + histname3;
00140          HFMeans[i] = new TProfile(histname.c_str(),histname1.c_str(), 10, -0.50, 9.5);
00141          histname1 = "HO ";
00142          histname = histname1 + histname2 + histname3;
00143          HOMeans[i] = new TProfile(histname.c_str(),histname1.c_str(), 10, -0.50, 9.5);
00144       }
00145 
00146       edm::ESHandle<HcalElectronicsMap> refEMap;
00147       iSetup.get<HcalElectronicsMapRcd>().get(refEMap);
00148       const HcalElectronicsMap* myRefEMap = refEMap.product();
00149       std::vector<HcalGenericDetId> listEMap = myRefEMap->allPrecisionId();
00150       for (std::vector<HcalGenericDetId>::const_iterator it = listEMap.begin(); it != listEMap.end(); it++)
00151       {     
00152          HcalGenericDetId mygenid(it->rawId());
00153          if(mygenid.isHcalDetId())
00154          {
00155             MCWidthsBunch a;
00156             HcalDetId chanid(mygenid.rawId());
00157             a.detid = chanid;
00158             a.usedflag = false;
00159             for(int i = 0; i != 4; i++)
00160             {
00161                for(int j = 0; j != 10; j++){
00162                for(int k = 0; k != 10; k++)
00163                {
00164                   a.sig[i][j][k] = 0;
00165                   a.prod[i][j][k] = 0;
00166                   a.num[i][j][k] = 0;
00167                }}
00168             }
00169             Bunches.push_back(a);
00170          }
00171       }
00172       firsttime = false;
00173    }
00174 
00175    std::vector<MCWidthsBunch>::iterator bunch_it;
00176 
00177    for(HBHEDigiCollection::const_iterator j = hbhe->begin(); j != hbhe->end(); j++)
00178    {
00179       const HBHEDataFrame digi = (const HBHEDataFrame)(*j);
00180       for(bunch_it = Bunches.begin(); bunch_it != Bunches.end(); bunch_it++)
00181          if(bunch_it->detid.rawId() == digi.id().rawId()) break;
00182       DetId searchid(digi.id().rawId());
00183       const HcalPedestal* peds = (myRefPeds->getValues( searchid ));
00184       bunch_it->usedflag = true;
00185       int firstcapid = digi.sample(0).capid();
00186       for(int ts = 0; ts != 10; ts++)
00187       {
00188          for(int j = ts; j != 10; j++){
00189          bunch_it->num[firstcapid][ts][j] += 1;
00190          bunch_it->prod[firstcapid][ts][j] += 
00191                                            ((digi.sample(j).adc()-peds->getValue(firstcapid))
00192                                            *(digi.sample(ts).adc()-peds->getValue(digi.sample(ts).capid())));
00193          }
00194       }
00195    }
00196 
00197    for(HODigiCollection::const_iterator j = ho->begin(); j != ho->end(); j++)
00198    {
00199       const HODataFrame digi = (const HODataFrame)(*j);
00200       for(bunch_it = Bunches.begin(); bunch_it != Bunches.end(); bunch_it++)
00201          if(bunch_it->detid.rawId() == digi.id().rawId()) break;
00202       DetId searchid(digi.id().rawId());
00203       const HcalPedestal* peds = (myRefPeds->getValues( searchid ));
00204       bunch_it->usedflag = true;
00205       int firstcapid = digi.sample(0).capid();
00206       for(int ts = 0; ts != 10; ts++)
00207       {
00208          for(int j = ts; j != 10; j++){
00209          bunch_it->num[firstcapid][ts][j] += 1;
00210          bunch_it->prod[firstcapid][ts][j] += 
00211                                            ((digi.sample(j).adc()-peds->getValue(firstcapid))
00212                                            *(digi.sample(ts).adc()-peds->getValue(digi.sample(ts).capid())));
00213          }
00214       }
00215    }
00216 
00217    for(HFDigiCollection::const_iterator j = hf->begin(); j != hf->end(); j++)
00218    {
00219       const HFDataFrame digi = (const HFDataFrame)(*j);
00220       for(bunch_it = Bunches.begin(); bunch_it != Bunches.end(); bunch_it++)
00221          if(bunch_it->detid.rawId() == digi.id().rawId()) break;
00222       DetId searchid(digi.id().rawId());
00223       const HcalPedestal* peds = (myRefPeds->getValues( searchid ));
00224       bunch_it->usedflag = true;
00225       int firstcapid = digi.sample(0).capid();
00226       for(int ts = 0; ts != 10; ts++)
00227       {
00228          for(int j = ts; j != 10; j++){
00229          bunch_it->num[firstcapid][ts][j] += 1;
00230          bunch_it->prod[firstcapid][ts][j] += 
00231                                            ((digi.sample(j).adc()-peds->getValue(firstcapid))
00232                                            *(digi.sample(ts).adc()-peds->getValue(digi.sample(ts).capid())));
00233          }
00234       }
00235 
00236    }
00237 //this is the last brace
00238 }
00239 
00240 //define this as a plug-in
00241 DEFINE_FWK_MODULE(HcalPedestalMCWidths);