CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_2_SLHC4_patch1/src/CalibMuon/DTCalibration/plugins/DTT0Calibration.cc

Go to the documentation of this file.
00001 /*
00002  *  See header file for a description of this class.
00003  *
00004  *  $Date: 2012/05/11 17:17:17 $
00005  *  $Revision: 1.13 $
00006  *  \author S. Bolognesi - INFN Torino
00007  */
00008 #include "CalibMuon/DTCalibration/plugins/DTT0Calibration.h"
00009 #include "CalibMuon/DTCalibration/interface/DTCalibDBUtils.h"
00010 
00011 #include "FWCore/Framework/interface/Event.h"
00012 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00013 
00014 #include "Geometry/Records/interface/MuonGeometryRecord.h"
00015 #include "Geometry/Records/interface/MuonNumberingRecord.h"
00016 
00017 #include "DataFormats/DTDigi/interface/DTDigiCollection.h"
00018 #include "CondFormats/DTObjects/interface/DTT0.h"
00019 
00020 #include "TH1I.h"
00021 #include "TFile.h"
00022 #include "TKey.h"
00023 
00024 using namespace std;
00025 using namespace edm;
00026 // using namespace cond;
00027 
00028 // Constructor
00029 DTT0Calibration::DTT0Calibration(const edm::ParameterSet& pset) {
00030   // Get the debug parameter for verbose output
00031   debug = pset.getUntrackedParameter<bool>("debug");
00032   if(debug) 
00033     cout << "[DTT0Calibration]Constructor called!" << endl;
00034 
00035   // Get the label to retrieve digis from the event
00036   digiLabel = pset.getUntrackedParameter<string>("digiLabel");
00037 
00038   // The root file which contain the histos per layer
00039   string rootFileName = pset.getUntrackedParameter<string>("rootFileName","DTT0PerLayer.root");
00040   theFile = new TFile(rootFileName.c_str(), "RECREATE");
00041  
00042   theCalibWheel =  pset.getUntrackedParameter<string>("calibWheel", "All"); //FIXME amke a vector of integer instead of a string
00043   if(theCalibWheel != "All") {
00044     stringstream linestr;
00045     int selWheel;
00046     linestr << theCalibWheel;
00047     linestr >> selWheel;
00048     cout << "[DTT0CalibrationPerLayer] chosen wheel " << selWheel << endl;
00049   }
00050 
00051   // Sector/s to calibrate
00052   theCalibSector =  pset.getUntrackedParameter<string>("calibSector", "All"); //FIXME amke a vector of integer instead of a string
00053   if(theCalibSector != "All") {
00054     stringstream linestr;
00055     int selSector;
00056     linestr << theCalibSector;
00057     linestr >> selSector;
00058     cout << "[DTT0CalibrationPerLayer] chosen sector " << selSector << endl;
00059   }
00060 
00061   vector<string> defaultCell;
00062   defaultCell.push_back("None");
00063   cellsWithHistos = pset.getUntrackedParameter<vector<string> >("cellsWithHisto", defaultCell);
00064   for(vector<string>::const_iterator cell = cellsWithHistos.begin(); cell != cellsWithHistos.end(); cell++){
00065     if((*cell) != "None"){
00066       stringstream linestr;
00067       int wheel,sector,station,sl,layer,wire;
00068       linestr << (*cell);
00069       linestr >> wheel >> sector >> station >> sl >> layer >> wire;
00070       wireIdWithHistos.push_back(DTWireId(wheel,station,sector,sl,layer,wire));
00071     }
00072   }
00073 
00074   hT0SectorHisto=0;
00075 
00076   nevents=0;
00077   eventsForLayerT0 = pset.getParameter<unsigned int>("eventsForLayerT0");
00078   eventsForWireT0 = pset.getParameter<unsigned int>("eventsForWireT0");
00079   rejectDigiFromPeak = pset.getParameter<unsigned int>("rejectDigiFromPeak");
00080   tpPeakWidth = pset.getParameter<double>("tpPeakWidth");
00081   //useReferenceWireInLayer_ = true;
00082   correctByChamberMean_ = pset.getParameter<bool>("correctByChamberMean");
00083 }
00084 
00085 // Destructor
00086 DTT0Calibration::~DTT0Calibration(){
00087   if(debug) 
00088     cout << "[DTT0Calibration]Destructor called!" << endl;
00089 
00090   theFile->Close();
00091 }
00092 
00094 void DTT0Calibration::analyze(const edm::Event & event, const edm::EventSetup& eventSetup) {
00095   if(debug || event.id().event() % 500==0)
00096     cout << "--- [DTT0Calibration] Analysing Event: #Run: " << event.id().run()
00097          << " #Event: " << event.id().event() << endl;
00098   nevents++;
00099 
00100   // Get the digis from the event
00101   Handle<DTDigiCollection> digis; 
00102   event.getByLabel(digiLabel, digis);
00103 
00104   // Get the DT Geometry
00105   eventSetup.get<MuonGeometryRecord>().get(dtGeom);
00106 
00107   // Iterate through all digi collections ordered by LayerId   
00108   DTDigiCollection::DigiRangeIterator dtLayerIt;
00109   for (dtLayerIt = digis->begin();
00110        dtLayerIt != digis->end();
00111        ++dtLayerIt){
00112     // Get the iterators over the digis associated with this LayerId
00113     const DTDigiCollection::Range& digiRange = (*dtLayerIt).second;
00114   
00115     // Get the layerId
00116     const DTLayerId layerId = (*dtLayerIt).first; //FIXME: check to be in the right sector
00117 
00118     if((theCalibWheel != "All") && (layerId.superlayerId().chamberId().wheel() != selWheel))
00119       continue;
00120     if((theCalibSector != "All") && (layerId.superlayerId().chamberId().sector() != selSector))
00121       continue;
00122  
00123     //if(debug) {
00124     //  cout << "Layer " << layerId<<" with "<<distance(digiRange.first, digiRange.second)<<" digi"<<endl;
00125     //}
00126 
00127     // Loop over all digis in the given layer
00128     for (DTDigiCollection::const_iterator digi = digiRange.first;
00129          digi != digiRange.second;
00130          digi++) {
00131       double t0 = (*digi).countsTDC();
00132 
00133       //Use first bunch of events to fill t0 per layer
00134       if(nevents < eventsForLayerT0){
00135         //Get the per-layer histo from the map
00136         TH1I *hT0LayerHisto = theHistoLayerMap[layerId];
00137         //If it doesn't exist, book it
00138         if(hT0LayerHisto == 0){
00139           theFile->cd();
00140           hT0LayerHisto = new TH1I(getHistoName(layerId).c_str(),
00141                                    "T0 from pulses by layer (TDC counts, 1 TDC count = 0.781 ns)",
00142                                    200, t0-100, t0+100);
00143           if(debug)
00144             cout << "  New T0 per Layer Histo: " << hT0LayerHisto->GetName() << endl;
00145           theHistoLayerMap[layerId] = hT0LayerHisto;
00146         }
00147     
00148         //Fill the histos
00149         theFile->cd();
00150         if(hT0LayerHisto != 0) {
00151           //  if(debug)
00152           // cout<<"Filling histo "<<hT0LayerHisto->GetName()<<" with digi "<<t0<<" TDC counts"<<endl;
00153           hT0LayerHisto->Fill(t0);
00154         }
00155       }
00156 
00157       //Use all the remaining events to compute t0 per wire
00158       if(nevents>eventsForLayerT0){
00159         // Get the wireId
00160         const DTWireId wireId(layerId, (*digi).wire());
00161         if(debug) {
00162           cout << "   Wire: " << wireId << endl
00163                << "       time (TDC counts): " << (*digi).countsTDC()<< endl;
00164         }   
00165 
00166         //Fill the histos per wire for the chosen cells
00167         vector<DTWireId>::iterator it_wire = find(wireIdWithHistos.begin(),wireIdWithHistos.end(),wireId);
00168         if(it_wire != wireIdWithHistos.end()){
00169           if(theHistoWireMap.find(wireId) == theHistoWireMap.end()){
00170             theHistoWireMap[wireId] = new TH1I(getHistoName(wireId).c_str(),"T0 from pulses by wire (TDC counts, 1 TDC count = 0.781 ns)",7000,0,7000);
00171             if(debug) cout << "  New T0 per wire Histo: " << (theHistoWireMap[wireId])->GetName() << endl;
00172           }
00173           if(theHistoWireMap_ref.find(wireId) == theHistoWireMap_ref.end()){
00174             theHistoWireMap_ref[wireId] = new TH1I((getHistoName(wireId) + "_ref").c_str(),"T0 from pulses by wire (TDC counts, 1 TDC count = 0.781 ns)",7000,0,7000);
00175             if(debug) cout << "  New T0 per wire Histo: " << (theHistoWireMap_ref[wireId])->GetName() << endl;
00176           }
00177 
00178           TH1I* hT0WireHisto = theHistoWireMap[wireId];
00179           //Fill the histos
00180           theFile->cd();
00181           if(hT0WireHisto) hT0WireHisto->Fill(t0);
00182         }
00183 
00184         //Check the tzero has reasonable value
00185         if(abs(hT0SectorHisto->GetBinCenter(hT0SectorHisto->GetMaximumBin()) - t0) > rejectDigiFromPeak){
00186           if(debug)
00187             cout<<"digi skipped because t0 per sector "<<hT0SectorHisto->GetBinCenter(hT0SectorHisto->GetMaximumBin())<<endl;
00188           continue;
00189         }
00190 
00191         //Use second bunch of events to compute a t0 reference per wire
00192         if(nevents< (eventsForLayerT0 + eventsForWireT0)){
00193           //Fill reference wire histos
00194           if(it_wire != wireIdWithHistos.end()){
00195             TH1I* hT0WireHisto_ref = theHistoWireMap_ref[wireId];
00196             theFile->cd();
00197             if(hT0WireHisto_ref) hT0WireHisto_ref->Fill(t0); 
00198           } 
00199           if(!nDigiPerWire_ref[wireId]){
00200             mK_ref[wireId] = 0;
00201           }
00202           nDigiPerWire_ref[wireId] = nDigiPerWire_ref[wireId] + 1;
00203           mK_ref[wireId] = mK_ref[wireId] + (t0-mK_ref[wireId])/nDigiPerWire_ref[wireId];
00204         }
00205         //Use last all the remaining events to compute the mean and sigma t0 per wire
00206         else if(nevents>(eventsForLayerT0 + eventsForWireT0)){
00207           if(abs(t0-mK_ref[wireId]) > tpPeakWidth)
00208             continue;
00209           if(!nDigiPerWire[wireId]){
00210             theAbsoluteT0PerWire[wireId] = 0;
00211             qK[wireId] = 0;
00212             mK[wireId] = 0;
00213           }
00214           nDigiPerWire[wireId] = nDigiPerWire[wireId] + 1;
00215           theAbsoluteT0PerWire[wireId] = theAbsoluteT0PerWire[wireId] + t0;
00216           //theSigmaT0PerWire[wireId] = theSigmaT0PerWire[wireId] + (t0*t0);
00217           qK[wireId] = qK[wireId] + ((nDigiPerWire[wireId]-1)*(t0-mK[wireId])*(t0-mK[wireId])/nDigiPerWire[wireId]);
00218           mK[wireId] = mK[wireId] + (t0-mK[wireId])/nDigiPerWire[wireId];
00219         }
00220       }//end if(nevents>1000)
00221     }//end loop on digi
00222   }//end loop on layer
00223 
00224   //Use the t0 per layer histos to have an indication about the t0 position 
00225   if(nevents == eventsForLayerT0){
00226     for(map<DTLayerId, TH1I*>::const_iterator lHisto = theHistoLayerMap.begin();
00227         lHisto != theHistoLayerMap.end();
00228         lHisto++){
00229       if(debug)
00230         cout<<"Reading histogram "<<(*lHisto).second->GetName()<<" with mean "<<(*lHisto).second->GetMean()<<" and RMS "<<(*lHisto).second->GetRMS();
00231 
00232       //Take the mean as a first t0 estimation
00233       if((*lHisto).second->GetRMS()<5.0){
00234         if(hT0SectorHisto == 0){
00235           hT0SectorHisto = new TH1D("hT0AllLayerOfSector","T0 from pulses per layer in sector", 
00236                                     //20, (*lHisto).second->GetMean()-100, (*lHisto).second->GetMean()+100);
00237                                     700, 0, 7000);
00238         }
00239         if(debug)
00240           cout<<" accepted"<<endl;
00241         hT0SectorHisto->Fill((*lHisto).second->GetMean());
00242       }
00243       //Take the mean of noise + 400ns as a first t0 estimation
00244       // if((*lHisto).second->GetRMS()>10.0 && ((*lHisto).second->GetRMS()<15.0)){
00245 //      double t0_estim = (*lHisto).second->GetMean() + 400;
00246 //      if(hT0SectorHisto == 0){
00247 //        hT0SectorHisto = new TH1D("hT0AllLayerOfSector","T0 from pulses per layer in sector", 
00248 //                                  //20, t0_estim-100, t0_estim+100);
00249 //                                  700, 0, 7000);
00250 //      }
00251 //      if(debug)
00252 //        cout<<" accepted + 400ns"<<endl;
00253 //      hT0SectorHisto->Fill((*lHisto).second->GetMean() + 400);
00254 //       }
00255       if(debug)
00256         cout<<endl;
00257 
00258       theT0LayerMap[(*lHisto).second->GetName()] = (*lHisto).second->GetMean();
00259       theSigmaT0LayerMap[(*lHisto).second->GetName()] = (*lHisto).second->GetRMS();
00260     }
00261     if(!hT0SectorHisto){
00262       cout<<"[DTT0Calibration]: All the t0 per layer are still uncorrect: trying with greater number of events"<<endl;
00263       eventsForLayerT0 = eventsForLayerT0*2;
00264       return;
00265     }
00266     if(debug)
00267       cout<<"[DTT0Calibration] t0 reference for this sector "<<
00268         hT0SectorHisto->GetBinCenter(hT0SectorHisto->GetMaximumBin())<<endl;
00269   } 
00270 }
00271 
00272 
00273 void DTT0Calibration::endJob() {
00274 
00275   DTT0* t0sAbsolute = new DTT0();
00276   DTT0* t0sRelative = new DTT0();
00277   DTT0* t0sWRTChamber = new DTT0();
00278 
00279   //if(debug) 
00280     cout << "[DTT0CalibrationPerLayer]Writing histos to file!" << endl;
00281 
00282   theFile->cd();
00283   hT0SectorHisto->Write();
00284   for(map<DTWireId, TH1I*>::const_iterator wHisto = theHistoWireMap.begin();
00285       wHisto != theHistoWireMap.end();
00286       wHisto++) {
00287     (*wHisto).second->Write(); 
00288   }
00289   for(map<DTWireId, TH1I*>::const_iterator wHisto = theHistoWireMap_ref.begin();
00290       wHisto != theHistoWireMap_ref.end();
00291       wHisto++) {
00292     (*wHisto).second->Write();
00293   }
00294   for(map<DTLayerId, TH1I*>::const_iterator lHisto = theHistoLayerMap.begin();
00295       lHisto != theHistoLayerMap.end();
00296       lHisto++) {
00297     (*lHisto).second->Write(); 
00298   }  
00299 
00300   //if(debug) 
00301     cout << "[DTT0Calibration] Compute and store t0 and sigma per wire" << endl;
00302 
00303   for(map<DTWireId, double>::const_iterator wiret0 = theAbsoluteT0PerWire.begin();
00304       wiret0 != theAbsoluteT0PerWire.end();
00305       wiret0++){
00306     if(nDigiPerWire[(*wiret0).first]){
00307       double t0 = (*wiret0).second/nDigiPerWire[(*wiret0).first];
00308 
00309       theRelativeT0PerWire[(*wiret0).first] = t0 - hT0SectorHisto->GetBinCenter(hT0SectorHisto->GetMaximumBin());
00310 
00311       //theSigmaT0PerWire[(*wiret0).first] = sqrt((theSigmaT0PerWire[(*wiret0).first] / nDigiPerWire[(*wiret0).first]) - t0*t0);
00312       theSigmaT0PerWire[(*wiret0).first] = sqrt(qK[(*wiret0).first]/nDigiPerWire[(*wiret0).first]);
00313 
00314       cout << "Wire " << (*wiret0).first << " has t0 " << t0 << "(absolute) "
00315                                                        << theRelativeT0PerWire[(*wiret0).first] << "(relative)"
00316                                          << "    sigma " << theSigmaT0PerWire[(*wiret0).first] << endl;
00317 
00318       t0sAbsolute->set((*wiret0).first, t0, theSigmaT0PerWire[(*wiret0).first],DTTimeUnits::counts); 
00319     }
00320     else{
00321       cout<<"[DTT0Calibration] ERROR: no digis in wire "<<(*wiret0).first<<endl;
00322       abort();
00323     }
00324   }
00325 
00326   if(correctByChamberMean_){
00328      // Get all the sls from the setup
00329      const vector<DTSuperLayer*> superLayers = dtGeom->superLayers();     
00330      // Loop over all SLs
00331      for(vector<DTSuperLayer*>::const_iterator  sl = superLayers.begin();
00332            sl != superLayers.end(); sl++) {
00333 
00334 
00335         //Compute mean for odd and even superlayers
00336         double oddLayersMean=0;
00337         double evenLayersMean=0; 
00338         double oddLayersDen=0;
00339         double evenLayersDen=0;
00340         for(map<DTWireId, double>::const_iterator wiret0 = theRelativeT0PerWire.begin();
00341               wiret0 != theRelativeT0PerWire.end();
00342               wiret0++){
00343            if((*wiret0).first.layerId().superlayerId() == (*sl)->id()){
00344               if(debug)
00345                  cout<<"[DTT0Calibration] Superlayer "<<(*sl)->id()
00346                     <<"layer " <<(*wiret0).first.layerId().layer()<<" with "<<(*wiret0).second<<endl;
00347               if(((*wiret0).first.layerId().layer()) % 2){
00348                  oddLayersMean = oddLayersMean + (*wiret0).second;
00349                  oddLayersDen++;
00350               }
00351               else{
00352                  evenLayersMean = evenLayersMean + (*wiret0).second;
00353                  evenLayersDen++;
00354               }
00355            }
00356         }
00357         oddLayersMean = oddLayersMean/oddLayersDen;
00358         evenLayersMean = evenLayersMean/evenLayersDen;
00359         //if(debug && oddLayersMean)
00360         cout<<"[DTT0Calibration] Relative T0 mean for  odd layers "<<oddLayersMean<<"  even layers"<<evenLayersMean<<endl;
00361 
00362         //Compute sigma for odd and even superlayers
00363         double oddLayersSigma=0;
00364         double evenLayersSigma=0;
00365         for(map<DTWireId, double>::const_iterator wiret0 = theRelativeT0PerWire.begin();
00366               wiret0 != theRelativeT0PerWire.end();
00367               wiret0++){
00368            if((*wiret0).first.layerId().superlayerId() == (*sl)->id()){
00369               if(((*wiret0).first.layerId().layer()) % 2){
00370                  oddLayersSigma = oddLayersSigma + ((*wiret0).second - oddLayersMean) * ((*wiret0).second - oddLayersMean);
00371               }
00372               else{
00373                  evenLayersSigma = evenLayersSigma + ((*wiret0).second - evenLayersMean) * ((*wiret0).second - evenLayersMean);
00374               }
00375            }
00376         }
00377         oddLayersSigma = oddLayersSigma/oddLayersDen;
00378         evenLayersSigma = evenLayersSigma/evenLayersDen;
00379         oddLayersSigma = sqrt(oddLayersSigma);
00380         evenLayersSigma = sqrt(evenLayersSigma);
00381 
00382         //if(debug && oddLayersMean)
00383         cout<<"[DTT0Calibration] Relative T0 sigma for  odd layers "<<oddLayersSigma<<"  even layers"<<evenLayersSigma<<endl;
00384 
00385         //Recompute the mean for odd and even superlayers discarding fluctations
00386         double oddLayersFinalMean=0; 
00387         double evenLayersFinalMean=0;
00388         for(map<DTWireId, double>::const_iterator wiret0 = theRelativeT0PerWire.begin();
00389               wiret0 != theRelativeT0PerWire.end();
00390               wiret0++){
00391            if((*wiret0).first.layerId().superlayerId() == (*sl)->id()){
00392               if(((*wiret0).first.layerId().layer()) % 2){
00393                  if(abs((*wiret0).second - oddLayersMean) < (2*oddLayersSigma))
00394                     oddLayersFinalMean = oddLayersFinalMean + (*wiret0).second;
00395               }
00396               else{
00397                  if(abs((*wiret0).second - evenLayersMean) < (2*evenLayersSigma))
00398                     evenLayersFinalMean = evenLayersFinalMean + (*wiret0).second;
00399               }
00400            }
00401         }
00402         oddLayersFinalMean = oddLayersFinalMean/oddLayersDen;
00403         evenLayersFinalMean = evenLayersFinalMean/evenLayersDen;
00404         //if(debug && oddLayersMean)
00405         cout<<"[DTT0Calibration] Final relative T0 mean for  odd layers "<<oddLayersFinalMean<<"  even layers"<<evenLayersFinalMean<<endl;
00406 
00407         for(map<DTWireId, double>::const_iterator wiret0 = theRelativeT0PerWire.begin();
00408               wiret0 != theRelativeT0PerWire.end();
00409               wiret0++){
00410            if((*wiret0).first.layerId().superlayerId() == (*sl)->id()){
00411               double t0=-999;
00412               if(((*wiret0).first.layerId().layer()) % 2)
00413                  t0 = (*wiret0).second + (evenLayersFinalMean - oddLayersFinalMean);
00414               else
00415                  t0 = (*wiret0).second;
00416 
00417               cout << "[DTT0Calibration] Wire " << (*wiret0).first << " has t0 " << (*wiret0).second 
00418                    << " (relative, after even-odd layer corrections)  "
00419                    << "    sigma " << theSigmaT0PerWire[(*wiret0).first] << endl;
00420 
00421               //Store the results into DB
00422               t0sRelative->set((*wiret0).first, t0, theSigmaT0PerWire[(*wiret0).first],DTTimeUnits::counts); 
00423            }
00424         }
00425      }
00426 
00428      //if(debug) 
00429      cout << "[DTT0Calibration]Computing relative t0 wrt to chamber average" << endl;
00430      //Compute the reference for each chamber
00431      map<DTChamberId,double> sumT0ByChamber;
00432      map<DTChamberId,int> countT0ByChamber;
00433      for(DTT0::const_iterator tzero = t0sRelative->begin();
00434            tzero != t0sRelative->end(); tzero++) {
00435         int channelId = tzero->channelId;
00436         if ( channelId == 0 ) continue;
00437         DTWireId wireId(channelId);
00438         DTChamberId chamberId(wireId.chamberId());
00439         //sumT0ByChamber[chamberId] = sumT0ByChamber[chamberId] + tzero->t0mean;
00440         // @@@ better DTT0 usage
00441         float t0mean_f;
00442         float t0rms_f;
00443         t0sRelative->get(wireId,t0mean_f,t0rms_f,DTTimeUnits::counts);
00444         sumT0ByChamber[chamberId] = sumT0ByChamber[chamberId] + t0mean_f;
00445         // @@@ NEW DTT0 END
00446         countT0ByChamber[chamberId]++;
00447      }
00448 
00449      //Change reference for each wire and store the new t0s in the new map
00450      for(DTT0::const_iterator tzero = t0sRelative->begin();
00451            tzero != t0sRelative->end(); tzero++) {
00452         int channelId = tzero->channelId;
00453         if ( channelId == 0 ) continue;
00454         DTWireId wireId(channelId);
00455         DTChamberId chamberId(wireId.chamberId());
00456         //double t0mean = (tzero->t0mean) - (sumT0ByChamber[chamberId]/countT0ByChamber[chamberId]);
00457         //double t0rms = tzero->t0rms;
00458         // @@@ better DTT0 usage
00459         float t0mean_f;
00460         float t0rms_f;
00461         t0sRelative->get(wireId,t0mean_f,t0rms_f,DTTimeUnits::counts);
00462         double t0mean = t0mean_f - (sumT0ByChamber[chamberId]/countT0ByChamber[chamberId]);
00463         double t0rms = t0rms_f;
00464         // @@@ NEW DTT0 END
00465         t0sWRTChamber->set(wireId,
00466               t0mean,
00467               t0rms,
00468               DTTimeUnits::counts);
00469         //if(debug)
00470         //cout<<"Chamber "<<chamberId<<" has reference "<<(sumT0ByChamber[chamberId]/countT0ByChamber[chamberId]);
00471         cout << "Changing t0 of wire " << wireId << " from " << t0mean_f
00472              << " to " << t0mean << endl;
00473      }
00474   }
00475   
00477   if(debug) 
00478    cout << "[DTT0Calibration]Writing values in DB!" << endl;
00479   // FIXME: to be read from cfg?
00480   string t0Record = "DTT0Rcd";
00481   // Write the t0 map to DB
00482   if( correctByChamberMean_ ) DTCalibDBUtils::writeToDB(t0Record, t0sWRTChamber);
00483   else                        DTCalibDBUtils::writeToDB(t0Record, t0sAbsolute);
00484 }
00485 
00486 string DTT0Calibration::getHistoName(const DTWireId& wId) const {
00487   string histoName;
00488   stringstream theStream;
00489   theStream << "Ch_" << wId.wheel() << "_" << wId.station() << "_" << wId.sector()
00490             << "_SL" << wId.superlayer() << "_L" << wId.layer() << "_W"<< wId.wire() <<"_hT0Histo";
00491   theStream >> histoName;
00492   return histoName;
00493 }
00494 
00495 string DTT0Calibration::getHistoName(const DTLayerId& lId) const {
00496   string histoName;
00497   stringstream theStream;
00498   theStream << "Ch_" << lId.wheel() << "_" << lId.station() << "_" << lId.sector()
00499             << "_SL" << lId.superlayer() << "_L" << lId.layer() <<"_hT0Histo";
00500   theStream >> histoName;
00501   return histoName;
00502 }
00503