CMS 3D CMS Logo

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