CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_2_SLHC2_patch1/src/Validation/HcalDigis/src/Digi2Raw2Digi.cc

Go to the documentation of this file.
00001 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00002 #include "FWCore/ServiceRegistry/interface/Service.h"
00003 #include "FWCore/PluginManager/interface/ModuleDef.h"
00004 #include "FWCore/Framework/interface/MakerMacros.h"
00005 
00006 #include "Validation/HcalDigis/interface/Digi2Raw2Digi.h"
00007 #include "DataFormats/HcalDigi/interface/HcalDigiCollections.h"
00008 #include "DataFormats/HcalDetId/interface/HcalSubdetector.h"
00009 #include "DataFormats/HcalDetId/interface/HcalElectronicsId.h"
00010 #include "DataFormats/HcalDetId/interface/HcalGenericDetId.h"
00011 #include "DataFormats/HcalDetId/interface/HcalZDCDetId.h"
00012 #include "DataFormats/HcalDetId/interface/HcalDetId.h"
00013 #include "DataFormats/HcalDigi/interface/HcalQIESample.h"
00014 
00015 #include "DQMServices/Core/interface/DQMStore.h"
00016 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00017 
00018 #include "DataFormats/HcalDigi/interface/HBHEDataFrame.h"
00019 #include "DataFormats/HcalDigi/interface/HFDataFrame.h"
00020 #include "DataFormats/HcalDigi/interface/HODataFrame.h"
00021 #include "DataFormats/HcalDigi/interface/ZDCDataFrame.h"
00022 #include "DataFormats/HcalDigi/interface/CastorDataFrame.h"
00023 #include <vector>
00024 #include <utility>
00025 #include <ostream>
00026 #include <string>
00027 #include <algorithm>
00028 #include <cmath>
00029 
00030 
00031 // Qunatities of interest in :
00032 // DataFormats/HcalDigi/test/HcalDigiDump.cc - example of Digi dumping...
00033 // DataFormats/HcalDigi/interface/HcalQIESample.h - adc, capid etc.
00034 // DataFormats/HcalDigi/interface/HBHEDataFrame.h -
00035 // zsMarkAndPass, zsUnsuppressed etc.
00036 
00037 template<class Digi>
00038 
00039 
00040 void Digi2Raw2Digi::compare(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
00041 
00042   typename edm::Handle<edm::SortedCollection<Digi> > digiCollection1;
00043   typename edm::SortedCollection<Digi>::const_iterator digiItr1;
00044   typename edm::Handle<edm::SortedCollection<Digi> > digiCollection2;
00045   typename edm::SortedCollection<Digi>::const_iterator digiItr2;
00046   
00047   if(unsuppressed) {  // ZDC
00048      iEvent.getByLabel ("simHcalUnsuppressedDigis", digiCollection1); 
00049   }
00050   else iEvent.getByLabel (inputTag1_, digiCollection1);
00051   
00052   iEvent.getByLabel (inputTag2_, digiCollection2);
00053   
00054   int size1 = 0;
00055   int size2 = 0;
00056   
00057   for (digiItr1=digiCollection1->begin();digiItr1!=digiCollection1->end();digiItr1++) {
00058     size1++;
00059   }
00060 
00061   for (digiItr2=digiCollection2->begin();digiItr2!=digiCollection2->end();digiItr2++) {
00062     size2++;
00063   }
00064 
00065   //std::cout << "Digi collections   size1 = "<< size1 
00066   //    << "   size2 = " << size2 << std::endl;
00067   
00068 
00069   // CYCLE over first DIGI collection ======================================
00070   
00071   for (digiItr1=digiCollection1->begin();digiItr1!=digiCollection1->end();digiItr1++) {
00072     HcalGenericDetId HcalGenDetId(digiItr1->id());
00073     int tsize =  (*digiItr1).size();
00074     int match = 0;
00075 
00076     if(HcalGenDetId.isHcalZDCDetId()){
00077       //for zdc
00078       HcalZDCDetId element(digiItr1->id());
00079       int zside =  element.zside();
00080       int section = element.section();
00081       int channel = element.channel();
00082       int gsub = HcalGenDetId.genericSubdet();
00083  
00084         if(section==3){// lumi section not reconstructed
00085           size2++;
00086           match = 1;
00087           goto lumi;
00088         }
00089       
00090 
00091         //std::cout<< " Zdc genSubdet="<< gsub << " zside=" <<zside
00092         //       << " section= "<< section << " channel " <<channel
00093         //      <<std::endl; 
00094 
00095       for (digiItr2=digiCollection2->begin();digiItr2!=digiCollection2->end();digiItr2++) {
00096         HcalZDCDetId element2(digiItr2->id());
00097         
00098         //int zside2 =  element2.zside();
00099         //int section2 = element2.section();
00100         //int channel2 = element2.channel();
00101         //int gsub2 = HcalGenDetId.genericSubdet();
00102 
00103         //std::cout<< " Zdc genSubdet="<<gsub2 
00104         //       << " zside=" <<zside2
00105         //       << " section= "<<section2
00106         //       <<" channel "<<channel2 
00107         //       <<std::endl; 
00108 
00109         if(element == element2) {
00110           match = 1;
00111           int identical = 1; 
00112           for (int i=0; i<tsize; i++) {
00113             double adc  =  (*digiItr1)[i].adc();     
00114             int capid   =  (*digiItr1)[i].capid();
00115             //    std::cout << std::endl << "     capid1=" << capid 
00116             //                << " adc1=" << adc 
00117             //              << std::endl;
00118             double adc2 =  (*digiItr2)[i].adc();     
00119             int capid2  =  (*digiItr2)[i].capid();
00120             //    std::cout << "     capid2=" << capid2 << " adc2=" << adc2 
00121             //              << std::endl;
00122             if( capid != capid2 || adc !=  adc2) {
00123               std::cout << "===> PROBLEM !!!  gebsubdet=" << gsub 
00124                         << " zside=" <<zside
00125                         << " section= "<< section << " channel " <<channel
00126                         << std::endl;
00127               std::cout << "     capid1["<< i << "]=" << capid 
00128                         << " adc1["<< i << "]=" << adc 
00129                         << "     capid2["<< i << "]=" << capid2 
00130                         << " adc2["<< i << "]=" << adc2
00131                         << std::endl;   
00132               identical = 0;
00133               meStatus->Fill(1.); 
00134               break;
00135             }
00136           } // end of DataFrames array
00137           if(identical) meStatus->Fill(0.); 
00138           break; // matched HcalZDCID  is processed,
00139           //  go to next (primary collection) cell  
00140         } 
00141       } // end of cycle over 2d DIGI collection 
00142     lumi:
00143       if (!match) {
00144         meStatus->Fill(2.); 
00145         std::cout << "===> PROBLEM !!!  gsubdet=" << gsub
00146                   << " zside=" <<zside
00147                   << " section= "<< section << " channel " <<channel
00148                   << " HcalZDCId match is not found !!!" 
00149                   << std::endl;
00150       }
00151   
00152     }
00153     else{
00154       //for Hcal subdetectors
00155       HcalDetId cell(digiItr1->id()); 
00156       int depth = cell.depth();
00157       int iphi  = cell.iphi()-1;
00158       int ieta  = cell.ieta();
00159       int sub   = cell.subdet();
00160       //    if(ieta > 0) ieta--;
00161       //    std::cout << " Cell subdet=" << sub << "  ieta=" << ieta 
00162       //              << "  inphi=" << iphi << "  depth=" << depth << std::endl;
00163       
00164       // CYCLE over second DIGI collection ======================================
00165       for (digiItr2=digiCollection2->begin();digiItr2!=digiCollection2->end();digiItr2++) {
00166   
00167         HcalDetId cell2(digiItr2->id());
00168         
00169         if( cell == cell2) {
00170           match = 1;
00171           int identical = 1; 
00172           for (int i=0; i<tsize; i++) {
00173             double adc  =  (*digiItr1)[i].adc();     
00174             int capid   =  (*digiItr1)[i].capid();
00175             //    std::cout << std::endl << "     capid1=" << capid 
00176             //                << " adc1=" << adc 
00177             //              << std::endl;
00178             double adc2 =  (*digiItr2)[i].adc();     
00179             int capid2  =  (*digiItr2)[i].capid();
00180             //    std::cout << "     capid2=" << capid2 << " adc2=" << adc2 
00181             //              << std::endl;
00182             if( capid != capid2 || adc !=  adc2) {
00183               std::cout << "===> PROBLEM !!!  subdet=" << sub << "  ieta="
00184                         << ieta  << "  inphi=" << iphi << "  depth=" << depth
00185                         << std::endl;
00186               std::cout << "     capid1["<< i << "]=" << capid 
00187                         << " adc1["<< i << "]=" << adc 
00188                         << "     capid2["<< i << "]=" << capid2 
00189                         << " adc2["<< i << "]=" << adc2
00190                         << std::endl;   
00191             identical = 0;
00192             meStatus->Fill(1.); 
00193             break;
00194             }
00195           } // end of DataFrames array
00196           if(identical) meStatus->Fill(0.); 
00197           break; // matched HcalID  is processed,
00198           //  go to next (primary collection) cell  
00199         } 
00200       } // end of cycle over 2d DIGI collection 
00201       if (!match) {
00202         meStatus->Fill(2.); 
00203         std::cout << "===> PROBLEM !!!  subdet=" << sub << "  ieta="
00204                   << ieta  << "  inphi=" << iphi << "  depth=" << depth
00205                   << " HcalID match is not found !!!" 
00206                   << std::endl;
00207       }
00208     }
00209   } // end of cycle over 1st DIGI collection    
00210   
00211   if (size1 != size2) {
00212     meStatus->Fill(3.); 
00213     std::cout << "===> PROBLEM !!!  Different size of Digi collections : "
00214               << size1  << " and " << size2 
00215               << std::endl;
00216   }
00217 }
00218 
00219 
00220 Digi2Raw2Digi::Digi2Raw2Digi(const edm::ParameterSet& iConfig)
00221   : inputTag1_(iConfig.getParameter<edm::InputTag>("digiLabel1")),
00222     inputTag2_(iConfig.getParameter<edm::InputTag>("digiLabel2")),
00223     outputFile_(iConfig.getUntrackedParameter<std::string>("outputFile")),
00224     dbe_(0)
00225 { 
00226 
00227   // DQM ROOT output
00228   if ( outputFile_.size() != 0 ) {
00229     edm::LogInfo("OutputInfo")
00230       << " Hcal RecHit Task histograms will be saved to '" 
00231       << outputFile_.c_str() << "'";
00232   } else {
00233     edm::LogInfo("OutputInfo") 
00234       << " Hcal RecHit Task histograms will NOT be saved";
00235   }
00236 
00237   // DQM service initialization
00238   dbe_ = edm::Service<DQMStore>().operator->();
00239    
00240   if ( dbe_ ) {
00241     dbe_->setCurrentFolder("Digi2Raw2DigiV/Digi2Raw2DigiTask");
00242   }
00243 
00244   // const char * sub = hcalselector_.c_str();
00245   char histo[100];
00246   sprintf (histo, "Digi2Raw2Digi_status") ;
00247   // bins: 1)full match 2)ID match, not content 3) no match 
00248   // 4) number of events with diff number of Digis
00249   meStatus    = dbe_->book1D(histo, histo, 5, 0., 5.);
00250 
00251 }
00252 
00253 void Digi2Raw2Digi::beginJob() {}
00254 void Digi2Raw2Digi::endJob() {
00255   if ( outputFile_.size() != 0 && dbe_ ) dbe_->save(outputFile_);
00256 }
00257 
00258 Digi2Raw2Digi::~Digi2Raw2Digi() {}
00259 
00260 
00261 void Digi2Raw2Digi::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup)
00262 {
00263   unsuppressed = 0;
00264   
00265   //  std::cout << "=== HBHE ==================" << std::endl; 
00266   compare<HBHEDataFrame>(iEvent,iSetup);
00267 
00268   //  std::cout << "=== HO ====================" << std::endl; 
00269   compare<HODataFrame>(iEvent,iSetup);
00270 
00271   //  std::cout << "=== HF ====================" << std::endl; 
00272   compare<HFDataFrame>(iEvent,iSetup);  
00273 
00274 
00275   //  std::cout << "=== ZDC ===================" << std::endl; 
00276   unsuppressed = 1;
00277   compare<ZDCDataFrame>(iEvent,iSetup);
00278   
00279   
00280   //  std::cout << "=== CASTOR ================" << std::endl; 
00281   //  compare<CastorDataFrame>(iEvent,iSetup); 
00282   
00283 }
00284 
00285 
00286 DEFINE_FWK_MODULE (Digi2Raw2Digi) ;