00001 #include <iostream>
00002 #include <fstream>
00003 #include <sstream>
00004 #include <vector>
00005 #include <string>
00006 #include <stdexcept>
00007 #include <cstring>
00008 #include <cstdlib>
00009 #include <boost/tuple/tuple.hpp>
00010 #include <boost/format.hpp>
00011
00012 #include "FWCore/Framework/interface/Event.h"
00013 #include "FWCore/Framework/interface/EventSetup.h"
00014 #include "FWCore/Framework/interface/ESHandle.h"
00015 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00016 #include "FWCore/Utilities/interface/EDMException.h"
00017
00018 #include "DetectorDescription/Core/interface/DDFilteredView.h"
00019 #include "DetectorDescription/Core/interface/DDCompactView.h"
00020 #include "DetectorDescription/Core/interface/DDMaterial.h"
00021 #include "Geometry/Records/interface/IdealGeometryRecord.h"
00022
00023 #include "SimDataFormats/ValidationFormats/interface/MaterialAccountingStep.h"
00024 #include "SimDataFormats/ValidationFormats/interface/MaterialAccountingTrack.h"
00025 #include "MaterialAccountingGroup.h"
00026 #include "TrackingMaterialAnalyser.h"
00027 #include "TrackingMaterialPlotter.h"
00028
00029
00030 TrackingMaterialAnalyser::TrackingMaterialAnalyser(const edm::ParameterSet& iPSet)
00031 {
00032 m_material = iPSet.getParameter<edm::InputTag>("MaterialAccounting");
00033 m_groupNames = iPSet.getParameter<std::vector<std::string> >("Groups");
00034 const std::string & splitmode = iPSet.getParameter<std::string>("SplitMode");
00035 if (strcasecmp(splitmode.c_str(), "NearestLayer") == 0) {
00036 m_splitMode = NEAREST_LAYER;
00037 } else if (strcasecmp(splitmode.c_str(), "InnerLayer") == 0) {
00038 m_splitMode = INNER_LAYER;
00039 } else if (strcasecmp(splitmode.c_str(), "OuterLayer") == 0) {
00040 m_splitMode = OUTER_LAYER;
00041 } else {
00042 m_splitMode = UNDEFINED;
00043 throw edm::Exception(edm::errors::LogicError) << "Invalid SplitMode \"" << splitmode << "\". Acceptable values are \"NearestLayer\", \"InnerLayer\", \"OuterLayer\".";
00044 }
00045 m_skipAfterLastDetector = iPSet.getParameter<bool>("SkipAfterLastDetector");
00046 m_skipBeforeFirstDetector = iPSet.getParameter<bool>("SkipBeforeFirstDetector");
00047 m_saveSummaryPlot = iPSet.getParameter<bool>("SaveSummaryPlot");
00048 m_saveDetailedPlots = iPSet.getParameter<bool>("SaveDetailedPlots");
00049 m_saveParameters = iPSet.getParameter<bool>("SaveParameters");
00050 m_saveXml = iPSet.getParameter<bool>("SaveXML");
00051 if (m_saveSummaryPlot)
00052 m_plotter = new TrackingMaterialPlotter( 300., 120., 10 );
00053 else
00054 m_plotter = NULL;
00055 }
00056
00057
00058 TrackingMaterialAnalyser::~TrackingMaterialAnalyser(void)
00059 {
00060 if (m_plotter)
00061 delete m_plotter;
00062 }
00063
00064
00065 void TrackingMaterialAnalyser::saveParameters(const char* name)
00066 {
00067 std::ofstream parameters(name);
00068 std::cout << std::endl;
00069 for (unsigned int i = 0; i < m_groups.size(); ++i) {
00070 MaterialAccountingGroup & layer = *(m_groups[i]);
00071 std::cout << layer.name() << std::endl;
00072 std::cout << boost::format("\tnumber of hits: %9d") % layer.tracks() << std::endl;
00073 std::cout << boost::format("\tnormalized segment length: %9.1f ± %9.1f cm") % layer.averageLength() % layer.sigmaLength() << std::endl;
00074 std::cout << boost::format("\tnormalized radiation lengths: %9.3f ± %9.3f") % layer.averageRadiationLengths() % layer.sigmaRadiationLengths() << std::endl;
00075 std::cout << boost::format("\tnormalized energy loss: %9.3f ± %9.3f MeV") % layer.averageEnergyLoss() % layer.sigmaEnergyLoss() << std::endl;
00076 parameters << boost::format("%-20s\t%7d\t%5.1f ± %5.1f cm\t%6.4f ± %6.4f \t%6.4fe-03 ± %6.4fe-03 GeV")
00077 % layer.name()
00078 % layer.tracks()
00079 % layer.averageLength() % layer.sigmaLength()
00080 % layer.averageRadiationLengths() % layer.sigmaRadiationLengths()
00081 % layer.averageEnergyLoss() % layer.sigmaEnergyLoss()
00082 << std::endl;
00083 }
00084 std::cout << std::endl;
00085
00086 parameters.close();
00087 }
00088
00089
00090 void TrackingMaterialAnalyser::saveXml(const char* name)
00091 {
00092 std::ofstream xml(name);
00093 xml << "<?xml version=\"1.0\" encoding=\"utf-8\"?>" << std::endl;
00094 xml << "<Groups>" << std::endl;
00095 for (unsigned int i = 0; i < m_groups.size(); ++i) {
00096 MaterialAccountingGroup & layer = *(m_groups[i]);
00097 xml << " <Group name=\"" << layer.name() << "\">\n"
00098 << " <Parameter name=\"TrackerRadLength\" value=\"" << layer.averageRadiationLengths() << "\"/>\n"
00099 << " <Parameter name=\"TrackerXi\" value=\"" << layer.averageEnergyLoss() << "\"/>\n"
00100 << " </Group>\n"
00101 << std::endl;
00102 }
00103 xml << "</Groudp>" << std::endl;
00104 }
00105
00106
00107 void TrackingMaterialAnalyser::saveLayerPlots()
00108 {
00109 for (unsigned int i = 0; i < m_groups.size(); ++i) {
00110 MaterialAccountingGroup & layer = *(m_groups[i]);
00111 layer.savePlots();
00112 }
00113 }
00114
00115
00116 void TrackingMaterialAnalyser::endJob(void)
00117 {
00118 if (m_saveParameters)
00119 saveParameters("parameters");
00120
00121 if (m_saveXml)
00122 saveXml("parameters.xml");
00123
00124 if (m_saveDetailedPlots)
00125 saveLayerPlots();
00126
00127 if (m_saveSummaryPlot and m_plotter) {
00128 m_plotter->normalize();
00129 m_plotter->draw();
00130 }
00131 }
00132
00133
00134 void TrackingMaterialAnalyser::beginJob(const edm::EventSetup & setup)
00135 {
00136 edm::ESHandle<DDCompactView> hDDD;
00137 setup.get<IdealGeometryRecord>().get( hDDD );
00138
00139 m_groups.reserve( m_groupNames.size() );
00140 for (unsigned int i = 0; i < m_groupNames.size(); ++i)
00141 m_groups.push_back( new MaterialAccountingGroup( m_groupNames[i], * hDDD) );
00142
00143
00144 std::cout << "TrackingMaterialAnalyser: List of the tracker groups: " << std::endl;
00145 for (unsigned int i = 0; i < m_groups.size(); ++i)
00146 std::cout << '\t' << m_groups[i]->info() << std::endl;
00147 std::cout << std::endl;
00148 }
00149
00150
00151 void TrackingMaterialAnalyser::analyze(const edm::Event& event, const edm::EventSetup& setup)
00152 {
00153 edm::Handle< std::vector<MaterialAccountingTrack> > h_tracks;
00154 event.getByLabel(m_material, h_tracks);
00155
00156 for (std::vector<MaterialAccountingTrack>::const_iterator t = h_tracks->begin(), end = h_tracks->end(); t != end; ++t) {
00157 MaterialAccountingTrack track(*t);
00158 split( track );
00159 }
00160 }
00161
00162
00163
00164
00165
00166
00167
00168
00169 void TrackingMaterialAnalyser::split( MaterialAccountingTrack & track )
00170 {
00171
00172 std::vector<int> group( track.m_detectors.size() );
00173 for (unsigned int i = 0; i < track.m_detectors.size(); ++i)
00174 group[i] = findLayer( track.m_detectors[i] );
00175
00176 unsigned int detectors = track.m_detectors.size();
00177 if (detectors == 0) {
00178
00179
00180 if (m_plotter)
00181 for (unsigned int i = 1; i < track.m_steps.size(); ++i)
00182 m_plotter->plotSegmentUnassigned( track.m_steps[i] );
00183 } else {
00184 const double TOLERANCE = 0.0001;
00185 std::vector<double> limits(detectors + 2);
00186
00187
00188 if (m_skipBeforeFirstDetector)
00189 limits[0] = track.m_detectors[0].m_curvilinearIn - TOLERANCE;
00190 else
00191 limits[0] = - TOLERANCE;
00192 if (m_skipAfterLastDetector)
00193 limits[detectors] = track.m_detectors[detectors-1].m_curvilinearOut + TOLERANCE;
00194 else
00195 limits[detectors] = track.m_total.length() + TOLERANCE;
00196 limits[detectors+1] = INFINITY;
00197
00198
00199 switch (m_splitMode) {
00200
00201
00202 case NEAREST_LAYER:
00203 for (unsigned int i = 1; i < detectors; ++i)
00204 limits[i] = (track.m_detectors[i-1].m_curvilinearOut + track.m_detectors[i].m_curvilinearIn) / 2.;
00205 break;
00206
00207
00208
00209 case INNER_LAYER:
00210 for (unsigned int i = 1; i < detectors; ++i)
00211 limits[i] = track.m_detectors[i].m_curvilinearIn - TOLERANCE;
00212 break;
00213
00214
00215
00216 case OUTER_LAYER:
00217 for (unsigned int i = 1; i < detectors; ++i)
00218 limits[i] = track.m_detectors[i-1].m_curvilinearOut + TOLERANCE;
00219 break;
00220
00221 case UNDEFINED:
00222 default:
00223
00224 throw edm::Exception(edm::errors::LogicError) << "Invalid SplitMode";
00225 }
00226
00227
00228
00229
00230 double begin = 0.;
00231 double end = 0.;
00232 unsigned int i = 1;
00233
00234
00235
00236 while (end < limits[0]) {
00237 const MaterialAccountingStep & step = track.m_steps[i++];
00238 end = begin + step.length();
00239
00240
00241 if (m_plotter)
00242 m_plotter->plotSegmentUnassigned( step );
00243
00244 begin = end;
00245
00246 }
00247
00248
00249
00250
00251 if (begin < limits[0] and end > limits[0]) {
00252 const MaterialAccountingStep & step = track.m_steps[i++];
00253 end = begin + step.length();
00254
00255 double fraction = (limits[0] - begin) / (end - begin);
00256 std::pair<MaterialAccountingStep, MaterialAccountingStep> parts = step.split(fraction);
00257
00258
00259 track.m_detectors[0].account( parts.second, limits[1], end );
00260
00261 if (m_plotter) {
00262
00263 m_plotter->plotSegmentUnassigned( parts.first );
00264
00265
00266 m_plotter->plotSegmentInLayer( parts.second, group[0] );
00267 }
00268 begin = end;
00269 }
00270
00271 unsigned int index = 0;
00272 while (i < track.m_steps.size()) {
00273 const MaterialAccountingStep & step = track.m_steps[i++];
00274
00275 end = begin + step.length();
00276
00277 if (begin > limits[detectors]) {
00278
00279 if (m_plotter)
00280 m_plotter->plotSegmentUnassigned( step );
00281 begin = end;
00282 continue;
00283 }
00284
00285
00286
00287
00288
00289 if (begin < limits[index] or end > limits[index+2]) {
00290
00291 std::cerr << "MaterialAccountingTrack::split(): ERROR: internal logic error, expected " << limits[index] << " < " << begin << " < " << limits[index+1] << std::endl;
00292 break;
00293 }
00294
00295
00296 if (limits[index] <= begin and end <= limits[index+1]) {
00297
00298 track.m_detectors[index].account( step, begin, end );
00299 if (m_plotter)
00300 m_plotter->plotSegmentInLayer( step, group[index] );
00301 } else {
00302
00303 double fraction = (limits[index+1] - begin) / (end - begin);
00304 std::pair<MaterialAccountingStep, MaterialAccountingStep> parts = step.split(fraction);
00305
00306 if (m_plotter) {
00307 if (index > 0)
00308 m_plotter->plotSegmentInLayer( parts.first, group[index] );
00309 else
00310
00311 m_plotter->plotSegmentUnassigned( parts.first );
00312
00313 if (index+1 < detectors)
00314 m_plotter->plotSegmentInLayer( parts.second, group[index+1] );
00315 else
00316
00317 m_plotter->plotSegmentUnassigned( parts.second );
00318 }
00319
00320 track.m_detectors[index].account( parts.first, begin, limits[index+1] );
00321 ++index;
00322
00323
00324 if (index < detectors)
00325 track.m_detectors[index].account( parts.second, limits[index+1], end );
00326 }
00327 begin = end;
00328 }
00329
00330 }
00331
00332
00333
00334 for (unsigned int i = 0; i < track.m_detectors.size(); ++i)
00335 if (group[i] != 0)
00336 m_groups[group[i]-1]->addDetector( track.m_detectors[i] );
00337
00338
00339 for (unsigned int i = 0; i < m_groups.size(); ++i)
00340 m_groups[i]->endOfTrack();
00341 }
00342
00343
00344
00345 int TrackingMaterialAnalyser::findLayer( const MaterialAccountingDetector & detector )
00346 {
00347 int index = 0;
00348 size_t inside = 0;
00349 for (size_t i = 0; i < m_groups.size(); ++i)
00350 if (m_groups[i]->inside(detector)) {
00351 ++inside;
00352 index = i+1;
00353 }
00354 if (inside == 0) {
00355 index = 0;
00356 std::cerr << "TrackingMaterialAnalyser::findLayer(...): ERROR: detector does not belong to any DetLayer" << std::endl;
00357 std::cerr << "TrackingMaterialAnalyser::findLayer(...): detector position: " << std::fixed
00358 << " (r: " << std::setprecision(1) << std::setw(5) << detector.position().perp()
00359 << ", z: " << std::setprecision(1) << std::setw(6) << detector.position().z()
00360 << ", phi: " << std::setprecision(3) << std::setw(6) << detector.position().phi() << ")"
00361 << std::endl;
00362 }
00363 if (inside > 1) {
00364 index = 0;
00365 std::cerr << "TrackingMaterialAnalyser::findLayer(...): ERROR: detector belongs to " << inside << "DetLayers" << std::endl;
00366 std::cerr << "TrackingMaterialAnalyser::findLayer(...): detector position: " << std::fixed
00367 << " (r: " << std::setprecision(1) << std::setw(5) << detector.position().perp()
00368 << ", z: " << std::setprecision(1) << std::setw(6) << detector.position().z()
00369 << ", phi: " << std::setprecision(3) << std::setw(6) << detector.position().phi() << ")"
00370 << std::endl;
00371 }
00372
00373 return index;
00374 }
00375
00376
00377
00378 #include "FWCore/PluginManager/interface/ModuleDef.h"
00379 #include "FWCore/Framework/interface/MakerMacros.h"
00380 DEFINE_FWK_MODULE(TrackingMaterialAnalyser);