CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_2_9_HLT1_bphpatch4/src/CalibTracker/SiStripLorentzAngle/plugins/SiStripLAProfileBooker.cc

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