CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_5/src/CalibTracker/SiStripLorentzAngle/plugins/SiStripLAProfileBooker.cc

Go to the documentation of this file.
00001 /*  VI Janurary 2012 
00002  * This file need to be migrated to the new interface of matched hit as the mono/stero are not there anymore
00003  * what is returned are hits w/o localpoistion, just cluster and id
00004  */
00005 #include <string>
00006 #include <iostream>
00007 #include <fstream>
00008 #include <functional>
00009 
00010 
00011 #include "CalibTracker/SiStripLorentzAngle/interface/SiStripLAProfileBooker.h"
00012 #include "CalibTracker/Records/interface/SiStripDetCablingRcd.h"
00013 #include "CalibFormats/SiStripObjects/interface/SiStripDetCabling.h"
00014 
00015 #include "FWCore/Framework/interface/ESHandle.h"
00016 #include "FWCore/Framework/interface/EventSetup.h"
00017 #include "FWCore/Framework/interface/Event.h"
00018 #include "FWCore/ServiceRegistry/interface/Service.h"
00019 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00020 
00021 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
00022 #include "Geometry/CommonDetUnit/interface/GeomDetUnit.h"
00023 #include "Geometry/TrackerGeometryBuilder/interface/GluedGeomDet.h"
00024 #include "Geometry/CommonTopologies/interface/StripTopology.h"
00025 
00026 #include "DataFormats/Common/interface/Handle.h"
00027 #include "DataFormats/TrackReco/interface/Track.h"
00028 #include "DataFormats/SiStripDetId/interface/StripSubdetector.h"
00029 #include "DataFormats/TrackerCommon/interface/TrackerTopology.h"
00030 #include "Geometry/Records/interface/IdealGeometryRecord.h"
00031 #include "DataFormats/TrackerRecHit2D/interface/SiStripMatchedRecHit2D.h"
00032 #include "DataFormats/TrackerRecHit2D/interface/ProjectedSiStripRecHit2D.h"
00033 #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
00034 #include "DataFormats/GeometryVector/interface/GlobalVector.h"
00035 
00036 #include "DQM/SiStripCommon/interface/ExtractTObject.h"
00037 #include "DQM/SiStripCommon/interface/SiStripHistoId.h"
00038 #include "DQM/SiStripCommon/interface/SiStripFolderOrganizer.h"
00039 #include "DQMServices/Core/interface/DQMStore.h"
00040 
00041 #include "TrackingTools/TransientTrackingRecHit/interface/TransientTrackingRecHit.h"
00042 #include "TrackingTools/PatternTools/interface/Trajectory.h"
00043 
00044 #include <TProfile.h>
00045 #include <TStyle.h>
00046 #include <TF1.h>
00047 
00048 #include<list>
00049 
00050 class DetIdLess 
00051   : public std::binary_function< const SiStripRecHit2D*,const SiStripRecHit2D*,bool> {
00052 public:
00053   
00054   DetIdLess() {}
00055   
00056   bool operator()( const SiStripRecHit2D* a, const SiStripRecHit2D* b) const {
00057     return *(a->cluster())<*(b->cluster());
00058   }
00059 };
00060   
00061 
00062   //Constructor
00063 
00064 SiStripLAProfileBooker::SiStripLAProfileBooker(edm::ParameterSet const& conf) : 
00065   conf_(conf)
00066 {
00067 }
00068 
00069   //BeginRun
00070 
00071 void SiStripLAProfileBooker::beginRun(const edm::EventSetup& c){
00072 
00073   //Retrieve tracker topology from geometry
00074   edm::ESHandle<TrackerTopology> tTopoHandle;
00075   c.get<IdealGeometryRecord>().get(tTopoHandle);
00076   const TrackerTopology* const tTopo = tTopoHandle.product();
00077  
00078   //get magnetic field and geometry from ES
00079   edm::ESHandle<MagneticField> esmagfield;
00080   c.get<IdealMagneticFieldRecord>().get(esmagfield);
00081   const MagneticField *  magfield=&(*esmagfield);
00082   
00083   edm::ESHandle<TrackerGeometry> estracker;
00084   c.get<TrackerDigiGeometryRecord>().get(estracker);
00085   tracker=&(*estracker); 
00086 
00087   std::vector<uint32_t> activeDets;
00088   edm::ESHandle<SiStripDetCabling> tkmechstruct=0;
00089   if (conf_.getParameter<bool>("UseStripCablingDB")){ 
00090     c.get<SiStripDetCablingRcd>().get(tkmechstruct);
00091     activeDets.clear();
00092     tkmechstruct->addActiveDetectorsRawIds(activeDets);
00093   }
00094   else {
00095     const TrackerGeometry::DetIdContainer& Id = estracker->detIds();
00096     TrackerGeometry::DetIdContainer::const_iterator Iditer;    
00097     activeDets.clear();
00098     for(Iditer=Id.begin();Iditer!=Id.end();Iditer++){
00099       activeDets.push_back(Iditer->rawId());
00100     }
00101   }
00102    
00103   edm::InputTag TkTag = conf_.getParameter<edm::InputTag>("Tracks");
00104   //Get Ids;
00105   double ModuleRangeMin=conf_.getParameter<double>("ModuleXMin");
00106   double ModuleRangeMax=conf_.getParameter<double>("ModuleXMax");
00107   double TIBRangeMin=conf_.getParameter<double>("TIBXMin");
00108   double TIBRangeMax=conf_.getParameter<double>("TIBXMax");
00109   double TOBRangeMin=conf_.getParameter<double>("TOBXMin");
00110   double TOBRangeMax=conf_.getParameter<double>("TOBXMax");
00111   int TIB_bin=conf_.getParameter<int>("TIB_bin");
00112   int TOB_bin=conf_.getParameter<int>("TOB_bin");
00113   int SUM_bin=conf_.getParameter<int>("SUM_bin");
00114     
00115   hFile = new TFile (conf_.getUntrackedParameter<std::string>("treeName").c_str(), "RECREATE" );
00116   
00117   Hit_Tree = hFile->mkdir("Hit_Tree");
00118   Track_Tree = hFile->mkdir("Track_Tree");
00119   Event_Tree = hFile->mkdir("Event_Tree");
00120   
00121   HitsTree = new TTree("HitsTree", "HitsTree");
00122   
00123   HitsTree->Branch("RunNumber", &RunNumber, "RunNumber/I");
00124   HitsTree->Branch("EventNumber", &EventNumber, "EventNumber/I");
00125   HitsTree->Branch("TanTrackAngle", &TanTrackAngle, "TanTrackAngle/F");
00126   HitsTree->Branch("TanTrackAngleParallel", &TanTrackAngleParallel, "TanTrackAngleParallel/F");
00127   HitsTree->Branch("ClSize", &ClSize, "ClSize/I");
00128   HitsTree->Branch("HitCharge", &HitCharge, "HitCharge/I");
00129   HitsTree->Branch("Hit_Std_Dev", &hit_std_dev, "hit_std_dev/F");
00130   HitsTree->Branch("Type", &Type, "Type/I");
00131   HitsTree->Branch("Layer", &Layer, "Layer/I");
00132   HitsTree->Branch("Wheel", &Wheel, "Wheel/I");
00133   HitsTree->Branch("bw_fw", &bw_fw, "bw_fw/I");
00134   HitsTree->Branch("Ext_Int", &Ext_Int, "Ext_Int/I");
00135   HitsTree->Branch("MonoStereo", &MonoStereo, "MonoStereo/I");
00136   HitsTree->Branch("MagField", &MagField, "MagField/F");
00137   HitsTree->Branch("SignCorrection", &SignCorrection, "SignCorrection/F");
00138   HitsTree->Branch("XGlobal", &XGlobal, "XGlobal/F");
00139   HitsTree->Branch("YGlobal", &YGlobal, "YGlobal/F");
00140   HitsTree->Branch("ZGlobal", &ZGlobal, "ZGlobal/F");
00141   HitsTree->Branch("ParticleCharge", &ParticleCharge, "ParticleCharge/I");
00142   HitsTree->Branch("Momentum", &Momentum, "Momentum/F");
00143   HitsTree->Branch("pt", &pt, "pt/F");
00144   HitsTree->Branch("chi2norm", &chi2norm, "chi2norm/F");
00145   HitsTree->Branch("EtaTrack", &EtaTrack, "EtaTrack/F");
00146   HitsTree->Branch("PhiTrack", &PhiTrack, "PhiTrack/F");
00147   HitsTree->Branch("TrajSize", &trajsize, "trajsize/I");
00148   HitsTree->Branch("HitNr", &HitNr, "HitNr/I");
00149   HitsTree->Branch("HitPerTrack", &HitPerTrack, "HitPerTrack/I");
00150   HitsTree->Branch("id_detector", &id_detector, "id_detector/I");
00151   HitsTree->Branch("thick_detector", &thick_detector, "thick_detector/F");
00152   HitsTree->Branch("pitch_detector", &pitch_detector, "pitch_detector/F");  
00153   HitsTree->Branch("Amplitudes", Amplitudes, "Amplitudes[ClSize]/I");
00154   
00155   HitsTree->SetDirectory(Hit_Tree);
00156   
00157   TrackTree = new TTree("TrackTree", "TrackTree");
00158   
00159   TrackTree->Branch("TrackCounter", &TrackCounter, "TrackCounter/I");
00160   
00161   TrackTree->SetDirectory(Track_Tree);
00162   
00163   EventTree = new TTree("EventTree", "EventTree");
00164   
00165   EventTree->Branch("EventCounter", &EventCounter, "EventCounter/I");
00166   
00167   EventTree->SetDirectory(Event_Tree);
00168   
00169       
00170   // use SistripHistoId for producing histogram id (and title)
00171   SiStripHistoId hidmanager;
00172 
00173   // create SiStripFolderOrganizer
00174   SiStripFolderOrganizer folder_organizer;
00175 
00176   dbe_ = edm::Service<DQMStore>().operator->();
00177   
00178   //get all detids
00179   
00180   for(std::vector<uint32_t>::const_iterator Id = activeDets.begin(); Id!=activeDets.end(); Id++){
00181     
00182     //  for(Iditer=Id.begin();Iditer!=Id.end();Iditer++){ //loop on detids
00183     DetId Iditero=DetId(*Id);
00184     DetId *Iditer=&Iditero;
00185     if((Iditer->subdetId() == int(StripSubdetector::TIB)) || (Iditer->subdetId() == int(StripSubdetector::TOB))){ //include only barrel
00186       
00187       int module_bin = 0;
00188       if(Iditer->subdetId() == int(StripSubdetector::TIB)){
00189         module_bin = TIB_bin;
00190       }else{
00191         module_bin = TOB_bin;
00192       }
00193       
00194       // create a TProfile for each module
00195       StripSubdetector subid(*Iditer);
00196       std::string hid;
00197       //Mono single sided detectors
00198       LocalPoint p;
00199       const GeomDetUnit * stripdet=dynamic_cast<const GeomDetUnit*>(tracker->idToDet(subid));
00200       if(stripdet==0)continue;
00201       const StripTopology& topol=(const StripTopology&)stripdet->topology();
00202       float thickness=stripdet->specificSurface().bounds().thickness();
00203       
00204       folder_organizer.setDetectorFolder(Iditer->rawId(), tTopo);
00205       hid = hidmanager.createHistoId(TkTag.label().c_str(),"det",Iditer->rawId());
00206       MonitorElement * profile=dbe_->bookProfile(hid,hid,module_bin,ModuleRangeMin,ModuleRangeMax,20,0,5,"");
00207       detparameters *param=new detparameters;
00208       histos[Iditer->rawId()] = profile;
00209       detmap[Iditer->rawId()] = param;
00210       param->thickness = thickness*10000;
00211       param->pitch = topol.localPitch(p)*10000;
00212       
00213       const GlobalPoint globalp = (stripdet->surface()).toGlobal(p);
00214       GlobalVector globalmagdir = magfield->inTesla(globalp);
00215       param->magfield=(stripdet->surface()).toLocal(globalmagdir);
00216       
00217       profile->setAxisTitle("tan(#theta_{t})",1);
00218       profile->setAxisTitle("Cluster size",2);
00219       
00220       // create a summary histo if it does not exist already
00221       std::string name;
00222       unsigned int layerid;
00223       getlayer(subid,tTopo,name,layerid);
00224       name+=TkTag.label().c_str();
00225       if(summaryhisto.find(layerid)==(summaryhisto.end())){
00226         folder_organizer.setSiStripFolder();
00227         MonitorElement * summaryprofile=0;
00228         if (subid.subdetId()==int (StripSubdetector::TIB)||subid.subdetId()==int (StripSubdetector::TID))
00229           summaryprofile=dbe_->bookProfile(name,name,SUM_bin,TIBRangeMin,TIBRangeMax,20,0,5,"");
00230         else if (subid.subdetId()==int (StripSubdetector::TOB)||subid.subdetId()==int (StripSubdetector::TEC))
00231           summaryprofile=dbe_->bookProfile(name,name,SUM_bin,TOBRangeMin,TOBRangeMax,20,0,5,"");
00232         if(summaryprofile){
00233           detparameters *summaryparam=new detparameters;
00234           summaryhisto[layerid] = summaryprofile;
00235           summarydetmap[layerid] = summaryparam;
00236           summaryparam->thickness = thickness*10000;
00237           summaryparam->pitch = topol.localPitch(p)*10000;
00238           summaryprofile->setAxisTitle("tan(#theta_{t})",1);
00239           summaryprofile->setAxisTitle("Cluster size",2);
00240         }
00241       }
00242       
00243     } 
00244   } 
00245   
00246   trackcollsize = 0;
00247   trajsize = 0;
00248   RunNumber = 0;
00249   EventNumber = 0;
00250   hitcounter = 0;
00251   hitcounter_2ndloop = 0;
00252   worse_double_hit = 0;
00253   better_double_hit = 0;
00254   eventcounter = 0;
00255   
00256   EventCounter = 1;
00257   TrackCounter = 1;
00258   
00259 }
00260 
00261 SiStripLAProfileBooker::~SiStripLAProfileBooker() {  
00262   detparmap::iterator detpariter;
00263   for( detpariter=detmap.begin(); detpariter!=detmap.end();++detpariter)delete detpariter->second;
00264   for( detpariter=summarydetmap.begin(); detpariter!=summarydetmap.end();++detpariter)delete detpariter->second;
00265   delete hFile;
00266 }  
00267 
00268 // Analyzer: Functions that gets called by framework every event
00269 
00270 void SiStripLAProfileBooker::analyze(const edm::Event& e, const edm::EventSetup& es)
00271 {
00272   //Retrieve tracker topology from geometry
00273   edm::ESHandle<TrackerTopology> tTopoHandle;
00274   es.get<IdealGeometryRecord>().get(tTopoHandle);
00275   const TrackerTopology* const tTopo = tTopoHandle.product();
00276   
00277   RunNumber = e.id().run();
00278   EventNumber = e.id().event();
00279   
00280   eventcounter++;
00281   
00282   EventTree->Fill();
00283   
00284   //Analysis of Trajectory-RecHits
00285   
00286   edm::InputTag TkTag = conf_.getParameter<edm::InputTag>("Tracks");
00287   
00288   edm::Handle<reco::TrackCollection> trackCollection;
00289   e.getByLabel(TkTag,trackCollection);
00290   
00291   edm::Handle<std::vector<Trajectory> > TrajectoryCollection;
00292   e.getByLabel(TkTag,TrajectoryCollection);
00293   
00294   edm::Handle<TrajTrackAssociationCollection> TrajTrackMap;
00295   e.getByLabel(TkTag, TrajTrackMap);
00296   
00297   const reco::TrackCollection *tracks=trackCollection.product();
00298   
00299   // FIXME this has to be changed to use pointers to clusters...
00300   std::map<const SiStripRecHit2D*,std::pair<float,float>,DetIdLess> hitangleassociation;
00301   std::list<SiStripRecHit2D> cache;  // ugly, inefficient, effective in making the above working
00302   
00303   trackcollsize=tracks->size();
00304   trajsize=TrajectoryCollection->size();
00305   
00306   edm::LogInfo("SiStripLAProfileBooker::analyze") <<" Number of tracks in event = "<<trackcollsize<<"\n";
00307   edm::LogInfo("SiStripLAProfileBooker::analyze") <<" Number of trajectories in event = "<<trajsize<<"\n";
00308   
00309   TrajTrackAssociationCollection::const_iterator TrajTrackIter;
00310   
00311   for(TrajTrackIter = TrajTrackMap->begin(); TrajTrackIter!= TrajTrackMap->end(); TrajTrackIter++){ //loop on trajectories
00312     
00313     if(TrajTrackIter->key->foundHits()>=5){
00314       
00315       TrackTree->Fill();
00316       
00317       ParticleCharge = -99;
00318       Momentum = -99;
00319       pt = -99;
00320       chi2norm = -99;
00321       HitPerTrack = -99;
00322       EtaTrack = -99;
00323       PhiTrack = -99;
00324       
00325       ParticleCharge = TrajTrackIter->val->charge();
00326       pt = TrajTrackIter->val->pt();
00327       Momentum = TrajTrackIter->val->p();
00328       chi2norm = TrajTrackIter->val->normalizedChi2();
00329       EtaTrack = TrajTrackIter->val->eta();
00330       PhiTrack = TrajTrackIter->val->phi();
00331       HitPerTrack = TrajTrackIter->key->foundHits();
00332       
00333       std::vector<TrajectoryMeasurement> TMeas=TrajTrackIter->key->measurements();
00334       std::vector<TrajectoryMeasurement>::iterator itm;
00335       
00336       for (itm=TMeas.begin();itm!=TMeas.end();itm++){ //loop on hits
00337         
00338         int i;
00339         for(i=0;i<100;i++){Amplitudes[i]=0;}
00340         
00341         TanTrackAngle = -99;
00342         TanTrackAngleParallel=-99;
00343         ClSize = -99;
00344         HitCharge = 0;
00345         Type = -99;
00346         Layer = -99;
00347         Wheel = -99;
00348         bw_fw = -99;
00349         Ext_Int = -99;
00350         MonoStereo = -99;
00351         MagField = -99;
00352         SignCorrection = -99;
00353         XGlobal = -99;
00354         YGlobal = -99;
00355         ZGlobal = -99;
00356         barycenter = -99;
00357         hit_std_dev = -99;
00358         sumx = 0;
00359         id_detector=-1;
00360         thick_detector=-1;
00361         pitch_detector=-1;       
00362         HitNr = 1;    
00363         
00364         TrajectoryStateOnSurface tsos=itm->updatedState();
00365         const TransientTrackingRecHit::ConstRecHitPointer thit=itm->recHit();
00366         if((thit->geographicalId().subdetId() == int(StripSubdetector::TIB)) ||  thit->geographicalId().subdetId()== int(StripSubdetector::TOB)){ //include only barrel
00367           const SiStripMatchedRecHit2D* matchedhit=dynamic_cast<const SiStripMatchedRecHit2D*>((*thit).hit());
00368           const ProjectedSiStripRecHit2D* phit=dynamic_cast<const ProjectedSiStripRecHit2D*>((*thit).hit());
00369           const SiStripRecHit2D* hit=dynamic_cast<const SiStripRecHit2D*>((*thit).hit());
00370           if(phit) hit=&(phit->originalHit());
00371           
00372           LocalVector trackdirection=tsos.localDirection();
00373           
00374           if(matchedhit){//if matched hit...
00375             
00376             GluedGeomDet * gdet=(GluedGeomDet *)tracker->idToDet(matchedhit->geographicalId());
00377             
00378             GlobalVector gtrkdir=gdet->toGlobal(trackdirection);        
00379             
00380             // THIS THE POINTER TO THE MONO HIT OF A MATCHED HIT 
00381             
00382             // top be migrated to the more direct interface of matchedhit
00383             cache.push_back(matchedhit->monoHit()); 
00384             const SiStripRecHit2D * monohit = &cache.back();   
00385             const SiStripRecHit2D::ClusterRef & monocluster=monohit->cluster();
00386             const GeomDetUnit * monodet=gdet->monoDet();
00387             // this does not exists anymore! either project the matched or use CPE
00388             const LocalPoint monoposition = monohit->localPosition();   
00389             
00390             StripSubdetector detid=(StripSubdetector)monohit->geographicalId();
00391             id_detector = detid.rawId();
00392             thick_detector=monodet->specificSurface().bounds().thickness();
00393             const StripTopology& mtopol=(StripTopology&)monodet->topology();
00394             pitch_detector = mtopol.localPitch(monoposition);
00395             const GlobalPoint monogposition = (monodet->surface()).toGlobal(monoposition);
00396             ClSize = (monocluster->amplitudes()).size();
00397             
00398             const std::vector<uint8_t> amplitudes = monocluster->amplitudes();
00399             
00400             barycenter = monocluster->barycenter()- 0.5; 
00401             uint16_t FirstStrip = monocluster->firstStrip();
00402             std::vector<uint8_t>::const_iterator idigi;
00403             std::vector<uint8_t>::const_iterator begin=amplitudes.begin();
00404             nstrip=0;
00405             for(idigi=begin; idigi!=amplitudes.end(); idigi++){
00406               Amplitudes[nstrip]=*idigi;
00407               sumx+=pow(((FirstStrip+idigi-begin)-barycenter),2)*(*idigi);
00408               HitCharge+=*idigi;
00409             }
00410             hit_std_dev = sqrt(sumx/HitCharge);
00411             
00412             
00413             XGlobal = monogposition.x();
00414             YGlobal = monogposition.y();
00415             ZGlobal = monogposition.z();
00416             
00417             Type = detid.subdetId();
00418             MonoStereo=detid.stereo();
00419             
00420             if(detid.subdetId() == int (StripSubdetector::TIB)){
00421               
00422               Layer = tTopo->tibLayer(detid);
00423               bw_fw = tTopo->tibStringInfo(detid)[0];
00424               Ext_Int = tTopo->tibStringInfo(detid)[1];
00425             }
00426             if(detid.subdetId() == int (StripSubdetector::TOB)){
00427               
00428               Layer = tTopo->tobLayer(detid);
00429               bw_fw = tTopo->tobRodInfo(detid)[0];
00430             }
00431             if(detid.subdetId() == int (StripSubdetector::TID)){
00432               
00433               Wheel = tTopo->tidWheel(detid);
00434               bw_fw = tTopo->tidModuleInfo(detid)[0];
00435             }
00436             if(detid.subdetId() == int (StripSubdetector::TEC)){
00437               
00438               Wheel = tTopo->tecWheel(detid);
00439               bw_fw = tTopo->tecPetalInfo(detid)[0];
00440             }
00441             
00442             
00443             LocalVector monotkdir=monodet->toLocal(gtrkdir);
00444             
00445             if(monotkdir.z()!=0){
00446               
00447               // THE LOCAL ANGLE (MONO)
00448               float tanangle = monotkdir.x()/monotkdir.z();
00449               TanTrackAngleParallel = monotkdir.y()/monotkdir.z();            
00450               TanTrackAngle = tanangle;
00451               detparmap::iterator TheDet=detmap.find(detid.rawId());
00452               LocalVector localmagdir;
00453               if(TheDet!=detmap.end())localmagdir=TheDet->second->magfield;
00454               MagField = localmagdir.mag();
00455               if(MagField != 0.){
00456                 LocalVector monoylocal(0,1,0);
00457                 float signcorrection = (localmagdir * monoylocal)/(MagField);
00458                 if(signcorrection!=0)SignCorrection=1/signcorrection;
00459               }
00460               
00461               std::map<const SiStripRecHit2D *,std::pair<float,float>,DetIdLess>::iterator alreadystored=hitangleassociation.find(monohit);
00462               
00463               
00464               if(alreadystored != hitangleassociation.end()){//decide which hit take
00465                 if(itm->estimate() >  alreadystored->second.first){
00466                   worse_double_hit++;}
00467                 if(itm->estimate() <  alreadystored->second.first){
00468                   better_double_hit++;
00469                   hitangleassociation.insert(std::make_pair(monohit, std::make_pair(itm->estimate(),tanangle)));
00470                   
00471                 }}
00472               else{
00473                 hitangleassociation.insert(make_pair(monohit, std::make_pair(itm->estimate(),tanangle)));
00474                 HitsTree->Fill();
00475                 hitcounter++;}
00476               
00477               // THIS THE POINTER TO THE STEREO HIT OF A MATCHED HIT 
00478               
00479               // top be migrated to the more direct interface of matchedhit
00480               cache.push_back(matchedhit->stereoHit());
00481               const SiStripRecHit2D * stereohit = &cache.back();
00482               const SiStripRecHit2D::ClusterRef & stereocluster=stereohit->cluster();
00483               const GeomDetUnit * stereodet=gdet->stereoDet();
00484               // this does not exists anymore! either project the matched or use CPE
00485               const LocalPoint stereoposition = stereohit->localPosition();    
00486               StripSubdetector detid=(StripSubdetector)stereohit->geographicalId();
00487               id_detector = detid.rawId();
00488               thick_detector=stereodet->specificSurface().bounds().thickness();
00489               const StripTopology& stopol=(StripTopology&)stereodet->topology();
00490               pitch_detector = stopol.localPitch(stereoposition);
00491               const GlobalPoint stereogposition = (stereodet->surface()).toGlobal(stereoposition);
00492               
00493               ClSize = (stereocluster->amplitudes()).size();
00494               
00495               const std::vector<uint8_t> amplitudes = stereocluster->amplitudes();
00496               
00497               barycenter = stereocluster->barycenter()- 0.5; 
00498               uint16_t FirstStrip = stereocluster->firstStrip();
00499               std::vector<uint8_t>::const_iterator idigi;
00500               std::vector<uint8_t>::const_iterator begin=amplitudes.begin();
00501               nstrip=0;
00502               for(idigi=begin; idigi!=amplitudes.end(); idigi++){
00503                 Amplitudes[nstrip]=*idigi;
00504                 sumx+=pow(((FirstStrip+idigi-begin)-barycenter),2)*(*idigi);
00505                 HitCharge+=*idigi;
00506               }
00507               hit_std_dev = sqrt(sumx/HitCharge);
00508               
00509               XGlobal = stereogposition.x();
00510               YGlobal = stereogposition.y();
00511               ZGlobal = stereogposition.z();
00512               
00513               Type = detid.subdetId();
00514               MonoStereo=detid.stereo();
00515               
00516               if(detid.subdetId() == int (StripSubdetector::TIB)){
00517                 
00518                 Layer = tTopo->tibLayer(detid);
00519                 bw_fw = tTopo->tibStringInfo(detid)[0];
00520                 Ext_Int = tTopo->tibStringInfo(detid)[1];
00521               }
00522               if(detid.subdetId() == int (StripSubdetector::TOB)){
00523                 
00524                 Layer = tTopo->tobLayer(detid);
00525                 bw_fw = tTopo->tobRodInfo(detid)[0];
00526               }
00527               if(detid.subdetId() == int (StripSubdetector::TID)){
00528                 
00529                 Wheel = tTopo->tidWheel(detid);
00530                 bw_fw = tTopo->tidModuleInfo(detid)[0];
00531               }
00532               if(detid.subdetId() == int (StripSubdetector::TEC)){
00533                 
00534                 Wheel = tTopo->tecWheel(detid);
00535                 bw_fw = tTopo->tecPetalInfo(detid)[0];
00536               }
00537               
00538               
00539               LocalVector stereotkdir=stereodet->toLocal(gtrkdir);
00540               
00541               if(stereotkdir.z()!=0){
00542                 
00543                 // THE LOCAL ANGLE (STEREO)
00544                 float tanangle = stereotkdir.x()/stereotkdir.z();
00545                 TanTrackAngleParallel = stereotkdir.y()/stereotkdir.z();
00546                 TanTrackAngle = tanangle;
00547                 detparmap::iterator TheDet=detmap.find(detid.rawId());
00548                 LocalVector localmagdir;
00549                 if(TheDet!=detmap.end())localmagdir=TheDet->second->magfield;
00550                 MagField = localmagdir.mag();
00551                 LocalVector stereoylocal(0,1,0);
00552                 if(MagField != 0.){
00553                   float signcorrection = (localmagdir * stereoylocal)/(MagField);
00554                   if(signcorrection!=0)SignCorrection=1/signcorrection;}
00555                 
00556                 std::map<const SiStripRecHit2D *,std::pair<float,float>,DetIdLess>::iterator alreadystored=hitangleassociation.find(stereohit);
00557                 
00558                 if(alreadystored != hitangleassociation.end()){//decide which hit take
00559                   if(itm->estimate() >  alreadystored->second.first){
00560                     worse_double_hit++;}
00561                   if(itm->estimate() <  alreadystored->second.first){
00562                     better_double_hit++;
00563                     hitangleassociation.insert(std::make_pair(stereohit, std::make_pair(itm->estimate(),tanangle)));
00564                     
00565                   }}
00566                 else{
00567                   hitangleassociation.insert(std::make_pair(stereohit, std::make_pair(itm->estimate(),tanangle)));
00568                   HitsTree->Fill();
00569                   hitcounter++;}
00570                 
00571               }
00572             }
00573           }
00574           else if(hit){
00575             
00576             
00577             //  hit= POINTER TO THE RECHIT
00578             
00579             const SiStripRecHit2D::ClusterRef & cluster=hit->cluster();
00580             
00581             GeomDetUnit * gdet=(GeomDetUnit *)tracker->idToDet(hit->geographicalId());
00582             const LocalPoint position = hit->localPosition();    
00583             StripSubdetector detid=(StripSubdetector)hit->geographicalId();
00584             id_detector = detid.rawId();
00585             thick_detector=gdet->specificSurface().bounds().thickness();
00586             const StripTopology& topol=(StripTopology&)gdet->topology();
00587             pitch_detector = topol.localPitch(position);
00588             const GlobalPoint gposition = (gdet->surface()).toGlobal(position);
00589             
00590             ClSize = (cluster->amplitudes()).size();
00591             
00592             const std::vector<uint8_t> amplitudes = cluster->amplitudes();
00593             
00594             barycenter = cluster->barycenter()- 0.5; 
00595             uint16_t FirstStrip = cluster->firstStrip();
00596             std::vector<uint8_t>::const_iterator idigi;
00597             std::vector<uint8_t>::const_iterator begin=amplitudes.begin();
00598             nstrip=0;
00599             for(idigi=begin; idigi!=amplitudes.end(); idigi++){
00600               Amplitudes[nstrip]=*idigi;
00601               sumx+=pow(((FirstStrip+idigi-begin)-barycenter),2)*(*idigi);
00602               HitCharge+=*idigi;
00603             }
00604             hit_std_dev = sqrt(sumx/HitCharge);
00605             
00606             XGlobal = gposition.x();
00607             YGlobal = gposition.y();
00608             ZGlobal = gposition.z();
00609             
00610             Type = detid.subdetId();
00611             MonoStereo=detid.stereo();
00612             
00613             if(detid.subdetId() == int (StripSubdetector::TIB)){
00614               
00615               Layer = tTopo->tibLayer(detid);
00616               bw_fw = tTopo->tibStringInfo(detid)[0];
00617               Ext_Int = tTopo->tibStringInfo(detid)[1];
00618             }
00619             if(detid.subdetId() == int (StripSubdetector::TOB)){
00620               
00621               Layer = tTopo->tobLayer(detid);
00622               bw_fw = tTopo->tobRodInfo(detid)[0];
00623             }
00624             if(detid.subdetId() == int (StripSubdetector::TID)){
00625               
00626               Wheel = tTopo->tidWheel(detid);
00627               bw_fw = tTopo->tidModuleInfo(detid)[0];
00628             }
00629             if(detid.subdetId() == int (StripSubdetector::TEC)){
00630               
00631               Wheel = tTopo->tecWheel(detid);
00632               bw_fw = tTopo->tecPetalInfo(detid)[0];
00633             }
00634             
00635             if(trackdirection.z()!=0){
00636               
00637               // THE LOCAL ANGLE 
00638               float tanangle = trackdirection.x()/trackdirection.z();
00639               TanTrackAngleParallel = trackdirection.y()/trackdirection.z();
00640               TanTrackAngle = tanangle;
00641               detparmap::iterator TheDet=detmap.find(detid.rawId());
00642               LocalVector localmagdir;
00643               if(TheDet!=detmap.end())localmagdir=TheDet->second->magfield;
00644               MagField = localmagdir.mag();
00645               if(MagField != 0.){
00646                 LocalVector ylocal(0,1,0);
00647                 float signcorrection = (localmagdir * ylocal)/(MagField);
00648                 if(signcorrection!=0)SignCorrection=1/signcorrection;}
00649               
00650               std::map<const SiStripRecHit2D *,std::pair<float,float>, DetIdLess>::iterator alreadystored=hitangleassociation.find(hit);
00651               
00652               if(alreadystored != hitangleassociation.end()){//decide which hit take
00653                 if(itm->estimate() >  alreadystored->second.first){
00654                   worse_double_hit++;}
00655                 if(itm->estimate() <  alreadystored->second.first){
00656                   better_double_hit++;
00657                   hitangleassociation.insert(std::make_pair(hit, std::make_pair(itm->estimate(),tanangle)));
00658                   
00659                 }}
00660               else{
00661                 hitangleassociation.insert(std::make_pair(hit,std::make_pair(itm->estimate(), tanangle) ) );
00662                 HitsTree->Fill();
00663                 hitcounter++;}
00664               
00665               
00666             }
00667           }
00668         }
00669       }
00670     }
00671   }
00672   std::map<const SiStripRecHit2D *,std::pair<float,float>,DetIdLess>::iterator hitsiter;
00673   
00674   
00675   for(hitsiter=hitangleassociation.begin();hitsiter!=hitangleassociation.end();hitsiter++){
00676     
00677     hitcounter_2ndloop++;
00678     
00679     const SiStripRecHit2D* hit=hitsiter->first;
00680     const SiStripRecHit2D::ClusterRef & cluster=hit->cluster();
00681     
00682     size=(cluster->amplitudes()).size();
00683     
00684     StripSubdetector detid=(StripSubdetector)hit->geographicalId();  
00685     
00686     float tangent = hitsiter->second.second;
00687     
00688     //Sign and XZ plane projection correction applied in TrackLocalAngle (TIB|TOB layers)
00689     
00690     detparmap::iterator thedet=detmap.find(detid.rawId());
00691     LocalVector localmagdir;
00692     if(thedet!=detmap.end())localmagdir=thedet->second->magfield;
00693     float localmagfield = localmagdir.mag();
00694     
00695     if(localmagfield != 0.){
00696       
00697       LocalVector ylocal(0,1,0);
00698       
00699       float normprojection = (localmagdir * ylocal)/(localmagfield);
00700       
00701       if(normprojection == 0.)LogDebug("SiStripLAProfileBooker::analyze")<<"Error: YBprojection = 0";
00702       
00703       else{
00704         float signprojcorrection = 1/normprojection;
00705         tangent*=signprojcorrection;
00706       }
00707     }
00708     
00709     //Filling histograms
00710     
00711     histomap::iterator thehisto=histos.find(detid.rawId());
00712     
00713     if(thehisto==histos.end())edm::LogError("SiStripLAProfileBooker::analyze")<<"Error: the profile associated to"<<detid.rawId()<<"does not exist! ";    
00714     else thehisto->second->Fill(tangent,size);
00715     
00716     //Summary histograms
00717     std::string name;
00718     unsigned int layerid;
00719     getlayer(detid,tTopo,name,layerid);
00720     histomap::iterator thesummaryhisto=summaryhisto.find(layerid);
00721     if(thesummaryhisto==summaryhisto.end())edm::LogError("SiStripLAProfileBooker::analyze")<<"Error: the profile associated to subdet "<<name<<"does not exist! ";   
00722     else thesummaryhisto->second->Fill(tangent,size);
00723     
00724   }
00725   
00726   
00727 }
00728 
00729 
00730 //Makename function
00731  
00732 void SiStripLAProfileBooker::getlayer(const DetId & detid, const TrackerTopology* tTopo, std::string &name,unsigned int &layerid){
00733     int layer=0;
00734     std::stringstream layernum;
00735 
00736     if(detid.subdetId() == int (StripSubdetector::TIB)){
00737       
00738       name+="TIB_Layer_";
00739       layer = tTopo->tibLayer(detid);
00740     }
00741 
00742     else if(detid.subdetId() == int (StripSubdetector::TID)){
00743       
00744       name+="TID_Ring_";
00745       layer = tTopo->tidRing(detid);
00746     }
00747 
00748     else if(detid.subdetId() == int (StripSubdetector::TOB)){
00749       
00750       name+="TOB_Layer_";
00751       layer = tTopo->tobLayer(detid);
00752 
00753     }
00754 
00755     else if(detid.subdetId() == int (StripSubdetector::TEC)){
00756       
00757       name+="TEC_Ring_";
00758       layer = tTopo->tecRing(detid);
00759     }
00760     layernum<<layer;
00761     name+=layernum.str();
00762     layerid=detid.subdetId()*10+layer;
00763   
00764 }
00765  
00766 void SiStripLAProfileBooker::endJob(){
00767 
00768   std::string outputFile_ =conf_.getUntrackedParameter<std::string>("fileName", "LorentzAngle.root");
00769   dbe_->save(outputFile_);
00770   
00771   hFile->Write();
00772   hFile->Close();
00773 }