ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/FastOpenGlDisplayer/src/Frog_Analyzer.cc
Revision: 1.9
Committed: Sat May 10 14:14:03 2008 UTC (16 years, 11 months ago) by querten
Content type: text/plain
Branch: MAIN
Changes since 1.8: +47 -11 lines
Log Message:
Allow Both RecoTrack and Trajectory in the .vis

File Contents

# User Rev Content
1 querten 1.1 // -*- C++ -*-
2     //
3     // Package: Frog_Analyzer
4     // Class: Frog_Analyzer
5     //
6     /**\class Frog_Analyzer Frog_Analyzer.cc Visualisation/Frog/src/Frog_Analyzer.cc
7    
8     Description: <one line class summary>
9    
10     Implementation:
11     <Notes on implementation>
12     */
13     //
14     // Original Author: Loic QUERTENMONT
15     // Created: Fri Oct 26 07:22:12 CEST 2007
16 querten 1.9 // $Id: Frog_Analyzer.cc,v 1.8 2008/05/05 08:49:58 querten Exp $
17 querten 1.1 //
18     //
19    
20    
21     // system include files
22     #include <memory>
23    
24     // user include files
25     #include "FWCore/Framework/interface/Frameworkfwd.h"
26     #include "FWCore/Framework/interface/EDAnalyzer.h"
27    
28     #include "FWCore/Framework/interface/Event.h"
29     #include "FWCore/Framework/interface/MakerMacros.h"
30    
31     #include "FWCore/ParameterSet/interface/ParameterSet.h"
32     #include "FWCore/ServiceRegistry/interface/Service.h"
33     #include "FWCore/Framework/interface/ESHandle.h"
34    
35     #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
36     #include "Geometry/CommonDetUnit/interface/GeomDetUnit.h"
37     #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
38     #include "Geometry/TrackerNumberingBuilder/interface/GeometricDet.h"
39     #include "Geometry/CommonTopologies/interface/PixelTopology.h"
40     #include "Geometry/CommonTopologies/interface/StripTopology.h"
41     #include "Geometry/TrackerGeometryBuilder/interface/PixelGeomDetType.h"
42     #include "Geometry/TrackerGeometryBuilder/interface/StripGeomDetType.h"
43     #include "Geometry/TrackerGeometryBuilder/interface/PixelGeomDetUnit.h"
44     #include "Geometry/TrackerGeometryBuilder/interface/StripGeomDetUnit.h"
45     #include "DataFormats/GeometrySurface/interface/BoundSurface.h"
46     #include "DataFormats/DetId/interface/DetId.h"
47    
48     #include "Geometry/CaloEventSetup/interface/CaloTopologyRecord.h"
49     #include "Geometry/Records/interface/IdealGeometryRecord.h"
50     #include "Geometry/CaloGeometry/interface/CaloSubdetectorGeometry.h"
51     #include "Geometry/CaloGeometry/interface/CaloGeometry.h"
52     #include "Geometry/EcalBarrelAlgo/interface/EcalBarrelGeometry.h"
53     #include "Geometry/EcalEndcapAlgo/interface/EcalEndcapGeometry.h"
54     #include "Geometry/CaloGeometry/interface/CaloCellGeometry.h"
55     #include "DataFormats/EcalDetId/interface/EcalSubdetector.h"
56    
57     #include "Geometry/CommonDetUnit/interface/TrackingGeometry.h"
58     #include "Geometry/Records/interface/MuonGeometryRecord.h"
59     #include "Geometry/DTGeometry/interface/DTGeometry.h"
60     #include "Geometry/CSCGeometry/interface/CSCGeometry.h"
61     #include "Geometry/RPCGeometry/interface/RPCGeometry.h"
62    
63    
64     #include "DataFormats/CaloRecHit/interface/CaloRecHit.h"
65     #include "DataFormats/EcalRecHit/interface/EcalRecHitCollections.h"
66     #include "DataFormats/HcalRecHit/interface/HcalRecHitCollections.h"
67    
68     #include "SimDataFormats/TrackingHit/interface/PSimHitContainer.h"
69     #include "SimDataFormats/TrackingHit/interface/PSimHit.h"
70     #include "SimDataFormats/CrossingFrame/interface/CrossingFrame.h"
71     #include "SimDataFormats/CrossingFrame/interface/MixCollection.h"
72     #include "SimDataFormats/CaloHit/interface/PCaloHit.h"
73     #include "SimDataFormats/Track/interface/SimTrack.h"
74     #include "SimDataFormats/Vertex/interface/SimVertex.h"
75    
76     #include "DataFormats/Provenance/interface/BranchDescription.h"
77     #include "DataFormats/Provenance/interface/Provenance.h"
78     #include "DataFormats/Candidate/interface/Candidate.h"
79     #include "DataFormats/HepMCCandidate/interface/GenParticleCandidate.h"
80     #include "DataFormats/RecoCandidate/interface/RecoChargedCandidate.h"
81     #include "DataFormats/TrackReco/interface/Track.h"
82 roberfro 1.5 #include "TrackingTools/PatternTools/interface/Trajectory.h"
83 querten 1.6 #include "TrackingTools/PatternTools/interface/TrajTrackAssociation.h"
84    
85 querten 1.1
86 querten 1.2 #include "DataFormats/DTRecHit/interface/DTRecSegment4DCollection.h"
87     #include "DataFormats/CSCRecHit/interface/CSCSegmentCollection.h"
88    
89     #include "DataFormats/DTRecHit/interface/DTRecHitCollection.h"
90     #include "DataFormats/CSCRecHit/interface/CSCRecHit2DCollection.h"
91     #include "DataFormats/RPCRecHit/interface/RPCRecHitCollection.h"
92    
93 querten 1.1 #include "DataFormats/Math/interface/Point3D.h"
94     #include "DataFormats/GeometrySurface/interface/TrapezoidalPlaneBounds.h"
95     #include "DataFormats/GeometrySurface/interface/RectangularPlaneBounds.h"
96    
97     #include "Visualisation/Frog/interface/SimEvent.h"
98     #include "Visualisation/Frog/interface/Geometry.h"
99    
100     using namespace edm;
101     using namespace std;
102    
103     //
104     // class decleration
105     //
106    
107     class Frog_Analyzer : public edm::EDAnalyzer {
108     public:
109     explicit Frog_Analyzer(const edm::ParameterSet&);
110     ~Frog_Analyzer();
111    
112    
113     private:
114     virtual void beginJob(const edm::EventSetup& iSetup);
115     virtual void analyze (const edm::Event&, const edm::EventSetup&);
116     virtual void endJob ();
117    
118    
119     MySimEvents* MyEvents;
120     Geometry* Geom_Tracker;
121     Geometry* Geom_ECAL;
122     Geometry* Geom_HCAL;
123     Geometry* Geom_Muon;
124 querten 1.7 Geometry* Geom_Fwd;
125 querten 1.1
126    
127     std::string OutputFile;
128    
129     std::vector<InputTag> SimTrackProducers;
130     std::vector<InputTag> SimVertexProducers;
131     std::vector<InputTag> SimHitProducers;
132    
133 querten 1.9 std::vector<InputTag> TrackProducers;
134 roberfro 1.5 std::vector<InputTag> TrajectoryProducers;
135 querten 1.1 std::vector<InputTag> EcalRecHitProducers;
136     std::vector<InputTag> HcalHBHERecHitProducers;
137     std::vector<InputTag> HcalHORecHitProducers;
138     std::vector<InputTag> HcalHFRecHitProducers;
139 querten 1.2 std::vector<InputTag> DTSegmentProducers;
140     std::vector<InputTag> CSCSegmentProducers;
141    
142     std::vector<InputTag> RPCHitsProducers;
143 querten 1.1
144    
145    
146    
147     // ----------member data ---------------------------
148     };
149    
150     //
151     // constructors and destructor
152     //
153     Frog_Analyzer::Frog_Analyzer(const edm::ParameterSet& iConfig)
154     {
155     OutputFile = iConfig.getParameter<std::string >("OutputFile");
156    
157     SimTrackProducers = iConfig.getParameter<std::vector<InputTag> >("SimTrackProducers");
158     SimVertexProducers = iConfig.getParameter<std::vector<InputTag> >("SimVertexProducers");
159     SimHitProducers = iConfig.getParameter<std::vector<InputTag> >("SimHitProducers");
160    
161 querten 1.9 TrackProducers = iConfig.getParameter<std::vector<InputTag> >("TrajectoryProducers");
162     TrajectoryProducers = iConfig.getParameter<std::vector<InputTag> >("TrajectoryProducers");
163 querten 1.1 EcalRecHitProducers = iConfig.getParameter<std::vector<InputTag> >("EcalRecHitProducers");
164     HcalHBHERecHitProducers = iConfig.getParameter<std::vector<InputTag> >("HcalHBHERecHitProducers");
165     HcalHORecHitProducers = iConfig.getParameter<std::vector<InputTag> >("HcalHORecHitProducers");
166     HcalHFRecHitProducers = iConfig.getParameter<std::vector<InputTag> >("HcalHFRecHitProducers");
167    
168 querten 1.2 DTSegmentProducers = iConfig.getParameter<std::vector<InputTag> >("DTSegmentProducers");
169     CSCSegmentProducers = iConfig.getParameter<std::vector<InputTag> >("CSCSegmentProducers");
170    
171     RPCHitsProducers = iConfig.getParameter<std::vector<InputTag> >("RPCHitsProducers");
172    
173 querten 1.1 }
174    
175    
176     Frog_Analyzer::~Frog_Analyzer()
177     {
178     }
179    
180     // ------------ method called once each job just before starting event loop ------------
181     void
182     Frog_Analyzer::beginJob(const edm::EventSetup& iSetup)
183     {
184     DetId Detid;
185     int SubDet;
186    
187     MyEvents = new MySimEvents();
188     Geom_Tracker = new Geometry();
189     Geom_ECAL = new Geometry();
190     Geom_HCAL = new Geometry();
191     Geom_Muon = new Geometry();
192 querten 1.7 Geom_Fwd = new Geometry();
193    
194 querten 1.1
195    
196     // ### TRACKER GEOMETRY ###
197    
198     edm::ESHandle<TrackerGeometry> tkGeom;
199     iSetup.get<TrackerDigiGeometryRecord>().get( tkGeom );
200     vector<GeomDet*> TrackerDets = tkGeom->dets();
201    
202     for(unsigned int i=0;i<TrackerDets.size();i++){
203     Detid = TrackerDets[i]->geographicalId();
204     // SubDet = Detid.subdetId();
205    
206     GeomDet* DetUnit = TrackerDets[i];
207     if(!DetUnit)continue;
208     const BoundPlane plane = DetUnit->surface();
209     const TrapezoidalPlaneBounds* trapezoidalBounds( dynamic_cast<const TrapezoidalPlaneBounds*>(&(plane.bounds())));
210     const RectangularPlaneBounds* rectangularBounds( dynamic_cast<const RectangularPlaneBounds*>(&(plane.bounds())));
211    
212     float width = 0;
213     float length = 0;
214     float thickness = 0;
215     float TrapezoidalParam = 0;
216    
217     if(trapezoidalBounds)
218     {
219     std::vector<float> const & parameters = (*trapezoidalBounds).parameters();
220     width = parameters[0]*2;
221     length = parameters[3]*2;
222     thickness = (*trapezoidalBounds).thickness();
223     TrapezoidalParam = parameters[1]/parameters[0];
224     }else if(rectangularBounds){
225     width = DetUnit->surface().bounds().width();
226     length = DetUnit->surface().bounds().length();
227     thickness = DetUnit->surface().bounds().thickness();
228     TrapezoidalParam = 1;
229     }
230    
231     Surface::GlobalPoint WidthVector = plane.toGlobal( LocalPoint(width/2, 0, 0) );
232     Surface::GlobalPoint LengthVector = plane.toGlobal( LocalPoint(0, length/2, 0) );
233     Surface::GlobalPoint ThickVector = plane.toGlobal( LocalPoint(0, 0, thickness/2) );
234    
235     GlobalVector Pos = GlobalVector(DetUnit->position().basicVector());
236    
237     Geom_Tracker->Add_TrackerDet(Detid.rawId(), TrapezoidalParam,
238     Pos.x(), Pos.y(), Pos.z(),
239     WidthVector.x() -Pos.x(), WidthVector.y() -Pos.y(), WidthVector.z() -Pos.z(),
240     LengthVector.x()-Pos.x(), LengthVector.y()-Pos.y(), LengthVector.z()-Pos.z(),
241     ThickVector.x() -Pos.x(), ThickVector.y() -Pos.y(), ThickVector.z() -Pos.z());
242    
243    
244     }
245    
246    
247     // ### CALO GEOMETRY ###
248    
249     edm::ESHandle<CaloGeometry> CaloGeom;
250     iSetup.get<IdealGeometryRecord>().get( CaloGeom );
251     const vector<DetId> CaloDets = CaloGeom->getValidDetIds();
252    
253     for(unsigned int i=0;i<CaloDets.size();i++)
254     {
255     Detid = CaloDets[i];
256     SubDet = Detid.subdetId();
257    
258     Geometry* Geom_temp = NULL;
259     if(Detid.det()==DetId::Ecal){
260     Geom_temp = Geom_ECAL;
261    
262     if(SubDet<1 || SubDet>3){
263     printf("Don't Save EcalTriggerTower or EcalLaserPnDiode\n");
264     continue;
265     }
266     }else if(Detid.det()==DetId::Hcal){
267     Geom_temp = Geom_HCAL;
268    
269     if(SubDet<1 || SubDet>4){
270     printf("Don't Save HcalEmpty, HcalTriggerTower, HcalOther\n");
271     continue;
272     }
273 querten 1.7 }else if(Detid.det() == DetId::Calo){
274     Geom_temp = NULL;
275    
276     // Calo Tower Geometry (SudDet=1)
277     // if(SubDet==1)
278     // //printf("Don't Save CaloTowers\n");
279     // continue;
280    
281     // Forward Calo Geometry --> Castor (SudDet=3) and ZDS (SubDet=2)
282     if(SubDet==HcalCastorDetId::SubdetectorId || SubDet==HcalZDCDetId::SubdetectorId)
283     Geom_temp = Geom_Fwd;
284     // printf("Save FWD CaloTowers\n");
285 querten 1.1 }
286    
287     const CaloCellGeometry* CellGeom = CaloGeom->getGeometry(Detid);
288     GlobalPoint CellPos = CellGeom->getPosition();
289     const CaloCellGeometry::CornersVec CellCorners = CellGeom->getCorners();
290    
291 querten 1.7 if(Geom_temp == Geom_ECAL || Geom_temp == Geom_HCAL){
292     float cX = (CellCorners[0].x() + CellCorners[2].x())/2;
293     float cY = (CellCorners[0].y() + CellCorners[2].y())/2;
294     float cZ = (CellCorners[0].z() + CellCorners[2].z())/2;
295    
296     float wX = (CellCorners[1].x() - CellCorners[0].x())/2;
297     float wY = (CellCorners[1].y() - CellCorners[0].y())/2;
298     float wZ = (CellCorners[1].z() - CellCorners[0].z())/2;
299    
300     float hX = (CellCorners[3].x() - CellCorners[0].x())/2;
301     float hY = (CellCorners[3].y() - CellCorners[0].y())/2;
302     float hZ = (CellCorners[3].z() - CellCorners[0].z())/2;
303    
304     float F = sqrt( pow(CellCorners[4].x()+CellCorners[6].x(),2) + pow(CellCorners[4].y()+CellCorners[6].y(),2) + pow(CellCorners[4].z()+CellCorners[6].z(),2) );
305     F /= sqrt( pow(CellCorners[0].x()+CellCorners[2].x(),2) + pow(CellCorners[0].y()+CellCorners[2].y(),2) + pow(CellCorners[0].z()+CellCorners[2].z(),2) );
306    
307     Geom_temp->Add_CaloDet(Detid.rawId(),
308     cX , cY , cZ ,
309     wX , wY , wZ ,
310     hX , hY , hZ , F );
311     }else if(Geom_temp == Geom_Fwd){
312     /*
313     Geom_temp->Add_FwdDet(Detid.rawId(),
314     CellCorners[0].x(), CellCorners[0].y(), CellCorners[0].z(),
315     CellCorners[1].x(), CellCorners[1].y(), CellCorners[1].z(),
316     CellCorners[2].x(), CellCorners[2].y(), CellCorners[2].z(),
317     CellCorners[3].x(), CellCorners[3].y(), CellCorners[3].z(),
318     CellCorners[4].x(), CellCorners[4].y(), CellCorners[4].z(),
319     CellCorners[5].x(), CellCorners[5].y(), CellCorners[5].z(),
320     CellCorners[6].x(), CellCorners[6].y(), CellCorners[6].z(),
321     CellCorners[7].x(), CellCorners[7].y(), CellCorners[7].z());
322     printf("SubDet = %i : (%6.2f,%6.2f,%6.2f) --> 7 (%6.2f,%6.2f,%6.2f) \n", SubDet,CellCorners[0].x(), CellCorners[0].y(), CellCorners[0].z(),CellCorners[6].x(), CellCorners[6].y(), CellCorners[6].z());
323     */
324     }
325 querten 1.1 }
326 querten 1.7
327 querten 1.1
328     // ### MUON GEOMETRY ###
329    
330     edm::ESHandle<DTGeometry> DtGeom;
331     iSetup.get<MuonGeometryRecord>().get( DtGeom );
332     const vector<GeomDet*> DtDets = DtGeom->dets();
333    
334     edm::ESHandle<CSCGeometry> CscGeom;
335     iSetup.get<MuonGeometryRecord>().get( CscGeom );
336     const vector<GeomDet*> CscDets = CscGeom->dets();
337    
338     edm::ESHandle<RPCGeometry> RpcGeom;
339     iSetup.get<MuonGeometryRecord>().get( RpcGeom );
340     const vector<GeomDet*> RpcDets = RpcGeom->dets();
341    
342     vector<GeomDet*> MuonDets;
343     for(unsigned int i=0;i<DtDets.size() ;i++){MuonDets.push_back(DtDets [i]);}
344     for(unsigned int i=0;i<CscDets.size();i++){MuonDets.push_back(CscDets[i]);}
345     for(unsigned int i=0;i<RpcDets.size();i++){MuonDets.push_back(RpcDets[i]);}
346    
347    
348     for(unsigned int i=0;i<MuonDets.size();i++)
349     {
350     Detid = DetId(MuonDets[i]->geographicalId());
351     SubDet = Detid.subdetId();
352    
353     GeomDet* DetUnit = MuonDets[i];
354     if(!DetUnit)continue;
355     const BoundPlane plane = DetUnit->surface();
356     const TrapezoidalPlaneBounds* trapezoidalBounds( dynamic_cast<const TrapezoidalPlaneBounds*>(&(plane.bounds())));
357     const RectangularPlaneBounds* rectangularBounds( dynamic_cast<const RectangularPlaneBounds*>(&(plane.bounds())));
358    
359     float width = 0;
360     float length = 0;
361     float thickness = 0;
362     float TrapezoidalParam = 0;
363    
364     if(trapezoidalBounds)
365     {
366     std::vector<float> const & parameters = (*trapezoidalBounds).parameters();
367     width = parameters[0]*2;
368     length = parameters[3]*2;
369     thickness = (*trapezoidalBounds).thickness();
370     TrapezoidalParam = parameters[1]/parameters[0];
371     }else if(rectangularBounds){
372     width = DetUnit->surface().bounds().width();
373     length = DetUnit->surface().bounds().length();
374     thickness = DetUnit->surface().bounds().thickness();
375     TrapezoidalParam = 1;
376     }
377    
378     Surface::GlobalPoint WidthVector = plane.toGlobal( LocalPoint(width/2, 0, 0) );
379     Surface::GlobalPoint LengthVector = plane.toGlobal( LocalPoint(0, length/2, 0) );
380     Surface::GlobalPoint ThickVector = plane.toGlobal( LocalPoint(0, 0, thickness/2) );
381    
382     GlobalVector Pos = GlobalVector(DetUnit->position().basicVector());
383    
384     Geom_Muon->Add_TrackerDet(Detid.rawId(), TrapezoidalParam,
385     Pos.x(), Pos.y(), Pos.z(),
386     WidthVector.x() -Pos.x(), WidthVector.y() -Pos.y(), WidthVector.z() -Pos.z(),
387     LengthVector.x()-Pos.x(), LengthVector.y()-Pos.y(), LengthVector.z()-Pos.z(),
388     ThickVector.x() -Pos.x(), ThickVector.y() -Pos.y(), ThickVector.z() -Pos.z());
389     }
390    
391    
392     // ### Save .geom ###
393    
394 querten 1.9 // printf("FWD Size = %i\n", Geom_Fwd->Det_ECAL_ALL.size());
395     // printf("FWD Size = %i\n", Geom_Fwd->Det_HCAL_ALL.size());
396     // printf("FWD Size = %i\n", Geom_Fwd->Det_FWD_ALL.size());
397 querten 1.7
398    
399    
400 querten 1.1 Geom_Tracker->Save("Tracker.geom");
401     Geom_ECAL ->Save("Ecal.geom");
402     Geom_HCAL ->Save("Hcal.geom");
403     Geom_Muon ->Save("Muon.geom");
404 querten 1.7 // Geom_Fwd ->Save("Fwd.geom");
405    
406    
407 querten 1.1
408     }
409    
410     // ------------ method called once each job just after ending the event loop ------------
411     void
412     Frog_Analyzer::endJob() {
413    
414     MyEvents->Save((char*) OutputFile.c_str());
415     // MyEvents->Load((char*) OutputFile.c_str());
416     }
417    
418    
419    
420     //
421     // member functions
422     //
423    
424     // ------------ method called to for each event ------------
425     void
426     Frog_Analyzer::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup)
427     {
428     MySimEvent* MyEvent = new MySimEvent;
429    
430 querten 1.2 // access the tracker
431     edm::ESHandle<TrackerGeometry> tkGeom;
432     iSetup.get<TrackerDigiGeometryRecord>().get( tkGeom );
433 querten 1.1
434 querten 1.2 edm::ESHandle<DTGeometry> DtGeom;
435     iSetup.get<MuonGeometryRecord>().get( DtGeom );
436    
437     edm::ESHandle<CSCGeometry> CscGeom;
438     iSetup.get<MuonGeometryRecord>().get( CscGeom );
439 querten 1.1
440 querten 1.2 edm::ESHandle<RPCGeometry> RpcGeom;
441     iSetup.get<MuonGeometryRecord>().get( RpcGeom );
442 querten 1.1
443    
444     // Save Data
445 querten 1.2 // ### SIMTRACK ###
446 querten 1.1 for(unsigned int i=0;i<SimTrackProducers.size();i++){
447     edm::Handle<std::vector< SimTrack > > h_SimTracks;
448     iEvent.getByLabel(SimTrackProducers[i], h_SimTracks);
449     std::vector< SimTrack > SimTrackColl = *h_SimTracks.product();
450    
451     for ( unsigned int a = 0; a < SimTrackColl.size(); ++a ) {
452 querten 1.4 MySimTrack* MysimTrack = new MySimTrack;
453 querten 1.1 SimTrack simTrack =SimTrackColl[a];
454    
455 querten 1.4 MysimTrack->track_id =simTrack.trackId();
456     MysimTrack->Type =simTrack.type();
457     MysimTrack->parent_vertex=simTrack.vertIndex();
458     MysimTrack->Px =simTrack.momentum().x();
459     MysimTrack->Py =simTrack.momentum().y();
460     MysimTrack->Pz =simTrack.momentum().z();
461     MysimTrack->E =simTrack.momentum().e();
462     MysimTrack->charge =simTrack.charge();
463 querten 1.1
464     MyEvent->MySimTrackCollection.push_back(MysimTrack);
465     }
466     }
467    
468 querten 1.2 // ### SIMVERTEX ###
469 querten 1.1 for(unsigned int i=0;i<SimVertexProducers.size();i++){
470     edm::Handle<std::vector< SimVertex > > h_Vertex;
471     iEvent.getByLabel(SimVertexProducers[i], h_Vertex);
472     std::vector< SimVertex > VertexColl = *h_Vertex.product();
473    
474     for (unsigned int b = 0; b < VertexColl.size(); ++b ) {
475 querten 1.4 MySimVertex* MyVertex = new MySimVertex;
476 querten 1.1 SimVertex Vertex = VertexColl[b];
477    
478 querten 1.4 MyVertex->parentTrack_id=Vertex.parentIndex ();
479     MyVertex->x =Vertex.position().x();
480     MyVertex->y =Vertex.position().y();
481     MyVertex->z =Vertex.position().z();
482 querten 1.1
483     MyEvent->MySimVertexCollection.push_back(MyVertex);
484     }
485     }
486    
487 querten 1.2 // ### SIMHIT ###
488 querten 1.1 for(unsigned int i=0;i<SimHitProducers.size();i++){
489     edm::Handle<std::vector< PSimHit > > h_Hits;
490     iEvent.getByLabel(SimHitProducers[i], h_Hits);
491     std::vector< PSimHit > Hits = *h_Hits.product();
492    
493     for(unsigned int h=0; h<Hits.size(); h++)
494     {
495 querten 1.2 DetId theDetUnitId(Hits[h].detUnitId());
496     const GeomDet* theDet = tkGeom->idToDet(theDetUnitId);
497 querten 1.1
498 querten 1.4 MyPSimHit* Hit = new MyPSimHit;
499     Hit->x = theDet->surface().toGlobal(Hits[h].localPosition()).x();
500     Hit->y = theDet->surface().toGlobal(Hits[h].localPosition()).y();
501     Hit->z = theDet->surface().toGlobal(Hits[h].localPosition()).z();
502     Hit->ProcessType = Hits[h].processType();
503     Hit->dEdX = Hits[h].energyLoss();
504 querten 1.1
505     MyEvent->MyPSimHitCollection.push_back(Hit);
506     }
507     }
508    
509    
510 querten 1.9 // ### TRAJECTORY ###
511 roberfro 1.5 for(unsigned int i=0;i<TrajectoryProducers.size();i++){
512 querten 1.6 Handle<TrajTrackAssociationCollection> trajTrackAssociationHandle;
513     iEvent.getByLabel(TrajectoryProducers[i], trajTrackAssociationHandle);
514     const TrajTrackAssociationCollection TrajToTrackMap = *trajTrackAssociationHandle.product();
515    
516     for(TrajTrackAssociationCollection::const_iterator it = TrajToTrackMap.begin(); it!=TrajToTrackMap.end(); ++it) {
517     Trajectory recoTraj = *it->key;
518     MyRecoTrack* MyrecoTrack = new MyRecoTrack;
519    
520     MyRecoTrackInfo* info = new MyRecoTrackInfo;
521     info->Coll = i;
522     info->P = it->val->p();
523     info->Pt = it->val->pt();
524 querten 1.8 info->chi2 = it->key->chiSquared();
525 querten 1.6 MyrecoTrack->Info = info;
526    
527 querten 1.1
528 roberfro 1.5 std::vector<TrajectoryMeasurement> measurements = recoTraj.measurements();
529     for(unsigned int h=0;h<measurements.size();h++){
530    
531     TrajectoryMeasurement::ConstRecHitPointer h_it = measurements[h].recHit();
532 querten 1.1 if(!h_it->isValid() )continue;
533     DetId detId = h_it->geographicalId();
534 querten 1.6 // const GeomDet * theDet = tkGeom->idToDet(detId);
535     // LocalPoint localPos = h_it->localPosition();
536 roberfro 1.5 GlobalPoint globalPos = measurements[h].updatedState().globalPosition();
537 querten 1.1
538 querten 1.4 MyRecoHit* hit = new MyRecoHit;
539 roberfro 1.5 /*
540 querten 1.4 hit->x = theDet->surface().toGlobal(localPos).x();
541     hit->y = theDet->surface().toGlobal(localPos).y();
542     hit->z = theDet->surface().toGlobal(localPos).z();
543 roberfro 1.5 */
544     hit->x = globalPos.x();
545     hit->y = globalPos.y();
546     hit->z = globalPos.z();
547 querten 1.4 hit->DetId = detId.rawId();
548     hit->Charge = -1;
549     MyrecoTrack->Hits.push_back(hit);
550 querten 1.1 }
551     MyEvent->MyRecoTrackCollection.push_back(MyrecoTrack);
552     }
553     }
554    
555 querten 1.9 // ### TRACK ###
556    
557     for(unsigned int i=0;i<TrackProducers.size();i++){
558     edm::Handle<std::vector< reco::Track > > h_Tracks;
559     iEvent.getByLabel(TrackProducers[i], h_Tracks);
560     std::vector< reco::Track > TrackColl = *h_Tracks.product();
561    
562     for ( unsigned int t = 0; t < TrackColl.size(); ++t ) {
563     reco::Track recoTrack =TrackColl[t];
564    
565     MyRecoTrack* MyrecoTrack = new MyRecoTrack;
566     MyRecoTrackInfo* info = new MyRecoTrackInfo;
567     info->Coll = TrajectoryProducers.size()+i;
568     info->P = recoTrack.p();
569     info->Pt = recoTrack.pt();
570     info->chi2 = recoTrack.chi2();
571     MyrecoTrack->Info = info;
572    
573    
574     for(unsigned int h=0;h<recoTrack.recHitsSize();h++){
575     TrackingRecHitRef h_it = recoTrack.recHit(h);
576     if(!h_it->isValid() )continue;
577     DetId detId = h_it->geographicalId();
578     const GeomDet * theDet = tkGeom->idToDet(detId);
579     LocalPoint localPos = h_it->localPosition();
580    
581     MyRecoHit* hit = new MyRecoHit;
582     hit->x = theDet->surface().toGlobal(localPos).x();
583     hit->y = theDet->surface().toGlobal(localPos).y();
584     hit->z = theDet->surface().toGlobal(localPos).z();
585     hit->DetId = detId.rawId();
586     hit->Charge = -1;
587     MyrecoTrack->Hits.push_back(hit);
588     }
589     MyEvent->MyRecoTrackCollection.push_back(MyrecoTrack);
590     }
591     }
592    
593    
594 querten 1.1
595 querten 1.2 // ### ECALRecHits ###
596 querten 1.1 for(unsigned int i=0;i<EcalRecHitProducers.size();i++){
597     edm::Handle<EcalRecHitCollection > h_Ecal_RecHits;
598     iEvent.getByLabel(EcalRecHitProducers[i], h_Ecal_RecHits);
599     EcalRecHitCollection Ecal_RecHits = *h_Ecal_RecHits.product();
600    
601     for(unsigned int eh=0;eh<Ecal_RecHits.size();eh++){
602 querten 1.4 MyCaloHit* temp_EcalHit = new MyCaloHit;
603     temp_EcalHit->E = Ecal_RecHits[eh].energy();
604     temp_EcalHit->t = Ecal_RecHits[eh].time();
605     temp_EcalHit->DetId = (Ecal_RecHits[eh].detid()).rawId();
606 querten 1.1
607     MyEvent->MyEcalCaloHitCollection.push_back(temp_EcalHit);
608     }
609     }
610    
611    
612 querten 1.2 // ### HCALRecHits ###
613 querten 1.1 for(unsigned int i=0;i<HcalHBHERecHitProducers.size();i++){
614     edm::Handle<HBHERecHitCollection > h_HcalHBHE_RecHits;
615     iEvent.getByLabel(HcalHBHERecHitProducers[i], h_HcalHBHE_RecHits);
616     HBHERecHitCollection HcalHBHE_RecHits = *h_HcalHBHE_RecHits.product();
617    
618     for(unsigned int hh=0;hh<HcalHBHE_RecHits.size();hh++){
619 querten 1.4 MyCaloHit* temp_HcalHBHEHit = new MyCaloHit;
620     temp_HcalHBHEHit->E = HcalHBHE_RecHits[hh].energy();
621     temp_HcalHBHEHit->t = HcalHBHE_RecHits[hh].time();
622     temp_HcalHBHEHit->DetId = (HcalHBHE_RecHits[hh].detid()).rawId();
623 querten 1.1
624     MyEvent->MyHcalCaloHitCollection.push_back(temp_HcalHBHEHit);
625     }
626     }
627    
628    
629     for(unsigned int i=0;i<HcalHORecHitProducers.size();i++){
630     edm::Handle<HORecHitCollection > h_HcalHO_RecHits;
631     iEvent.getByLabel(HcalHORecHitProducers[i], h_HcalHO_RecHits);
632     HORecHitCollection HcalHO_RecHits = *h_HcalHO_RecHits.product();
633    
634     for(unsigned int hh=0;hh<HcalHO_RecHits.size();hh++){
635 querten 1.4 MyCaloHit* temp_HcalHOHit = new MyCaloHit;
636     temp_HcalHOHit->E = HcalHO_RecHits[hh].energy();
637     temp_HcalHOHit->t = HcalHO_RecHits[hh].time();
638     temp_HcalHOHit->DetId = (HcalHO_RecHits[hh].detid()).rawId();
639 querten 1.1
640     MyEvent->MyHcalCaloHitCollection.push_back(temp_HcalHOHit);
641     }
642     }
643    
644    
645     for(unsigned int i=0;i<HcalHFRecHitProducers.size();i++){
646     edm::Handle<HFRecHitCollection > h_HcalHF_RecHits;
647     iEvent.getByLabel(HcalHFRecHitProducers[i], h_HcalHF_RecHits);
648     HFRecHitCollection HcalHF_RecHits = *h_HcalHF_RecHits.product();
649    
650     for(unsigned int hh=0;hh<HcalHF_RecHits.size();hh++){
651 querten 1.4 MyCaloHit* temp_HcalHFHit = new MyCaloHit;
652     temp_HcalHFHit->E = HcalHF_RecHits[hh].energy();
653     temp_HcalHFHit->t = HcalHF_RecHits[hh].time();
654     temp_HcalHFHit->DetId = (HcalHF_RecHits[hh].detid()).rawId();
655 querten 1.1
656     MyEvent->MyHcalCaloHitCollection.push_back(temp_HcalHFHit);
657     }
658     }
659    
660 querten 1.2 // ### Muon Segments ###
661     for(unsigned int i=0;i<CSCSegmentProducers.size();i++){
662     edm::Handle<CSCSegmentCollection > h_CSC_Segments;
663     iEvent.getByLabel(CSCSegmentProducers[i], h_CSC_Segments);
664     CSCSegmentCollection CSC_Segments = *h_CSC_Segments.product();
665    
666     for(unsigned int s=0;s<CSC_Segments.size();s++){
667     DetId theDetUnitId = CSC_Segments[s].geographicalId();
668     const GeomDet* theDet = CscGeom->idToDet(theDetUnitId);
669    
670 querten 1.4 MyMuonSegment* temp_CscSeg = new MyMuonSegment;
671     temp_CscSeg->DetId = theDetUnitId.rawId();
672     temp_CscSeg->PosX = theDet->surface().toGlobal(CSC_Segments[s].localPosition()).x();
673     temp_CscSeg->PosY = theDet->surface().toGlobal(CSC_Segments[s].localPosition()).y();
674     temp_CscSeg->PosZ = theDet->surface().toGlobal(CSC_Segments[s].localPosition()).z();
675     temp_CscSeg->DirX = theDet->surface().toGlobal(CSC_Segments[s].localDirection()).x();
676     temp_CscSeg->DirY = theDet->surface().toGlobal(CSC_Segments[s].localDirection()).y();
677     temp_CscSeg->DirZ = theDet->surface().toGlobal(CSC_Segments[s].localDirection()).z();
678 querten 1.2
679     MyEvent->MyMuonSegmentCollection.push_back(temp_CscSeg);
680     }
681     }
682    
683     for(unsigned int i=0;i<DTSegmentProducers.size();i++){
684     edm::Handle<DTRecSegment4DCollection > h_DT_Segments;
685     iEvent.getByLabel(DTSegmentProducers[i], h_DT_Segments);
686     DTRecSegment4DCollection DT_Segments = *h_DT_Segments.product();
687    
688     for(unsigned int s=0;s<DT_Segments.size();s++){
689     DetId theDetUnitId = DT_Segments[s].geographicalId();
690     const GeomDet* theDet = DtGeom->idToDet(theDetUnitId);
691    
692 querten 1.4 MyMuonSegment* temp_DtSeg = new MyMuonSegment;
693     temp_DtSeg->DetId = theDetUnitId.rawId();
694     temp_DtSeg->PosX = theDet->surface().toGlobal(DT_Segments[s].localPosition()).x();
695     temp_DtSeg->PosY = theDet->surface().toGlobal(DT_Segments[s].localPosition()).y();
696     temp_DtSeg->PosZ = theDet->surface().toGlobal(DT_Segments[s].localPosition()).z();
697     temp_DtSeg->DirX = theDet->surface().toGlobal(DT_Segments[s].localDirection()).x();
698     temp_DtSeg->DirY = theDet->surface().toGlobal(DT_Segments[s].localDirection()).y();
699     temp_DtSeg->DirZ = theDet->surface().toGlobal(DT_Segments[s].localDirection()).z();
700 querten 1.2
701     MyEvent->MyMuonSegmentCollection.push_back(temp_DtSeg);
702     }
703     }
704    
705     // ### Muon Hits ###
706     for(unsigned int i=0;i<RPCHitsProducers.size();i++){
707     edm::Handle<RPCRecHitCollection > h_RPC_Hits;
708     iEvent.getByLabel(RPCHitsProducers[i], h_RPC_Hits);
709     RPCRecHitCollection RPC_Hits = *h_RPC_Hits.product();
710    
711     for(unsigned int h=0;h<RPC_Hits.size();h++){
712     DetId theDetUnitId = RPC_Hits[h].geographicalId();
713     const GeomDet* theDet = RpcGeom->idToDet(theDetUnitId);
714    
715 querten 1.4 MyMuonHit* temp_RpcHit = new MyMuonHit;
716     temp_RpcHit->DetId = theDetUnitId.rawId();
717     temp_RpcHit->x = theDet->surface().toGlobal(RPC_Hits[h].localPosition()).x();
718     temp_RpcHit->y = theDet->surface().toGlobal(RPC_Hits[h].localPosition()).y();
719     temp_RpcHit->z = theDet->surface().toGlobal(RPC_Hits[h].localPosition()).z();
720 querten 1.2
721     MyEvent->MyMuonHitCollection.push_back(temp_RpcHit);
722     }
723     }
724    
725 querten 1.1 MyEvents->Events.push_back(MyEvent);
726     }
727    
728    
729     //define this as a plug-in
730     DEFINE_FWK_MODULE(Frog_Analyzer);
731    
732