ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/FastOpenGlDisplayer/src/Frog_Analyzer.cc
Revision: 1.5
Committed: Tue Apr 29 16:42:16 2008 UTC (17 years ago) by roberfro
Content type: text/plain
Branch: MAIN
Changes since 1.4: +21 -11 lines
Log Message:
Use of TrajectoryMeasurement instead of rechits

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 roberfro 1.5 // $Id: Frog_Analyzer.cc,v 1.4 2008/04/27 07:27:30 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.1
84 querten 1.2 #include "DataFormats/DTRecHit/interface/DTRecSegment4DCollection.h"
85     #include "DataFormats/CSCRecHit/interface/CSCSegmentCollection.h"
86    
87     #include "DataFormats/DTRecHit/interface/DTRecHitCollection.h"
88     #include "DataFormats/CSCRecHit/interface/CSCRecHit2DCollection.h"
89     #include "DataFormats/RPCRecHit/interface/RPCRecHitCollection.h"
90    
91 querten 1.1 #include "DataFormats/Math/interface/Point3D.h"
92     #include "DataFormats/GeometrySurface/interface/TrapezoidalPlaneBounds.h"
93     #include "DataFormats/GeometrySurface/interface/RectangularPlaneBounds.h"
94    
95     #include "Visualisation/Frog/interface/SimEvent.h"
96     #include "Visualisation/Frog/interface/Geometry.h"
97    
98     using namespace edm;
99     using namespace std;
100    
101     //
102     // class decleration
103     //
104    
105     class Frog_Analyzer : public edm::EDAnalyzer {
106     public:
107     explicit Frog_Analyzer(const edm::ParameterSet&);
108     ~Frog_Analyzer();
109    
110    
111     private:
112     virtual void beginJob(const edm::EventSetup& iSetup);
113     virtual void analyze (const edm::Event&, const edm::EventSetup&);
114     virtual void endJob ();
115    
116    
117     MySimEvents* MyEvents;
118     Geometry* Geom_Tracker;
119     Geometry* Geom_ECAL;
120     Geometry* Geom_HCAL;
121     Geometry* Geom_Muon;
122    
123    
124     std::string OutputFile;
125    
126     std::vector<InputTag> SimTrackProducers;
127     std::vector<InputTag> SimVertexProducers;
128     std::vector<InputTag> SimHitProducers;
129    
130 roberfro 1.5 std::vector<InputTag> TrajectoryProducers;
131 querten 1.1 std::vector<InputTag> EcalRecHitProducers;
132     std::vector<InputTag> HcalHBHERecHitProducers;
133     std::vector<InputTag> HcalHORecHitProducers;
134     std::vector<InputTag> HcalHFRecHitProducers;
135 querten 1.2 std::vector<InputTag> DTSegmentProducers;
136     std::vector<InputTag> CSCSegmentProducers;
137    
138     std::vector<InputTag> RPCHitsProducers;
139 querten 1.1
140    
141    
142    
143     // ----------member data ---------------------------
144     };
145    
146     //
147     // constructors and destructor
148     //
149     Frog_Analyzer::Frog_Analyzer(const edm::ParameterSet& iConfig)
150     {
151     OutputFile = iConfig.getParameter<std::string >("OutputFile");
152    
153     SimTrackProducers = iConfig.getParameter<std::vector<InputTag> >("SimTrackProducers");
154     SimVertexProducers = iConfig.getParameter<std::vector<InputTag> >("SimVertexProducers");
155     SimHitProducers = iConfig.getParameter<std::vector<InputTag> >("SimHitProducers");
156    
157 roberfro 1.5 TrajectoryProducers = iConfig.getParameter<std::vector<InputTag> >("TrajectoryProducers");
158 querten 1.1 EcalRecHitProducers = iConfig.getParameter<std::vector<InputTag> >("EcalRecHitProducers");
159     HcalHBHERecHitProducers = iConfig.getParameter<std::vector<InputTag> >("HcalHBHERecHitProducers");
160     HcalHORecHitProducers = iConfig.getParameter<std::vector<InputTag> >("HcalHORecHitProducers");
161     HcalHFRecHitProducers = iConfig.getParameter<std::vector<InputTag> >("HcalHFRecHitProducers");
162    
163 querten 1.2 DTSegmentProducers = iConfig.getParameter<std::vector<InputTag> >("DTSegmentProducers");
164     CSCSegmentProducers = iConfig.getParameter<std::vector<InputTag> >("CSCSegmentProducers");
165    
166     RPCHitsProducers = iConfig.getParameter<std::vector<InputTag> >("RPCHitsProducers");
167    
168 querten 1.1 }
169    
170    
171     Frog_Analyzer::~Frog_Analyzer()
172     {
173     }
174    
175     // ------------ method called once each job just before starting event loop ------------
176     void
177     Frog_Analyzer::beginJob(const edm::EventSetup& iSetup)
178     {
179     DetId Detid;
180     int SubDet;
181    
182     MyEvents = new MySimEvents();
183     Geom_Tracker = new Geometry();
184     Geom_ECAL = new Geometry();
185     Geom_HCAL = new Geometry();
186     Geom_Muon = new Geometry();
187    
188    
189     // ### TRACKER GEOMETRY ###
190    
191     edm::ESHandle<TrackerGeometry> tkGeom;
192     iSetup.get<TrackerDigiGeometryRecord>().get( tkGeom );
193     vector<GeomDet*> TrackerDets = tkGeom->dets();
194    
195     for(unsigned int i=0;i<TrackerDets.size();i++){
196     Detid = TrackerDets[i]->geographicalId();
197     // SubDet = Detid.subdetId();
198    
199     GeomDet* DetUnit = TrackerDets[i];
200     if(!DetUnit)continue;
201     const BoundPlane plane = DetUnit->surface();
202     const TrapezoidalPlaneBounds* trapezoidalBounds( dynamic_cast<const TrapezoidalPlaneBounds*>(&(plane.bounds())));
203     const RectangularPlaneBounds* rectangularBounds( dynamic_cast<const RectangularPlaneBounds*>(&(plane.bounds())));
204    
205     float width = 0;
206     float length = 0;
207     float thickness = 0;
208     float TrapezoidalParam = 0;
209    
210     if(trapezoidalBounds)
211     {
212     std::vector<float> const & parameters = (*trapezoidalBounds).parameters();
213     width = parameters[0]*2;
214     length = parameters[3]*2;
215     thickness = (*trapezoidalBounds).thickness();
216     TrapezoidalParam = parameters[1]/parameters[0];
217     }else if(rectangularBounds){
218     width = DetUnit->surface().bounds().width();
219     length = DetUnit->surface().bounds().length();
220     thickness = DetUnit->surface().bounds().thickness();
221     TrapezoidalParam = 1;
222     }
223    
224     Surface::GlobalPoint WidthVector = plane.toGlobal( LocalPoint(width/2, 0, 0) );
225     Surface::GlobalPoint LengthVector = plane.toGlobal( LocalPoint(0, length/2, 0) );
226     Surface::GlobalPoint ThickVector = plane.toGlobal( LocalPoint(0, 0, thickness/2) );
227    
228     GlobalVector Pos = GlobalVector(DetUnit->position().basicVector());
229    
230     Geom_Tracker->Add_TrackerDet(Detid.rawId(), TrapezoidalParam,
231     Pos.x(), Pos.y(), Pos.z(),
232     WidthVector.x() -Pos.x(), WidthVector.y() -Pos.y(), WidthVector.z() -Pos.z(),
233     LengthVector.x()-Pos.x(), LengthVector.y()-Pos.y(), LengthVector.z()-Pos.z(),
234     ThickVector.x() -Pos.x(), ThickVector.y() -Pos.y(), ThickVector.z() -Pos.z());
235    
236    
237     }
238    
239    
240     // ### CALO GEOMETRY ###
241    
242     edm::ESHandle<CaloGeometry> CaloGeom;
243     iSetup.get<IdealGeometryRecord>().get( CaloGeom );
244     const vector<DetId> CaloDets = CaloGeom->getValidDetIds();
245    
246     for(unsigned int i=0;i<CaloDets.size();i++)
247     {
248     Detid = CaloDets[i];
249     SubDet = Detid.subdetId();
250    
251     Geometry* Geom_temp = NULL;
252     if(Detid.det()==DetId::Ecal){
253     Geom_temp = Geom_ECAL;
254    
255     if(SubDet<1 || SubDet>3){
256     printf("Don't Save EcalTriggerTower or EcalLaserPnDiode\n");
257     continue;
258     }
259     }else if(Detid.det()==DetId::Hcal){
260     Geom_temp = Geom_HCAL;
261    
262     if(SubDet<1 || SubDet>4){
263     printf("Don't Save HcalEmpty, HcalTriggerTower, HcalOther\n");
264     continue;
265     }
266     }else{
267     continue;
268     }
269    
270     const CaloCellGeometry* CellGeom = CaloGeom->getGeometry(Detid);
271     GlobalPoint CellPos = CellGeom->getPosition();
272     const CaloCellGeometry::CornersVec CellCorners = CellGeom->getCorners();
273    
274 querten 1.4 float cX = (CellCorners[0].x() + CellCorners[2].x())/2;
275     float cY = (CellCorners[0].y() + CellCorners[2].y())/2;
276     float cZ = (CellCorners[0].z() + CellCorners[2].z())/2;
277    
278     float wX = (CellCorners[1].x() - CellCorners[0].x())/2;
279     float wY = (CellCorners[1].y() - CellCorners[0].y())/2;
280     float wZ = (CellCorners[1].z() - CellCorners[0].z())/2;
281    
282     float hX = (CellCorners[3].x() - CellCorners[0].x())/2;
283     float hY = (CellCorners[3].y() - CellCorners[0].y())/2;
284     float hZ = (CellCorners[3].z() - CellCorners[0].z())/2;
285    
286     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) );
287     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) );
288    
289 querten 1.1 Geom_temp->Add_CaloDet(Detid.rawId(),
290 querten 1.4 cX , cY , cZ ,
291     wX , wY , wZ ,
292     hX , hY , hZ , F );
293 querten 1.1 }
294    
295     // ### MUON GEOMETRY ###
296    
297     edm::ESHandle<DTGeometry> DtGeom;
298     iSetup.get<MuonGeometryRecord>().get( DtGeom );
299     const vector<GeomDet*> DtDets = DtGeom->dets();
300    
301     edm::ESHandle<CSCGeometry> CscGeom;
302     iSetup.get<MuonGeometryRecord>().get( CscGeom );
303     const vector<GeomDet*> CscDets = CscGeom->dets();
304    
305     edm::ESHandle<RPCGeometry> RpcGeom;
306     iSetup.get<MuonGeometryRecord>().get( RpcGeom );
307     const vector<GeomDet*> RpcDets = RpcGeom->dets();
308    
309     vector<GeomDet*> MuonDets;
310     for(unsigned int i=0;i<DtDets.size() ;i++){MuonDets.push_back(DtDets [i]);}
311     for(unsigned int i=0;i<CscDets.size();i++){MuonDets.push_back(CscDets[i]);}
312     for(unsigned int i=0;i<RpcDets.size();i++){MuonDets.push_back(RpcDets[i]);}
313    
314    
315     for(unsigned int i=0;i<MuonDets.size();i++)
316     {
317     Detid = DetId(MuonDets[i]->geographicalId());
318     SubDet = Detid.subdetId();
319    
320     GeomDet* DetUnit = MuonDets[i];
321     if(!DetUnit)continue;
322     const BoundPlane plane = DetUnit->surface();
323     const TrapezoidalPlaneBounds* trapezoidalBounds( dynamic_cast<const TrapezoidalPlaneBounds*>(&(plane.bounds())));
324     const RectangularPlaneBounds* rectangularBounds( dynamic_cast<const RectangularPlaneBounds*>(&(plane.bounds())));
325    
326     float width = 0;
327     float length = 0;
328     float thickness = 0;
329     float TrapezoidalParam = 0;
330    
331     if(trapezoidalBounds)
332     {
333     std::vector<float> const & parameters = (*trapezoidalBounds).parameters();
334     width = parameters[0]*2;
335     length = parameters[3]*2;
336     thickness = (*trapezoidalBounds).thickness();
337     TrapezoidalParam = parameters[1]/parameters[0];
338     }else if(rectangularBounds){
339     width = DetUnit->surface().bounds().width();
340     length = DetUnit->surface().bounds().length();
341     thickness = DetUnit->surface().bounds().thickness();
342     TrapezoidalParam = 1;
343     }
344    
345     Surface::GlobalPoint WidthVector = plane.toGlobal( LocalPoint(width/2, 0, 0) );
346     Surface::GlobalPoint LengthVector = plane.toGlobal( LocalPoint(0, length/2, 0) );
347     Surface::GlobalPoint ThickVector = plane.toGlobal( LocalPoint(0, 0, thickness/2) );
348    
349     GlobalVector Pos = GlobalVector(DetUnit->position().basicVector());
350    
351     Geom_Muon->Add_TrackerDet(Detid.rawId(), TrapezoidalParam,
352     Pos.x(), Pos.y(), Pos.z(),
353     WidthVector.x() -Pos.x(), WidthVector.y() -Pos.y(), WidthVector.z() -Pos.z(),
354     LengthVector.x()-Pos.x(), LengthVector.y()-Pos.y(), LengthVector.z()-Pos.z(),
355     ThickVector.x() -Pos.x(), ThickVector.y() -Pos.y(), ThickVector.z() -Pos.z());
356     }
357    
358    
359     // ### Save .geom ###
360    
361     Geom_Tracker->Save("Tracker.geom");
362     Geom_ECAL ->Save("Ecal.geom");
363     Geom_HCAL ->Save("Hcal.geom");
364     Geom_Muon ->Save("Muon.geom");
365    
366     }
367    
368     // ------------ method called once each job just after ending the event loop ------------
369     void
370     Frog_Analyzer::endJob() {
371    
372     MyEvents->Save((char*) OutputFile.c_str());
373     // MyEvents->Load((char*) OutputFile.c_str());
374     }
375    
376    
377    
378     //
379     // member functions
380     //
381    
382     // ------------ method called to for each event ------------
383     void
384     Frog_Analyzer::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup)
385     {
386     MySimEvent* MyEvent = new MySimEvent;
387    
388 querten 1.2 // access the tracker
389     edm::ESHandle<TrackerGeometry> tkGeom;
390     iSetup.get<TrackerDigiGeometryRecord>().get( tkGeom );
391 querten 1.1
392 querten 1.2 edm::ESHandle<DTGeometry> DtGeom;
393     iSetup.get<MuonGeometryRecord>().get( DtGeom );
394    
395     edm::ESHandle<CSCGeometry> CscGeom;
396     iSetup.get<MuonGeometryRecord>().get( CscGeom );
397 querten 1.1
398 querten 1.2 edm::ESHandle<RPCGeometry> RpcGeom;
399     iSetup.get<MuonGeometryRecord>().get( RpcGeom );
400 querten 1.1
401    
402     // Save Data
403 querten 1.2 // ### SIMTRACK ###
404 querten 1.1 for(unsigned int i=0;i<SimTrackProducers.size();i++){
405     edm::Handle<std::vector< SimTrack > > h_SimTracks;
406     iEvent.getByLabel(SimTrackProducers[i], h_SimTracks);
407     std::vector< SimTrack > SimTrackColl = *h_SimTracks.product();
408    
409     for ( unsigned int a = 0; a < SimTrackColl.size(); ++a ) {
410 querten 1.4 MySimTrack* MysimTrack = new MySimTrack;
411 querten 1.1 SimTrack simTrack =SimTrackColl[a];
412    
413 querten 1.4 MysimTrack->track_id =simTrack.trackId();
414     MysimTrack->Type =simTrack.type();
415     MysimTrack->parent_vertex=simTrack.vertIndex();
416     MysimTrack->Px =simTrack.momentum().x();
417     MysimTrack->Py =simTrack.momentum().y();
418     MysimTrack->Pz =simTrack.momentum().z();
419     MysimTrack->E =simTrack.momentum().e();
420     MysimTrack->charge =simTrack.charge();
421 querten 1.1
422     MyEvent->MySimTrackCollection.push_back(MysimTrack);
423     }
424     }
425    
426 querten 1.2 // ### SIMVERTEX ###
427 querten 1.1 for(unsigned int i=0;i<SimVertexProducers.size();i++){
428     edm::Handle<std::vector< SimVertex > > h_Vertex;
429     iEvent.getByLabel(SimVertexProducers[i], h_Vertex);
430     std::vector< SimVertex > VertexColl = *h_Vertex.product();
431    
432     for (unsigned int b = 0; b < VertexColl.size(); ++b ) {
433 querten 1.4 MySimVertex* MyVertex = new MySimVertex;
434 querten 1.1 SimVertex Vertex = VertexColl[b];
435    
436 querten 1.4 MyVertex->parentTrack_id=Vertex.parentIndex ();
437     MyVertex->x =Vertex.position().x();
438     MyVertex->y =Vertex.position().y();
439     MyVertex->z =Vertex.position().z();
440 querten 1.1
441     MyEvent->MySimVertexCollection.push_back(MyVertex);
442     }
443     }
444    
445 querten 1.2 // ### SIMHIT ###
446 querten 1.1 for(unsigned int i=0;i<SimHitProducers.size();i++){
447     edm::Handle<std::vector< PSimHit > > h_Hits;
448     iEvent.getByLabel(SimHitProducers[i], h_Hits);
449     std::vector< PSimHit > Hits = *h_Hits.product();
450    
451     for(unsigned int h=0; h<Hits.size(); h++)
452     {
453 querten 1.2 DetId theDetUnitId(Hits[h].detUnitId());
454     const GeomDet* theDet = tkGeom->idToDet(theDetUnitId);
455 querten 1.1
456 querten 1.4 MyPSimHit* Hit = new MyPSimHit;
457     Hit->x = theDet->surface().toGlobal(Hits[h].localPosition()).x();
458     Hit->y = theDet->surface().toGlobal(Hits[h].localPosition()).y();
459     Hit->z = theDet->surface().toGlobal(Hits[h].localPosition()).z();
460     Hit->ProcessType = Hits[h].processType();
461     Hit->dEdX = Hits[h].energyLoss();
462 querten 1.1
463     MyEvent->MyPSimHitCollection.push_back(Hit);
464     }
465     }
466    
467    
468 querten 1.2 // ### TRACKS ###
469 roberfro 1.5 for(unsigned int i=0;i<TrajectoryProducers.size();i++){
470     edm::Handle<std::vector< Trajectory > > h_Trajs;
471     iEvent.getByLabel(TrajectoryProducers[i], h_Trajs);
472     std::vector< Trajectory > TrajColl = *h_Trajs.product();
473 querten 1.1
474 roberfro 1.5 for ( unsigned int t = 0; t < TrajColl.size(); ++t ) {
475 querten 1.4 MyRecoTrack* MyrecoTrack = new MyRecoTrack;
476 roberfro 1.5 Trajectory recoTraj =TrajColl[t];
477 querten 1.1
478 roberfro 1.5 std::vector<TrajectoryMeasurement> measurements = recoTraj.measurements();
479    
480     for(unsigned int h=0;h<measurements.size();h++){
481    
482     TrajectoryMeasurement::ConstRecHitPointer h_it = measurements[h].recHit();
483 querten 1.1 if(!h_it->isValid() )continue;
484     DetId detId = h_it->geographicalId();
485 querten 1.2 const GeomDet * theDet = tkGeom->idToDet(detId);
486 querten 1.1 LocalPoint localPos = h_it->localPosition();
487 roberfro 1.5 GlobalPoint globalPos = measurements[h].updatedState().globalPosition();
488 querten 1.1
489 querten 1.4 MyRecoHit* hit = new MyRecoHit;
490 roberfro 1.5 /*
491 querten 1.4 hit->x = theDet->surface().toGlobal(localPos).x();
492     hit->y = theDet->surface().toGlobal(localPos).y();
493     hit->z = theDet->surface().toGlobal(localPos).z();
494 roberfro 1.5 */
495     hit->x = globalPos.x();
496     hit->y = globalPos.y();
497     hit->z = globalPos.z();
498 querten 1.4 hit->DetId = detId.rawId();
499     hit->Charge = -1;
500     MyrecoTrack->Hits.push_back(hit);
501 querten 1.1 }
502     MyEvent->MyRecoTrackCollection.push_back(MyrecoTrack);
503     }
504     }
505    
506    
507 querten 1.2 // ### ECALRecHits ###
508 querten 1.1 for(unsigned int i=0;i<EcalRecHitProducers.size();i++){
509     edm::Handle<EcalRecHitCollection > h_Ecal_RecHits;
510     iEvent.getByLabel(EcalRecHitProducers[i], h_Ecal_RecHits);
511     EcalRecHitCollection Ecal_RecHits = *h_Ecal_RecHits.product();
512    
513     for(unsigned int eh=0;eh<Ecal_RecHits.size();eh++){
514 querten 1.4 MyCaloHit* temp_EcalHit = new MyCaloHit;
515     temp_EcalHit->E = Ecal_RecHits[eh].energy();
516     temp_EcalHit->t = Ecal_RecHits[eh].time();
517     temp_EcalHit->DetId = (Ecal_RecHits[eh].detid()).rawId();
518 querten 1.1
519     MyEvent->MyEcalCaloHitCollection.push_back(temp_EcalHit);
520     }
521     }
522    
523    
524 querten 1.2 // ### HCALRecHits ###
525 querten 1.1 for(unsigned int i=0;i<HcalHBHERecHitProducers.size();i++){
526     edm::Handle<HBHERecHitCollection > h_HcalHBHE_RecHits;
527     iEvent.getByLabel(HcalHBHERecHitProducers[i], h_HcalHBHE_RecHits);
528     HBHERecHitCollection HcalHBHE_RecHits = *h_HcalHBHE_RecHits.product();
529    
530     for(unsigned int hh=0;hh<HcalHBHE_RecHits.size();hh++){
531 querten 1.4 MyCaloHit* temp_HcalHBHEHit = new MyCaloHit;
532     temp_HcalHBHEHit->E = HcalHBHE_RecHits[hh].energy();
533     temp_HcalHBHEHit->t = HcalHBHE_RecHits[hh].time();
534     temp_HcalHBHEHit->DetId = (HcalHBHE_RecHits[hh].detid()).rawId();
535 querten 1.1
536     MyEvent->MyHcalCaloHitCollection.push_back(temp_HcalHBHEHit);
537     }
538     }
539    
540    
541     for(unsigned int i=0;i<HcalHORecHitProducers.size();i++){
542     edm::Handle<HORecHitCollection > h_HcalHO_RecHits;
543     iEvent.getByLabel(HcalHORecHitProducers[i], h_HcalHO_RecHits);
544     HORecHitCollection HcalHO_RecHits = *h_HcalHO_RecHits.product();
545    
546     for(unsigned int hh=0;hh<HcalHO_RecHits.size();hh++){
547 querten 1.4 MyCaloHit* temp_HcalHOHit = new MyCaloHit;
548     temp_HcalHOHit->E = HcalHO_RecHits[hh].energy();
549     temp_HcalHOHit->t = HcalHO_RecHits[hh].time();
550     temp_HcalHOHit->DetId = (HcalHO_RecHits[hh].detid()).rawId();
551 querten 1.1
552     MyEvent->MyHcalCaloHitCollection.push_back(temp_HcalHOHit);
553     }
554     }
555    
556    
557     for(unsigned int i=0;i<HcalHFRecHitProducers.size();i++){
558     edm::Handle<HFRecHitCollection > h_HcalHF_RecHits;
559     iEvent.getByLabel(HcalHFRecHitProducers[i], h_HcalHF_RecHits);
560     HFRecHitCollection HcalHF_RecHits = *h_HcalHF_RecHits.product();
561    
562     for(unsigned int hh=0;hh<HcalHF_RecHits.size();hh++){
563 querten 1.4 MyCaloHit* temp_HcalHFHit = new MyCaloHit;
564     temp_HcalHFHit->E = HcalHF_RecHits[hh].energy();
565     temp_HcalHFHit->t = HcalHF_RecHits[hh].time();
566     temp_HcalHFHit->DetId = (HcalHF_RecHits[hh].detid()).rawId();
567 querten 1.1
568     MyEvent->MyHcalCaloHitCollection.push_back(temp_HcalHFHit);
569     }
570     }
571    
572 querten 1.2 // ### Muon Segments ###
573     for(unsigned int i=0;i<CSCSegmentProducers.size();i++){
574     edm::Handle<CSCSegmentCollection > h_CSC_Segments;
575     iEvent.getByLabel(CSCSegmentProducers[i], h_CSC_Segments);
576     CSCSegmentCollection CSC_Segments = *h_CSC_Segments.product();
577    
578     for(unsigned int s=0;s<CSC_Segments.size();s++){
579     DetId theDetUnitId = CSC_Segments[s].geographicalId();
580     const GeomDet* theDet = CscGeom->idToDet(theDetUnitId);
581    
582 querten 1.4 MyMuonSegment* temp_CscSeg = new MyMuonSegment;
583     temp_CscSeg->DetId = theDetUnitId.rawId();
584     temp_CscSeg->PosX = theDet->surface().toGlobal(CSC_Segments[s].localPosition()).x();
585     temp_CscSeg->PosY = theDet->surface().toGlobal(CSC_Segments[s].localPosition()).y();
586     temp_CscSeg->PosZ = theDet->surface().toGlobal(CSC_Segments[s].localPosition()).z();
587     temp_CscSeg->DirX = theDet->surface().toGlobal(CSC_Segments[s].localDirection()).x();
588     temp_CscSeg->DirY = theDet->surface().toGlobal(CSC_Segments[s].localDirection()).y();
589     temp_CscSeg->DirZ = theDet->surface().toGlobal(CSC_Segments[s].localDirection()).z();
590 querten 1.2
591     MyEvent->MyMuonSegmentCollection.push_back(temp_CscSeg);
592     }
593     }
594    
595     for(unsigned int i=0;i<DTSegmentProducers.size();i++){
596     edm::Handle<DTRecSegment4DCollection > h_DT_Segments;
597     iEvent.getByLabel(DTSegmentProducers[i], h_DT_Segments);
598     DTRecSegment4DCollection DT_Segments = *h_DT_Segments.product();
599    
600     for(unsigned int s=0;s<DT_Segments.size();s++){
601     DetId theDetUnitId = DT_Segments[s].geographicalId();
602     const GeomDet* theDet = DtGeom->idToDet(theDetUnitId);
603    
604 querten 1.4 MyMuonSegment* temp_DtSeg = new MyMuonSegment;
605     temp_DtSeg->DetId = theDetUnitId.rawId();
606     temp_DtSeg->PosX = theDet->surface().toGlobal(DT_Segments[s].localPosition()).x();
607     temp_DtSeg->PosY = theDet->surface().toGlobal(DT_Segments[s].localPosition()).y();
608     temp_DtSeg->PosZ = theDet->surface().toGlobal(DT_Segments[s].localPosition()).z();
609     temp_DtSeg->DirX = theDet->surface().toGlobal(DT_Segments[s].localDirection()).x();
610     temp_DtSeg->DirY = theDet->surface().toGlobal(DT_Segments[s].localDirection()).y();
611     temp_DtSeg->DirZ = theDet->surface().toGlobal(DT_Segments[s].localDirection()).z();
612 querten 1.2
613     MyEvent->MyMuonSegmentCollection.push_back(temp_DtSeg);
614     }
615     }
616    
617     // ### Muon Hits ###
618     for(unsigned int i=0;i<RPCHitsProducers.size();i++){
619     edm::Handle<RPCRecHitCollection > h_RPC_Hits;
620     iEvent.getByLabel(RPCHitsProducers[i], h_RPC_Hits);
621     RPCRecHitCollection RPC_Hits = *h_RPC_Hits.product();
622    
623     for(unsigned int h=0;h<RPC_Hits.size();h++){
624     DetId theDetUnitId = RPC_Hits[h].geographicalId();
625     const GeomDet* theDet = RpcGeom->idToDet(theDetUnitId);
626    
627 querten 1.4 MyMuonHit* temp_RpcHit = new MyMuonHit;
628     temp_RpcHit->DetId = theDetUnitId.rawId();
629     temp_RpcHit->x = theDet->surface().toGlobal(RPC_Hits[h].localPosition()).x();
630     temp_RpcHit->y = theDet->surface().toGlobal(RPC_Hits[h].localPosition()).y();
631     temp_RpcHit->z = theDet->surface().toGlobal(RPC_Hits[h].localPosition()).z();
632 querten 1.2
633     MyEvent->MyMuonHitCollection.push_back(temp_RpcHit);
634     }
635     }
636    
637 querten 1.1 MyEvents->Events.push_back(MyEvent);
638     }
639    
640    
641     //define this as a plug-in
642     DEFINE_FWK_MODULE(Frog_Analyzer);
643    
644