ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/FastOpenGlDisplayer/src/Frog_Analyzer.cc
Revision: 1.6
Committed: Thu May 1 06:52:20 2008 UTC (17 years ago) by querten
Content type: text/plain
Branch: MAIN
CVS Tags: Version_0_26, Version_0_25
Changes since 1.5: +19 -11 lines
Log Message:
Add Some Track Info

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