CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_2_SLHC4_patch1/src/CalibCalorimetry/EcalPedestalOffsets/src/testChannel.cc

Go to the documentation of this file.
00001 
00011 #include "DataFormats/EcalRawData/interface/EcalDCCHeaderBlock.h"
00012 #include "DataFormats/EcalRawData/interface/EcalRawDataCollections.h"
00013 
00014 #include "TFile.h"
00015 
00016 #include "CalibCalorimetry/EcalPedestalOffsets/interface/testChannel.h"
00017 
00019 testChannel::testChannel (const edm::ParameterSet& paramSet) :
00020   m_digiCollection (paramSet.getParameter<std::string> ("digiCollection")) ,
00021   m_digiProducer (paramSet.getParameter<std::string> ("digiProducer")) ,
00022   m_headerProducer (paramSet.getParameter<std::string> ("headerProducer")) ,
00023   m_xmlFile (paramSet.getParameter<std::string> ("xmlFile")) ,
00024   m_DACmin (paramSet.getParameter<int> ("DACmin")) ,
00025   m_DACmax (paramSet.getParameter<int> ("DACmax")) ,
00026   m_RMSmax (paramSet.getParameter<double> ("RMSmax")) ,
00027   m_bestPed (paramSet.getParameter<int> ("bestPed")) ,
00028   m_xtal (paramSet.getParameter<int> ("xtal")) ,
00029   m_pedVSDAC ("pedVSDAC","pedVSDAC",100,150,250,m_DACmax-m_DACmin,m_DACmin,m_DACmax) ,
00030   m_singlePedVSDAC_1 ("singlePedVSDAC_1","pedVSDAC (g1) for xtal "+m_xtal,100,150,250,m_DACmax-m_DACmin,m_DACmin,m_DACmax) ,
00031   m_singlePedVSDAC_2 ("singlePedVSDAC_2","pedVSDAC (g2) for xtal "+m_xtal,100,150,250,m_DACmax-m_DACmin,m_DACmin,m_DACmax) ,
00032   m_singlePedVSDAC_3 ("singlePedVSDAC_3","pedVSDAC (g3) for xtal "+m_xtal,100,150,250,m_DACmax-m_DACmin,m_DACmin,m_DACmax)
00033 {
00034   edm::LogInfo ("testChannel") << " reading "
00035                                << " m_DACmin: " << m_DACmin
00036                                << " m_DACmax: " << m_DACmax
00037                                << " m_RMSmax: " << m_RMSmax
00038                                << " m_bestPed: " << m_bestPed ;
00039 }
00040 
00041 
00043 testChannel::~testChannel ()
00044 {
00045 }
00046 
00047 
00049 void testChannel::beginJob ()
00050 {
00051    LogDebug ("testChannel") << "entering beginJob ..." ;
00052 }
00053 
00054 
00056 void testChannel::analyze (edm::Event const& event, 
00057                            edm::EventSetup const& eventSetup) 
00058 {
00059    LogDebug ("testChannel") << "entering analyze ..." ;
00060 
00061    // get the headers
00062    // (one header for each supermodule)
00063    edm::Handle<EcalRawDataCollection> DCCHeaders ;
00064    event.getByLabel (m_headerProducer, DCCHeaders) ;
00065    if(!DCCHeaders.isValid())
00066    {
00067      edm::LogError ("testChannel") << "Error! can't get the product " 
00068                                    << m_headerProducer.c_str () ;
00069    }
00070 
00071    std::map <int,int> DACvalues ;
00072    
00073    // loop over the headers
00074    for ( EcalRawDataCollection::const_iterator headerItr= DCCHeaders->begin () ;
00075          headerItr != DCCHeaders->end () ; 
00076              ++headerItr ) 
00077      {
00078        EcalDCCHeaderBlock::EcalDCCEventSettings settings = headerItr->getEventSettings () ;
00079        DACvalues[getHeaderSMId (headerItr->id ())] = settings.ped_offset ;
00080 //       std::cout << "DCCid: " << headerItr->id () << "" ;
00081 //       std::cout << "Ped offset DAC: " << settings.ped_offset << "" ;
00082      } 
00083 
00084    // get the digis
00085    // (one digi for each crystal)
00086    edm::Handle<EBDigiCollection> pDigis;
00087    event.getByLabel (m_digiProducer, pDigis) ;
00088    if(!pDigis.isValid())
00089    {
00090      edm::LogError ("testChannel") << "Error! can't get the product " 
00091                                    << m_digiCollection.c_str () ;
00092    }
00093    
00094    // loop over the digis
00095    for (EBDigiCollection::const_iterator itdigi = pDigis->begin () ; 
00096         itdigi != pDigis->end () ; 
00097         ++itdigi) 
00098     {    
00099        EBDataFrame df( *itdigi );
00100        int gainId = df.sample (0).gainId () ;
00101        int crystalId = EBDetId(itdigi->id ()).ic () ;
00102        int smId = EBDetId(itdigi->id ()).ism () ;
00103 
00104        edm::LogInfo ("testChannel") << "channel " << event.id ()  
00105                                     << "\tcry: " << crystalId 
00106                                     << "\tG: " << gainId 
00107                                     << "\tDAC: " << DACvalues[smId] ;
00108 
00109        // loop over the samples
00110        for (int iSample = 0; iSample < EBDataFrame::MAXSAMPLES ; ++iSample) 
00111          {
00112             edm::LogInfo ("testChannel") << "\t`-->" << df.sample (iSample).adc () ;
00113             m_pedVSDAC.Fill (df.sample (iSample).adc (),DACvalues[smId]) ;
00114             if (crystalId == m_xtal)
00115               { 
00116                 if (gainId == 1) m_singlePedVSDAC_1.Fill (df.sample (iSample).adc (),DACvalues[smId]) ;
00117                 if (gainId == 2) m_singlePedVSDAC_2.Fill (df.sample (iSample).adc (),DACvalues[smId]) ;
00118                 if (gainId == 3) m_singlePedVSDAC_3.Fill (df.sample (iSample).adc (),DACvalues[smId]) ;
00119               }
00120          } // loop over the samples
00121     } // loop over the digis
00122 
00123 }
00124 
00126 void testChannel::endJob () 
00127 {
00128   char ccout[80] ;
00129   sprintf (ccout,"out_%d.root",m_xtal) ;
00130   TFile out (ccout,"RECREATE") ;
00131   out.cd () ;
00132   m_pedVSDAC.Write () ;
00133   m_singlePedVSDAC_1.Write () ;
00134   m_singlePedVSDAC_2.Write () ;
00135   m_singlePedVSDAC_3.Write () ;
00136   TProfile * profilo1 = m_singlePedVSDAC_1.ProfileX () ;
00137   TProfile * profilo2 = m_singlePedVSDAC_2.ProfileX () ;
00138   TProfile * profilo3 = m_singlePedVSDAC_3.ProfileX () ;
00139   profilo1->Write ("singleProfile_1") ;
00140   profilo2->Write ("singleProfile_2") ;
00141   profilo3->Write ("singleProfile_3") ;
00142   out.Close () ;
00143 }
00144 
00145 /*
00146 void testChannel::writeDb (EcalCondDBInterface* econn, 
00147                            MonRunIOV* moniov) 
00148 {}
00149 */
00150 
00151 int testChannel::getHeaderSMId (const int headerId) 
00152 {
00153   //PG FIXME temporary solution
00154   //PG FIXME check it is consistent with the TB!
00155   return 1 ;
00156 }
00157 
00158 
00159 
00160 
00161 void testChannel::subscribe ()
00162 {}
00163 
00164 void testChannel::subscribeNew ()
00165 {}
00166 
00167 void testChannel::unsubscribe ()
00168 {}
00169 
00170