ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/FastOpenGlDisplayer/src/Frog_Analyzer.cc
Revision: 1.4
Committed: Sun Apr 27 07:27:30 2008 UTC (17 years ago) by querten
Content type: text/plain
Branch: MAIN
CVS Tags: Version_0_24
Changes since 1.3: +84 -75 lines
Log Message:
Improvements in Loading and File Size

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