#include <CalibMuon/DTCalibration/plugins/DTVDriftWriter.h>
Public Member Functions | |
void | analyze (const edm::Event &event, const edm::EventSetup &eventSetup) |
DTVDriftWriter (const edm::ParameterSet &pset) | |
Constructor. | |
void | endJob () |
virtual | ~DTVDriftWriter () |
Destructor. | |
Private Attributes | |
bool | debug |
edm::ParameterSet | theCalibFilePar |
TFile * | theFile |
DTMeanTimerFitter * | theFitter |
std::string | theGranularity |
DTMtime * | theMTime |
std::string | theRootInputFile |
std::string | theVDriftOutputFile |
Definition at line 32 of file DTVDriftWriter.h.
DTVDriftWriter::DTVDriftWriter | ( | const edm::ParameterSet & | pset | ) |
Constructor.
Definition at line 42 of file DTVDriftWriter.cc.
References GenMuonPlsPt100GeV_cfg::cout, debug, lat::endl(), edm::ParameterSet::getUntrackedParameter(), theCalibFilePar, theFile, theFitter, theGranularity, theMTime, theRootInputFile, and theVDriftOutputFile.
00042 { 00043 // get selected debug option 00044 debug = pset.getUntrackedParameter<bool>("debug", "false"); 00045 00046 // Open the root file which contains the histos 00047 theRootInputFile = pset.getUntrackedParameter<string>("rootFileName"); 00048 theFile = new TFile(theRootInputFile.c_str(), "READ"); 00049 00050 theFitter = new DTMeanTimerFitter(theFile); 00051 if(debug) 00052 theFitter->setVerbosity(1); 00053 00054 // the text file which will contain the histos 00055 theVDriftOutputFile = pset.getUntrackedParameter<string>("vDriftFileName"); 00056 00057 // get parameter set for DTCalibrationMap constructor 00058 theCalibFilePar = pset.getUntrackedParameter<ParameterSet>("calibFileConfig"); 00059 00060 // the granularity to be used for calib constants evaluation 00061 theGranularity = pset.getUntrackedParameter<string>("calibGranularity","bySL"); 00062 00063 theMTime = new DTMtime(); 00064 00065 if(debug) 00066 cout << "[DTVDriftWriter]Constructor called!" << endl; 00067 }
DTVDriftWriter::~DTVDriftWriter | ( | ) | [virtual] |
Destructor.
Definition at line 70 of file DTVDriftWriter.cc.
References GenMuonPlsPt100GeV_cfg::cout, debug, lat::endl(), theFile, and theFitter.
00070 { 00071 if(debug) 00072 cout << "[DTVDriftWriter]Destructor called!" << endl; 00073 theFile->Close(); 00074 delete theFitter; 00075 }
void DTVDriftWriter::analyze | ( | const edm::Event & | event, | |
const edm::EventSetup & | eventSetup | |||
) | [virtual] |
Implements edm::EDAnalyzer.
Definition at line 77 of file DTVDriftWriter.cc.
References DTCalibrationMap::addCell(), GenMuonPlsPt100GeV_cfg::cout, debug, lat::endl(), DTMeanTimerFitter::evaluateVDriftAndReso(), edm::EventSetup::get(), DTCalibrationMap::getConsts(), DTCalibrationMap::getKey(), N, DTTimeUnits::ns, DTChamberId::sector(), DTMtime::set(), DTChamberId::station(), DTSuperLayerId::superLayer(), theCalibFilePar, theFitter, theGranularity, theMTime, theVDriftOutputFile, DTChamberId::wheel(), and DTCalibrationMap::writeConsts().
00077 { 00078 if(debug) 00079 cout << "[DTVDriftWriter]Analyzer called!" << endl; 00080 00081 // Instantiate a DTCalibrationMap object 00082 DTCalibrationMap calibValuesFile(theCalibFilePar); 00083 00084 // Get the DT Geometry 00085 ESHandle<DTGeometry> dtGeom; 00086 eventSetup.get<MuonGeometryRecord>().get(dtGeom); 00087 00088 if(theGranularity == "bySL") { 00089 // Get all the sls from the setup 00090 const vector<DTSuperLayer*> superLayers = dtGeom->superLayers(); 00091 00092 // Loop over all SLs 00093 for(vector<DTSuperLayer*>::const_iterator slCell = superLayers.begin(); 00094 slCell != superLayers.end(); slCell++) { 00095 00096 DTSuperLayerId slId = (*slCell)->id(); 00097 // evaluate v_drift and sigma from the TMax histograms 00098 DTWireId wireId(slId, 0, 0); 00099 vector<float> newConstants; 00100 TString N=(((((TString) "TMax"+(long) wireId.wheel()) +(long) wireId.station()) 00101 +(long) wireId.sector())+(long) wireId.superLayer()); 00102 vector<float> vDriftAndReso = theFitter->evaluateVDriftAndReso(N); 00103 00104 // Don't write the constants for the SL if the vdrift was not computed 00105 if(vDriftAndReso.front() == -1) 00106 continue; 00107 00108 const DTCalibrationMap::CalibConsts* oldConstants = calibValuesFile.getConsts(wireId); 00109 00110 if(oldConstants != 0) { 00111 newConstants.push_back((*oldConstants)[0]); 00112 newConstants.push_back((*oldConstants)[1]); 00113 } else { 00114 newConstants.push_back(-1); 00115 newConstants.push_back(-1); 00116 } 00117 for(int ivd=0; ivd<=5;ivd++) { 00118 // 0=vdrift, 1=reso, 2=(3deltat0-2deltat0), 3=(2deltat0-1deltat0), 00119 // 4=(1deltat0-0deltat0), 5=deltat0 from hists with max entries, 00120 newConstants.push_back(vDriftAndReso[ivd]); 00121 } 00122 calibValuesFile.addCell(calibValuesFile.getKey(wireId), newConstants); 00123 00124 theMTime->set(slId, 00125 vDriftAndReso[0], 00126 vDriftAndReso[1], 00127 DTTimeUnits::ns); 00128 if(debug) { 00129 cout << " SL: " << slId 00130 << " vDrift = " << vDriftAndReso[0] 00131 << " reso = " << vDriftAndReso[1] << endl; 00132 } 00133 } 00134 } 00135 // to be implemented: granularity different from bySL 00136 00137 // if(theGranularity == "byChamber") { 00138 // const vector<DTChamber*> chambers = dMap.chambers(); 00139 00140 // // Loop over all chambers 00141 // for(vector<MuBarChamber*>::const_iterator chamber = chambers.begin(); 00142 // chamber != chambers.end(); chamber ++) { 00143 // MuBarChamberId chamber_id = (*chamber)->id(); 00144 // MuBarDigiParameters::Key wire_id(chamber_id, 0, 0, 0); 00145 // vector<float> newConstants; 00146 // vector<float> vDriftAndReso = evaluateVDriftAndReso(wire_id, f); 00147 // const CalibConsts* oldConstants = digiParams.getConsts(wire_id); 00148 // if(oldConstants !=0) { 00149 // newConstants = *oldConstants; 00150 // newConstants.push_back(vDriftAndReso[0]); 00151 // newConstants.push_back(vDriftAndReso[1]); 00152 // newConstants.push_back(vDriftAndReso[2]); 00153 // newConstants.push_back(vDriftAndReso[3]); 00154 // } else { 00155 // newConstants.push_back(-1); 00156 // newConstants.push_back(-1); 00157 // newConstants.push_back(vDriftAndReso[0]); 00158 // newConstants.push_back(vDriftAndReso[1]); 00159 // newConstants.push_back(vDriftAndReso[2]); 00160 // newConstants.push_back(vDriftAndReso[3]); 00161 // } 00162 // digiParams.addCell(wire_id, newConstants); 00163 // } 00164 // } 00165 //write values to a table 00166 calibValuesFile.writeConsts(theVDriftOutputFile); 00167 }
Reimplemented from edm::EDAnalyzer.
Definition at line 171 of file DTVDriftWriter.cc.
References GenMuonPlsPt100GeV_cfg::cout, debug, lat::endl(), ecalRecalibSequence_cff::record, and theMTime.
00171 { 00172 00173 if(debug) 00174 cout << "[DTVDriftWriter]Writing vdrift object to DB!" << endl; 00175 00176 // Write the ttrig object to DB 00177 string record = "DTMtimeRcd"; 00178 DTCalibDBUtils::writeToDB<DTMtime>(record, theMTime); 00179 }
bool DTVDriftWriter::debug [private] |
Definition at line 49 of file DTVDriftWriter.h.
Referenced by analyze(), DTVDriftWriter(), endJob(), and ~DTVDriftWriter().
TFile* DTVDriftWriter::theFile [private] |
Definition at line 52 of file DTVDriftWriter.h.
Referenced by DTVDriftWriter(), and ~DTVDriftWriter().
DTMeanTimerFitter* DTVDriftWriter::theFitter [private] |
Definition at line 67 of file DTVDriftWriter.h.
Referenced by analyze(), DTVDriftWriter(), and ~DTVDriftWriter().
std::string DTVDriftWriter::theGranularity [private] |
DTMtime* DTVDriftWriter::theMTime [private] |
Definition at line 70 of file DTVDriftWriter.h.
Referenced by analyze(), DTVDriftWriter(), and endJob().
std::string DTVDriftWriter::theRootInputFile [private] |
std::string DTVDriftWriter::theVDriftOutputFile [private] |