CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_3/src/CalibCalorimetry/EcalLaserAnalyzer/plugins/EcalABAnalyzer.cc

Go to the documentation of this file.
00001 /* 
00002  *  \class EcalABAnalyzer
00003  *
00004  *  $Date: 2012/02/09 10:07:36 $
00005  *  primary author: Julie Malcles - CEA/Saclay
00006  *  author: Gautier Hamel De Monchenault - CEA/Saclay
00007  */
00008 
00009 #include <TAxis.h>
00010 #include <TH1.h>
00011 #include <TProfile.h>
00012 #include <TTree.h>
00013 #include <TChain.h>
00014 #include <TFile.h>
00015 #include <TMath.h>
00016 
00017 #include "CalibCalorimetry/EcalLaserAnalyzer/plugins/EcalABAnalyzer.h"
00018 
00019 #include <sstream>
00020 #include <iostream>
00021 #include <fstream>
00022 #include <iomanip>
00023 #include <vector>
00024 
00025 #include <FWCore/MessageLogger/interface/MessageLogger.h>
00026 #include <FWCore/Utilities/interface/Exception.h>
00027 
00028 #include <FWCore/Framework/interface/ESHandle.h>
00029 #include <FWCore/Framework/interface/EventSetup.h>
00030 
00031 #include <Geometry/EcalMapping/interface/EcalElectronicsMapping.h>
00032 #include <Geometry/EcalMapping/interface/EcalMappingRcd.h>
00033 
00034 #include <FWCore/Framework/interface/Event.h>
00035 #include <FWCore/Framework/interface/MakerMacros.h>
00036 #include <FWCore/ParameterSet/interface/ParameterSet.h>
00037 #include <DataFormats/EcalDigi/interface/EcalDigiCollections.h>
00038 #include <DataFormats/EcalDetId/interface/EcalElectronicsId.h>
00039 #include <DataFormats/EcalDetId/interface/EcalDetIdCollections.h>
00040 #include <DataFormats/EcalRawData/interface/EcalRawDataCollections.h>
00041 
00042 #include <CalibCalorimetry/EcalLaserAnalyzer/interface/TShapeAnalysis.h>
00043 #include <CalibCalorimetry/EcalLaserAnalyzer/interface/TMom.h>
00044 #include <CalibCalorimetry/EcalLaserAnalyzer/interface/TAPD.h>
00045 #include <CalibCalorimetry/EcalLaserAnalyzer/interface/TAPDPulse.h>
00046 #include <CalibCalorimetry/EcalLaserAnalyzer/interface/PulseFitWithFunction.h>
00047 #include <CalibCalorimetry/EcalLaserAnalyzer/interface/ME.h>
00048 #include <CalibCalorimetry/EcalLaserAnalyzer/interface/MEGeom.h>
00049 
00050 using namespace std;
00051 
00052 
00053 //========================================================================
00054 EcalABAnalyzer::EcalABAnalyzer(const edm::ParameterSet& iConfig)
00055 //========================================================================
00056   :
00057 iEvent(0),
00058 
00059 // framework parameters with default values
00060 _nsamples(        iConfig.getUntrackedParameter< unsigned int >( "nSamples",           10 ) ),
00061 _presample(       iConfig.getUntrackedParameter< unsigned int >( "nPresamples",         2 ) ),
00062 _firstsample(     iConfig.getUntrackedParameter< unsigned int >( "firstSample",         1 ) ),
00063 _lastsample(      iConfig.getUntrackedParameter< unsigned int >( "lastSample",          2 ) ),
00064 _timingcutlow(    iConfig.getUntrackedParameter< unsigned int >( "timingCutLow",        2 ) ),
00065 _timingcuthigh(   iConfig.getUntrackedParameter< unsigned int >( "timingCutHigh",       9 ) ),
00066 _timingquallow(   iConfig.getUntrackedParameter< unsigned int >( "timingQualLow",       3 ) ),
00067 _timingqualhigh(  iConfig.getUntrackedParameter< unsigned int >( "timingQualHigh",      8 ) ),
00068 _ratiomincutlow(  iConfig.getUntrackedParameter< double       >( "ratioMinCutLow",    0.4 ) ),
00069 _ratiomincuthigh( iConfig.getUntrackedParameter< double       >( "ratioMinCutHigh",  0.95 ) ),
00070 _ratiomaxcutlow(  iConfig.getUntrackedParameter< double       >( "ratioMaxCutLow",    0.8 ) ),
00071 _presamplecut(    iConfig.getUntrackedParameter< double       >( "presampleCut",      5.0 ) ),
00072 _niter(           iConfig.getUntrackedParameter< unsigned int >( "nIter",               3 ) ),
00073 _alpha(           iConfig.getUntrackedParameter< double       >( "alpha",       1.5076494 ) ),
00074 _beta(            iConfig.getUntrackedParameter< double       >( "beta",        1.5136036 ) ),
00075 _nevtmax(         iConfig.getUntrackedParameter< unsigned int >( "nEvtMax",           200 ) ),
00076 _noise(           iConfig.getUntrackedParameter< double       >( "noise",             2.0 ) ),
00077 _chi2cut(         iConfig.getUntrackedParameter< double       >( "chi2cut",         100.0 ) ), 
00078 _ecalPart(        iConfig.getUntrackedParameter< std::string  >( "ecalPart",         "EB" ) ),
00079 _qualpercent(     iConfig.getUntrackedParameter< double       >( "qualPercent",       0.2 ) ), 
00080 _debug(           iConfig.getUntrackedParameter< int          >( "debug",              0  ) ), 
00081 nCrys(                                                                               NCRYSEB),
00082 runType(-1), runNum(0), fedID(-1), dccID(-1), side(2), lightside(2), iZ(1),
00083 phi(-1), eta(-1), event(0), color(-1), channelIteratorEE(0)
00084 
00085   //========================================================================
00086 
00087 {
00088 
00089 
00090   // Initialization from cfg file
00091 
00092   resdir_                 = iConfig.getUntrackedParameter<std::string>("resDir");
00093 
00094   digiCollection_         = iConfig.getParameter<std::string>("digiCollection");
00095   digiProducer_           = iConfig.getParameter<std::string>("digiProducer");
00096 
00097   eventHeaderCollection_  = iConfig.getParameter<std::string>("eventHeaderCollection");
00098   eventHeaderProducer_    = iConfig.getParameter<std::string>("eventHeaderProducer");
00099 
00100 
00101 
00102   // Geometrical constants initialization 
00103 
00104 
00105   if (_ecalPart == "EB") {
00106     nCrys    = NCRYSEB;
00107   } else {
00108     nCrys    = NCRYSEE;
00109   }
00110   iZ         =  1;
00111   if(_fedid <= 609 ) 
00112     iZ       = -1;
00113   
00114   for(unsigned int j=0;j<nCrys;j++){
00115     iEta[j]=-1;
00116     iPhi[j]=-1;
00117     iTowerID[j]=-1;
00118     iChannelID[j]=-1;
00119     idccID[j]=-1;
00120     iside[j]=-1;
00121     wasTimingOK[j]=true;
00122     wasGainOK[j]=true; 
00123     nevtAB[j]=0 ;
00124   }
00125   
00126   // Quality check flags
00127   
00128   isGainOK=true;
00129   isTimingOK=true;
00130 
00131   // Objects dealing with pulses
00132 
00133   APDPulse = new TAPDPulse(_nsamples, _presample, _firstsample, _lastsample, 
00134                            _timingcutlow, _timingcuthigh, _timingquallow, _timingqualhigh,
00135                            _ratiomincutlow,_ratiomincuthigh, _ratiomaxcutlow);
00136 
00137 
00138   // Objects needed for npresample calculation
00139 
00140   Delta01=new TMom();
00141   Delta12=new TMom();
00142   _fitab=true;
00143 
00144 }
00145 
00146 //========================================================================
00147 EcalABAnalyzer::~EcalABAnalyzer(){
00148   //========================================================================
00149 
00150 
00151   // do anything here that needs to be done at desctruction time
00152   // (e.g. close files, deallocate resources etc.)
00153 
00154 }
00155 
00156 
00157 
00158 //========================================================================
00159 void EcalABAnalyzer::beginJob() {
00160   //========================================================================
00161 
00162   
00163   //Calculate alpha and beta
00164 
00165 
00166   // Define output results filenames and shape analyzer object (alpha,beta) 
00167   //=====================================================================
00168   
00169   // 1) AlphaBeta files
00170   
00171   doesABTreeExist=true;
00172   
00173   std::stringstream nameabinitfile;
00174   nameabinitfile << resdir_ <<"/ABInit.root";
00175   alphainitfile=nameabinitfile.str();
00176   
00177   std::stringstream nameabfile;  
00178   std::stringstream link;  
00179   nameabfile << resdir_ <<"/AB.root";
00180   
00181   FILE *test;
00182   test = fopen(nameabinitfile.str().c_str(),"r"); 
00183   if(test == NULL) {
00184     doesABTreeExist=false;
00185     _fitab=true;
00186   };
00187   delete test;
00188 
00189   
00190   TFile *fAB=0; TTree *ABInit=0;
00191   if(doesABTreeExist){
00192     fAB=new TFile(nameabinitfile.str().c_str());
00193   }
00194   if(doesABTreeExist && fAB){
00195     ABInit = (TTree*) fAB->Get("ABCol0");
00196   }
00197   
00198   // 2) Shape analyzer
00199   
00200   if(doesABTreeExist && fAB && ABInit && ABInit->GetEntries()!=0){
00201     shapana= new TShapeAnalysis(ABInit, _alpha, _beta, 5.5, 1.0);
00202     doesABTreeExist=true;
00203   }else{
00204     shapana= new TShapeAnalysis(_alpha, _beta, 5.5, 1.0);
00205     doesABTreeExist=false;
00206     _fitab=true;  
00207   }  
00208   shapana -> set_const(_nsamples,_firstsample,_lastsample,
00209                        _presample, _nevtmax, _noise, _chi2cut);
00210   
00211   if(doesABTreeExist && fAB ) fAB->Close();
00212   
00213   if(_fitab){
00214     alphafile=nameabfile.str();
00215   }else{
00216     alphafile=alphainitfile;
00217     link<< "ln -s "<<resdir_<<"/ABInit.root "<< resdir_<<"/AB.root";  
00218     system(link.str().c_str());
00219   }
00220 
00221   // Define output results files' names
00222   
00223   std::stringstream namefile;
00224   namefile << resdir_ <<"/AB.root";
00225   alphafile=namefile.str(); 
00226 
00227 }
00228 
00229 //========================================================================
00230 void EcalABAnalyzer:: analyze( const edm::Event & e, const  edm::EventSetup& c){
00231 //========================================================================
00232 
00233   ++iEvent;
00234   
00235   // retrieving DCC header
00236   edm::Handle<EcalRawDataCollection> pDCCHeader;
00237   const  EcalRawDataCollection* DCCHeader=0;
00238   try {
00239     e.getByLabel(eventHeaderProducer_,eventHeaderCollection_, pDCCHeader);
00240     DCCHeader=pDCCHeader.product();
00241   }catch ( std::exception& ex ) {
00242     std::cerr << "Error! can't get the product retrieving DCC header" << eventHeaderCollection_.c_str() <<" "<< eventHeaderProducer_.c_str() << std::endl;
00243   }
00244   
00245   //retrieving crystal data from Event
00246   edm::Handle<EBDigiCollection>  pEBDigi;
00247   const  EBDigiCollection* EBDigi=0;
00248   edm::Handle<EEDigiCollection>  pEEDigi;
00249   const  EEDigiCollection* EEDigi=0;
00250   if (_ecalPart == "EB") {
00251     try {
00252       e.getByLabel(digiProducer_,digiCollection_, pEBDigi); 
00253       EBDigi=pEBDigi.product(); 
00254     }catch ( std::exception& ex ) {
00255       std::cerr << "Error! can't get the product retrieving EB crystal data " << digiCollection_.c_str() << std::endl;
00256     } 
00257   } else if (_ecalPart == "EE") {
00258     try {
00259       e.getByLabel(digiProducer_,digiCollection_, pEEDigi); 
00260       EEDigi=pEEDigi.product(); 
00261     }catch ( std::exception& ex ) {
00262       std::cerr << "Error! can't get the product retrieving EE crystal data " << digiCollection_.c_str() << std::endl;
00263     } 
00264   } else {
00265     std::cout <<" Wrong ecalPart in cfg file " << std::endl;
00266     return;
00267   }
00268 
00269   // retrieving electronics mapping
00270   edm::ESHandle< EcalElectronicsMapping > ecalmapping;
00271   const EcalElectronicsMapping* TheMapping=0; 
00272   try{
00273     c.get< EcalMappingRcd >().get(ecalmapping);
00274     TheMapping = ecalmapping.product();
00275   }catch ( std::exception& ex ) {
00276     std::cerr << "Error! can't get the product EcalMappingRcd"<< std::endl;
00277   }
00278 
00279   // =============================
00280   // Decode DCCHeader Information 
00281   // =============================
00282   
00283 
00284   for ( EcalRawDataCollection::const_iterator headerItr= DCCHeader->begin();headerItr != DCCHeader->end(); 
00285         ++headerItr ) {
00286     
00287     // Get run type and run number 
00288 
00289     int fed = headerItr->fedId();  
00290     if(fed!=_fedid && _fedid!=-999) continue; 
00291     
00292     runType=headerItr->getRunType();
00293     runNum=headerItr->getRunNumber();
00294     event=headerItr->getLV1();
00295 
00296     dccID=headerItr->getDccInTCCCommand();
00297     fedID=headerItr->fedId();  
00298     lightside=headerItr->getRtHalf();
00299 
00300     // Check fed corresponds to the DCC in TCC
00301     
00302     if( 600+dccID != fedID ) continue;
00303 
00304     // Cut on runType
00305 
00306     if(runType!=EcalDCCHeaderBlock::LASER_STD 
00307        && runType!=EcalDCCHeaderBlock::LASER_GAP 
00308        && runType!=EcalDCCHeaderBlock::LASER_POWER_SCAN 
00309        && runType!=EcalDCCHeaderBlock::LASER_DELAY_SCAN ) return; 
00310     
00311     // Retrieve laser color and event number
00312     
00313     EcalDCCHeaderBlock::EcalDCCEventSettings settings = headerItr->getEventSettings();
00314     color = settings.wavelength;
00315     if( color <0 ) return;
00316 
00317     std::vector<int>::iterator iter = find( colors.begin(), colors.end(), color );
00318     if( iter==colors.end() ){
00319       colors.push_back( color );
00320     }
00321   }
00322   
00323 
00324   // Cut on fedID
00325 
00326   if(fedID!=_fedid && _fedid!=-999) return; 
00327   
00328 
00329   // ===========================
00330   // Decode EBDigis Information
00331   // ===========================
00332   
00333   int  adcGain=0;
00334   
00335   if (EBDigi){
00336     
00337     // Loop on crystals
00338     //=================== 
00339 
00340     for ( EBDigiCollection::const_iterator digiItr= EBDigi->begin(); digiItr != EBDigi->end(); ++digiItr ) {  // Loop on crystals
00341       
00342       // Retrieve geometry
00343       //=================== 
00344 
00345       EBDetId id_crystal(digiItr->id()) ;
00346       EBDataFrame df( *digiItr );
00347 
00348       int etaG = id_crystal.ieta() ;  // global
00349       int phiG = id_crystal.iphi() ;  // global
00350 
00351       int etaL ;  // local
00352       int phiL ;  // local
00353       std::pair<int, int> LocalCoord=MEEBGeom::localCoord( etaG , phiG );
00354       
00355       etaL=LocalCoord.first ;
00356       phiL=LocalCoord.second ;
00357       
00358       eta = etaG;
00359       phi = phiG;
00360 
00361       side=MEEBGeom::side(etaG,phiG);
00362    
00363       // Recover the TT id and the electronic crystal numbering from EcalElectronicsMapping
00364 
00365       EcalElectronicsId elecid_crystal = TheMapping->getElectronicsId(id_crystal);
00366 
00367       int towerID=elecid_crystal.towerId();
00368       int strip=elecid_crystal.stripId();
00369       int xtal=elecid_crystal.xtalId();
00370       int channelID=  5*(strip-1) + xtal-1; 
00371       
00372       unsigned int channel=MEEBGeom::electronic_channel( etaL, phiL );
00373 
00374       assert( channel < nCrys );
00375       
00376       iEta[channel]=eta;
00377       iPhi[channel]=phi;   
00378       iTowerID[channel]=towerID;
00379       iChannelID[channel]=channelID;
00380       idccID[channel]=dccID;
00381       iside[channel]=side;
00382 
00383       // APD Pulse
00384       //=========== 
00385       
00386       // Loop on adc samples 
00387       
00388       for (unsigned int i=0; i< (*digiItr).size() ; ++i ) {   
00389 
00390         EcalMGPASample samp_crystal(df.sample(i));
00391         adc[i]=samp_crystal.adc() ;    
00392         adcG[i]=samp_crystal.gainId();   
00393         adc[i]*=adcG[i];
00394         if (i==0) adcGain=adcG[i];
00395         if (i>0) adcGain=TMath::Max(adcG[i],adcGain);  
00396       }
00397 
00398       APDPulse->setPulse(adc);
00399       
00400 
00401       // Quality checks
00402       //================
00403       
00404       if(adcGain!=1) nEvtBadGain[channel]++;   
00405       if(!APDPulse->isTimingQualOK()) nEvtBadTiming[channel]++;
00406       nEvtTot[channel]++;
00407 
00408 
00409       // Fill if Pulse is fine
00410       //=======================
00411       
00412       if( APDPulse->isPulseOK() && lightside==side){
00413         
00414         Delta01->addEntry(APDPulse->getDelta(0,1));
00415         Delta12->addEntry(APDPulse->getDelta(1,2));
00416         
00417         if( nevtAB[channel] < _nevtmax && _fitab ){
00418           if(doesABTreeExist)  shapana -> putAllVals(channel, adc, eta, phi);   
00419           else shapana -> putAllVals(channel, adc, eta, phi, dccID, side, towerID, channelID);
00420           nevtAB[channel]++ ;
00421         }
00422       }
00423     }
00424     
00425   } else if (EEDigi) {
00426     
00427     // Loop on crystals
00428     //===================
00429 
00430     for ( EEDigiCollection::const_iterator digiItr= EEDigi->begin(); digiItr != EEDigi->end(); ++digiItr ) {  // Loop on crystals
00431       
00432       // Retrieve geometry
00433       //=================== 
00434 
00435       EEDetId id_crystal(digiItr->id()) ;
00436       EEDataFrame df( *digiItr );
00437       
00438       phi = id_crystal.ix() ; 
00439       eta = id_crystal.iy() ; 
00440 
00441       int iX = (phi-1)/5+1;
00442       int iY = (eta-1)/5+1;
00443   
00444       side=MEEEGeom::side( iX, iY ,iZ);       
00445       EcalElectronicsId elecid_crystal = TheMapping->getElectronicsId(id_crystal);
00446       
00447       int towerID=elecid_crystal.towerId();
00448       int channelID=elecid_crystal.channelId()-1;  
00449       
00450       int hashedIndex=100000*eta+phi;
00451       
00452       if( channelMapEE.count(hashedIndex) == 0 ){
00453         channelMapEE[hashedIndex]=channelIteratorEE;
00454         channelIteratorEE++;
00455       }
00456       
00457       unsigned int channel=channelMapEE[hashedIndex];
00458      
00459       assert ( channel < nCrys );
00460 
00461       iEta[channel]=eta;
00462       iPhi[channel]=phi;
00463       iTowerID[channel]=towerID;
00464       iChannelID[channel]=channelID;
00465       idccID[channel]=dccID;
00466       iside[channel]=side;
00467 
00468       // APD Pulse
00469       //=========== 
00470 
00471       if( (*digiItr).size()>10) std::cout <<"SAMPLES SIZE > 10!" <<  (*digiItr).size()<< std::endl;
00472  
00473       // Loop on adc samples  
00474 
00475       for (unsigned int i=0; i< (*digiItr).size() ; ++i ) { 
00476 
00477         EcalMGPASample samp_crystal(df.sample(i));
00478         adc[i]=samp_crystal.adc() ;    
00479         adcG[i]=samp_crystal.gainId();   
00480         adc[i]*=adcG[i];
00481         
00482         if (i==0) adcGain=adcG[i];
00483         if (i>0) adcGain=TMath::Max(adcG[i],adcGain);  
00484       }
00485       
00486       APDPulse->setPulse(adc);
00487 
00488       // Quality checks
00489       //================
00490       
00491       if(adcGain!=1) nEvtBadGain[channel]++;   
00492       if(!APDPulse->isTimingQualOK()) nEvtBadTiming[channel]++;
00493       nEvtTot[channel]++;
00494       
00495       // Fill if Pulse is fine
00496       //=======================
00497       
00498       if( APDPulse->isPulseOK() && lightside==side){
00499         
00500         Delta01->addEntry(APDPulse->getDelta(0,1));
00501         Delta12->addEntry(APDPulse->getDelta(1,2));
00502         
00503         if( nevtAB[channel] < _nevtmax && _fitab ){
00504           if(doesABTreeExist)  shapana -> putAllVals(channel, adc, eta, phi);   
00505           else shapana -> putAllVals(channel, adc, eta, phi, dccID, side, towerID, channelID);  
00506           nevtAB[channel]++ ;
00507         }
00508       }
00509     }
00510   }
00511 }// analyze
00512 
00513 
00514 //========================================================================
00515 void EcalABAnalyzer::endJob() {
00516 //========================================================================
00517 
00518   std::cout <<  "\n\t+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+" << std::endl;
00519   std::cout <<   "\t+=+     Analyzing data: getting (alpha, beta)     +=+" << std::endl;
00520 
00521   // Adjust channel numbers for EE 
00522   //===============================
00523 
00524   if( _ecalPart == "EE" ) {
00525     nCrys=channelMapEE.size();
00526     shapana->set_nch(nCrys);
00527   }
00528     
00529 
00530   // Set presamples number 
00531   //======================
00532   
00533   double delta01=Delta01->getMean();
00534   double delta12=Delta12->getMean();
00535   if(delta12>_presamplecut) {
00536     _presample=2;
00537     if(delta01>_presamplecut) _presample=1;
00538   }
00539   
00540   APDPulse->setPresamples(_presample);
00541   shapana->set_presample(_presample);
00542 
00543   //  Get alpha and beta
00544   //======================
00545 
00546   if(_fitab){
00547     std::cout <<  "\n\t+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+" << std::endl;
00548     std::cout <<   "\t+=+     Analyzing data: getting (alpha, beta)     +=+" << std::endl;
00549     TFile *fAB=0; TTree *ABInit=0;
00550     if(doesABTreeExist){
00551       fAB=new TFile(alphainitfile.c_str());
00552     }
00553     if(doesABTreeExist && fAB){
00554       ABInit = (TTree*) fAB->Get("ABCol0");
00555     } 
00556     shapana->computeShape(alphafile, ABInit);
00557 
00558   // Set quality flags for gains and timing
00559 
00560   double BadGainEvtPercentage=0.0;
00561   double BadTimingEvtPercentage=0.0;
00562   
00563   int nChanBadGain=0;
00564   int nChanBadTiming=0;
00565   
00566   for (unsigned int i=0;i<nCrys;i++){
00567     if(nEvtTot[i]!=0){
00568       BadGainEvtPercentage=double(nEvtBadGain[i])/double(nEvtTot[i]);
00569       BadTimingEvtPercentage=double(nEvtBadTiming[i])/double(nEvtTot[i]);
00570     }
00571     if(BadGainEvtPercentage>_qualpercent) {
00572       wasGainOK[i]=false;
00573       nChanBadGain++;
00574     }
00575     if(BadTimingEvtPercentage>_qualpercent){
00576       wasTimingOK[i]=false;
00577       nChanBadTiming++;
00578     }
00579   }
00580   
00581   double BadGainChanPercentage=double(nChanBadGain)/double(nCrys);
00582   double BadTimingChanPercentage=double(nChanBadTiming)/double(nCrys);
00583   
00584   if(BadGainChanPercentage>_qualpercent) isGainOK = false;
00585   if(BadTimingChanPercentage>_qualpercent) isTimingOK = false;
00586 
00587   
00588   if( !isGainOK )
00589   std::cout <<    "\t+=+ ............................ WARNING! APD GAIN WAS NOT 1    +=+" << std::endl;
00590   if( !isTimingOK )
00591   std::cout <<    "\t+=+ ............................ WARNING! TIMING WAS BAD        +=+" << std::endl;
00592 
00593   std::cout <<    "\t+=+    .................................... done  +=+" << std::endl;
00594   std::cout <<    "\t+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+" << std::endl;
00595   }
00596   
00597 
00598 }
00599 
00600 DEFINE_FWK_MODULE(EcalABAnalyzer);
00601