CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_1/src/CalibMuon/DTCalibration/plugins/DTTTrigWriter.cc

Go to the documentation of this file.
00001 /*
00002  *  See header file for a description of this class.
00003  *
00004  *  $Date: 2010/08/02 16:11:35 $
00005  *  $Revision: 1.5 $
00006  *  \author S. Bolognesi
00007  */
00008 
00009 #include "CalibMuon/DTCalibration/plugins/DTTTrigWriter.h"
00010 #include "CalibMuon/DTCalibration/interface/DTTimeBoxFitter.h"
00011 #include "CalibMuon/DTCalibration/interface/DTCalibDBUtils.h"
00012 
00013 
00014 
00015 #include "FWCore/Framework/interface/EventSetup.h"
00016 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00017 #include "FWCore/Framework/interface/Event.h"
00018 #include "FWCore/Framework/interface/ESHandle.h"
00019 #include "DataFormats/MuonDetId/interface/DTSuperLayerId.h"
00020 #include "Geometry/DTGeometry/interface/DTGeometry.h"
00021 #include "Geometry/Records/interface/MuonGeometryRecord.h"
00022 
00023 #include "CondFormats/DTObjects/interface/DTTtrig.h"
00024 #include "CondFormats/DataRecord/interface/DTStatusFlagRcd.h"
00025 #include "CondFormats/DTObjects/interface/DTStatusFlag.h"
00026 
00027 /* C++ Headers */
00028 #include <vector> 
00029 #include <iostream>
00030 #include <fstream>
00031 #include <string>
00032 #include <sstream>
00033 #include "TFile.h"
00034 #include "TH1.h"
00035 
00036 using namespace std;
00037 using namespace edm;
00038 
00039 
00040 // Constructor
00041 DTTTrigWriter::DTTTrigWriter(const ParameterSet& pset) {
00042   // get selected debug option
00043   debug = pset.getUntrackedParameter<bool>("debug",false);
00044 
00045   // Open the root file which contains the histos
00046   theRootInputFile = pset.getUntrackedParameter<string>("rootFileName");
00047   theFile = new TFile(theRootInputFile.c_str(), "READ");
00048   theFile->cd();
00049   theFitter = new DTTimeBoxFitter();
00050   if(debug)
00051     theFitter->setVerbosity(1);
00052 
00053   double sigmaFit = pset.getUntrackedParameter<double>("sigmaTTrigFit",10.);
00054   theFitter->setFitSigma(sigmaFit);
00055 
00056   // the kfactor to be uploaded in the ttrig DB
00057   kFactor = pset.getUntrackedParameter<double>("kFactor",-0.7);
00058 
00059   // Create the object to be written to DB
00060   tTrig = new DTTtrig();
00061   
00062   if(debug)
00063     cout << "[DTTTrigWriter]Constructor called!" << endl;
00064 }
00065 
00066 
00067 
00068 // Destructor
00069 DTTTrigWriter::~DTTTrigWriter(){
00070   if(debug)
00071     cout << "[DTTTrigWriter]Destructor called!" << endl;
00072   theFile->Close();
00073   delete theFitter;
00074 }
00075 
00076 
00077 
00078 // Do the job
00079 void DTTTrigWriter::analyze(const Event & event, const EventSetup& eventSetup) {
00080   if(debug)
00081     cout << "[DTTTrigWriter]Analyzer called!" << endl;
00082 
00083   // Get the DT Geometry
00084   ESHandle<DTGeometry> dtGeom;
00085   eventSetup.get<MuonGeometryRecord>().get(dtGeom);
00086 
00087   // Get all the sls from the setup
00088   const vector<DTSuperLayer*> superLayers = dtGeom->superLayers(); 
00089     
00090   // Loop over all SLs
00091   for(vector<DTSuperLayer*>::const_iterator  sl = superLayers.begin();
00092       sl != superLayers.end(); sl++) {
00093       
00094     // Get the histo from file
00095     DTSuperLayerId slId = (*sl)->id();
00096     TH1F* histo = (TH1F*)theFile->Get((getTBoxName(slId)).c_str());
00097     if(histo) { // Check that the histo exists
00098       // Compute mean and sigma of the rising edge
00099       pair<double, double> meanAndSigma = theFitter->fitTimeBox(histo);
00100 
00101       // Write them in DB object
00102       tTrig->set(slId,
00103                  meanAndSigma.first,
00104                  meanAndSigma.second,
00105                  kFactor,
00106                  DTTimeUnits::ns);
00107       if(debug) {
00108         cout << " SL: " << slId
00109              << " mean = " << meanAndSigma.first
00110              << " sigma = " << meanAndSigma.second << endl;
00111       }
00112     }
00113   }
00114 }
00115 
00116 
00117 
00118 // Write objects to DB
00119 void DTTTrigWriter::endJob() {
00120   if(debug) 
00121         cout << "[DTTTrigWriter]Writing ttrig object to DB!" << endl;
00122 
00123   // FIXME: to be read from cfg?
00124   string tTrigRecord = "DTTtrigRcd";
00125   
00126   // Write the object to DB
00127   DTCalibDBUtils::writeToDB(tTrigRecord, tTrig);
00128 
00129 }  
00130 
00131 
00132 
00133 // Compute the name of the time box histo
00134 string DTTTrigWriter::getTBoxName(const DTSuperLayerId& slId) const {
00135   string histoName;
00136   stringstream theStream;
00137   theStream << "Ch_" << slId.wheel() << "_" << slId.station() << "_" << slId.sector()
00138             << "_SL" << slId.superlayer() << "_hTimeBox";
00139   theStream >> histoName;
00140   return histoName;
00141 }