CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/FastSimulation/CaloHitMakers/src/EcalHitMaker.cc

Go to the documentation of this file.
00001 #include "DataFormats/DetId/interface/DetId.h"
00002 #include "DataFormats/EcalDetId/interface/EcalSubdetector.h"
00003 //#include "DataFormats/EcalDetId/interface/EBDetId.h"
00004 //#include "DataFormats/EcalDetId/interface/EEDetId.h"
00005 
00006 //#include "Geometry/EcalAlgo/interface/EcalBarrelGeometry.h"
00007 //#include "Geometry/EcalAlgo/interface/EcalEndcapGeometry.h"
00008 //#include "Geometry/CaloGeometry/interface/CaloCellGeometry.h"
00009 
00010 #include "FastSimulation/CaloHitMakers/interface/EcalHitMaker.h"
00011 #include "FastSimulation/CaloGeometryTools/interface/CaloGeometryHelper.h"
00012 #include "FastSimulation/CaloGeometryTools/interface/CrystalWindowMap.h"
00013 #include "FastSimulation/CaloGeometryTools/interface/CaloDirectionOperations.h"
00014 #include "FastSimulation/Event/interface/FSimVertex.h"
00015 #include "FastSimulation/Event/interface/FSimTrack.h"
00016 #include "FastSimulation/CalorimeterProperties/interface/PreshowerLayer1Properties.h"
00017 #include "FastSimulation/CalorimeterProperties/interface/PreshowerLayer2Properties.h"
00018 #include "FastSimulation/CalorimeterProperties/interface/HCALProperties.h"
00019 #include "FastSimulation/Utilities/interface/RandomEngine.h"
00020 //#include "FastSimulation/Utilities/interface/Histos.h"
00021 
00022 //#include "Math/GenVector/Transform3D.h"
00023 #include "FastSimulation/CaloGeometryTools/interface/Transform3DPJ.h"
00024 
00025 #include <algorithm>
00026 #include <cmath>
00027 
00028 typedef ROOT::Math::Plane3D::Vector Vector;
00029 typedef ROOT::Math::Plane3D::Point Point;
00030 typedef ROOT::Math::Transform3DPJ Transform3DR;
00031 
00032 EcalHitMaker::EcalHitMaker(CaloGeometryHelper * theCalo,
00033                            const XYZPoint& ecalentrance, 
00034                            const DetId& cell, int onEcal,
00035                            unsigned size, unsigned showertype,
00036                            const RandomEngine* engine):
00037   CaloHitMaker(theCalo,DetId::Ecal,((onEcal==1)?EcalBarrel:EcalEndcap),onEcal,showertype),
00038   EcalEntrance_(ecalentrance),
00039   onEcal_(onEcal),
00040   myTrack_(NULL),
00041   random(engine)
00042 {
00043 #ifdef FAMOSDEBUG
00044   myHistos = Histos::instance();
00045 #endif
00046   //  myHistos->debug("Constructeur EcalHitMaker");
00047   simulatePreshower_ = true;
00048   X0depthoffset_ = 0. ;
00049   X0PS1_ = 0.;
00050   X0PS2_ = 0.; 
00051   X0PS2EE_ = 0.;
00052   X0ECAL_ = 0.;
00053   X0EHGAP_ = 0.;
00054   X0HCAL_ = 0.;
00055   L0PS1_ = 0.;
00056   L0PS2_ = 0.;
00057   L0PS2EE_ = 0.;
00058   L0ECAL_ = 0.;
00059   L0EHGAP_ = 0.;
00060   L0HCAL_ = 0.;
00061   maxX0_ = 0.;
00062   totalX0_ = 0;
00063   totalL0_ = 0.;
00064   pulledPadProbability_ = 1.;
00065   outsideWindowEnergy_ = 0.;
00066   rearleakage_ = 0.;
00067   bfactor_ = 1.; 
00068   ncrystals_ = 0;
00069 
00070   doreorg_ = !showertype;
00071 
00072   hitmaphasbeencalculated_ = false;
00073 
00074   if(onEcal) 
00075     myCalorimeter->buildCrystal(cell,pivot_);
00076   else
00077     pivot_=Crystal();
00078   central_=onEcal==1;
00079   ecalFirstSegment_=-1;
00080   
00081   myCrystalWindowMap_ = 0; 
00082   // In some cases, a "dummy" grid, not based on a cell, can be built. The previous variables
00083   // should however be initialized. In such a case onEcal=0
00084   if(!onEcal) return;
00085 
00086   // Same size in eta-phi
00087   etasize_ = size;
00088   phisize_ = size;
00089 
00090   // Build the grid
00091   // The result is put in CellsWindow and is ordered by distance to the pivot
00092   myCalorimeter->getWindow(pivot_.getDetId(),size,size,CellsWindow_);
00093 
00094   buildGeometry();
00095   //  std::cout << " Geometry built " << regionOfInterest_.size() << std::endl;
00096 
00097   truncatedGrid_ = CellsWindow_.size()!=(etasize_*phisize_);
00098 
00099   // A local vector of corners
00100   mycorners.resize(4);
00101   corners.resize(4);
00102   
00103 #ifdef DEBUGGW
00104   myHistos->fill("h10",EcalEntrance_.eta(),CellsWindow_.size());
00105   if(onEcal==2) 
00106     {
00107       myHistos->fill("h20",EcalEntrance_.perp(),CellsWindow_.size());
00108       if(EcalEntrance_.perp()>70&&EcalEntrance_.perp()<80&&CellsWindow_.size()<35)
00109         {
00110           std::cout << " Truncated grid " << CellsWindow_.size() << " " << EcalEntrance_.perp() << std::endl;
00111           std::cout << " Pivot " << myCalorimeter->getEcalEndcapGeometry()->getGeometry(pivot_.getDetId())->getPosition().perp();
00112           std::cout << EEDetId(pivot_.getDetId()) << std::endl;
00113 
00114           std::cout << " Test getClosestCell " << EcalEntrance_ << std::endl;
00115           DetId testcell = myCalorimeter->getClosestCell(EcalEntrance_, true, false);
00116           std::cout << " Result "<< EEDetId(testcell) << std::endl;
00117           std::cout << " Position " << myCalorimeter->getEcalEndcapGeometry()->getGeometry(testcell)->getPosition() << std::endl;
00118         }
00119     }
00120 
00121 #endif
00122 
00123 }
00124 
00125 EcalHitMaker::~EcalHitMaker()
00126 {
00127   if (myCrystalWindowMap_ != 0)
00128     {
00129       delete myCrystalWindowMap_;
00130     }
00131 }
00132 
00133 bool 
00134 EcalHitMaker::addHitDepth(double r,double phi,double depth)
00135 {
00136 //  std::cout << " Add hit depth called; Current deph is "  << currentdepth_;
00137 //  std::cout << " Required depth is " << depth << std::endl;
00138   depth+=X0depthoffset_;
00139   double sp(1.);
00140   r*=radiusFactor_;
00141   CLHEP::Hep2Vector point(r*std::cos(phi),r*std::sin(phi));
00142 
00143   unsigned xtal=fastInsideCell(point,sp);  
00144   //  if(cellid.isZero()) std::cout << " cell is Zero " << std::endl;
00145 //  if(xtal<1000) 
00146 //    {
00147 //      std::cout << "Result " << regionOfInterest_[xtal].getX0Back() << " " ;
00148 //      std::cout << depth << std::endl;
00149 //    }
00150 //  myHistos->fill("h5000",depth);
00151   if(xtal<1000)
00152     {
00153       //      myHistos->fill("h5002",regionOfInterest_[xtal].getX0Back(),depth);
00154       //      myHistos->fill("h5003",ecalentrance_.eta(),maxX0_);
00155       if(regionOfInterest_[xtal].getX0Back()>depth) 
00156         {
00157           hits_[xtal]+=spotEnergy;      
00158           //      myHistos->fill("h5005",r);
00159           return true;
00160         }
00161       else
00162         {
00163           rearleakage_+=spotEnergy;
00164         }
00165     }
00166   else
00167     {
00168 //      std::cout << " Return false " << std::endl;
00169 //      std::cout << " Add hit depth called; Current deph is "  << currentdepth_;
00170 //      std::cout << " Required depth is " << depth << std::endl;
00171 //      std::cout << " R = " << r << " " << radiusFactor_ << std::endl; 
00172     }
00173 
00174   outsideWindowEnergy_+=spotEnergy;
00175   return false;
00176 }
00177 
00178 
00179 //bool EcalHitMaker::addHitDepth(double r,double phi,double depth)
00180 //{
00181 //  std::map<unsigned,float>::iterator itcheck=hitMap_.find(pivot_.getDetId().rawId());
00182 //  if(itcheck==hitMap_.end())
00183 //    {
00184 //      hitMap_.insert(std::pair<uint32_t,float>(pivot_.getDetId().rawId(),spotEnergy));
00185 //    }
00186 //  else
00187 //    {
00188 //      itcheck->second+=spotEnergy;
00189 //    }
00190 //  return true;
00191 //}
00192 
00193 bool 
00194 EcalHitMaker::addHit(double r,double phi,unsigned layer)
00195 {
00196   //  std::cout <<" Addhit " << std::endl;
00197   //  std::cout << " Before insideCell " << std::endl;
00198   double sp(1.);
00199   //  std::cout << " Trying to add " << r << " " << phi << " " << radiusFactor_ << std::endl; 
00200   r*=radiusFactor_;
00201   CLHEP::Hep2Vector point(r*std::cos(phi),r*std::sin(phi));
00202   //  std::cout << "point " << point << std::endl;
00203   //  CellID cellid=insideCell(point,sp);
00204   unsigned xtal=fastInsideCell(point,sp);  
00205   //  if(cellid.isZero()) std::cout << " cell is Zero " << std::endl;
00206   if(xtal<1000) 
00207     {
00208       if(sp==1.)
00209         hits_[xtal]+=spotEnergy;        
00210       else      
00211         hits_[xtal]+=(random->flatShoot()<sp)*spotEnergy;
00212       return true;
00213     }
00214 
00215   outsideWindowEnergy_+=spotEnergy;
00216 //  std::cout << " This hit ; r= " << point << " hasn't been added "<<std::endl;
00217 //  std::cout << " Xtal " << xtal << std::endl;
00218 //  for(unsigned ip=0;ip<npadsatdepth_;++ip)
00219 //    {
00220 //      std::cout << padsatdepth_[ip] << std::endl;
00221 //    }
00222   
00223   return false;
00224 }
00225 
00226 // Temporary solution
00227 //bool EcalHitMaker::addHit(double r,double phi,unsigned layer)
00228 //{
00229 //  std::map<unsigned,float>::iterator itcheck=hitMap_.find(pivot_.getDetId().rawId());
00230 //  if(itcheck==hitMap_.end())
00231 //    {
00232 //      hitMap_.insert(std::pair<uint32_t,float>(pivot_.getDetId().rawId(),spotEnergy));
00233 //    }
00234 //  else
00235 //    {
00236 //      itcheck->second+=spotEnergy;
00237 //    }
00238 //  return true;
00239 //}
00240 
00241 
00242 unsigned 
00243 EcalHitMaker::fastInsideCell(const CLHEP::Hep2Vector & point,double & sp,bool debug) 
00244 {
00245 
00246   //  debug = true;
00247   bool found=false;
00248   unsigned niter=0;
00249   // something clever has to be implemented here
00250   unsigned d1,d2;
00251   convertIntegerCoordinates(point.x(),point.y(),d1,d2);
00252   //  std::cout << "Fastinside cell " << point.x() << " " <<  point.y() << " " << d1 << " "<< d2 << " " << nx_ << " " << ny_ << std::endl;
00253   if(d1>=nx_||d2>=ny_) 
00254     {
00255       //      std::cout << " Not in the map " <<std::endl;
00256       return 9999;
00257     }
00258   unsigned cell=myCrystalNumberArray_[d1][d2];
00259   // We are likely to be lucky
00260   //  std::cout << " Got the cell " << cell << std::endl;
00261   if (validPads_[cell]&&padsatdepth_[cell].inside(point)) 
00262     {
00263       //      std::cout << " We are lucky " << cell << std::endl;
00264       sp = padsatdepth_[cell].survivalProbability();
00265       return cell;
00266     }
00267   
00268   //  std::cout << "Starting the loop " << std::endl;
00269   bool status(true);
00270   const std::vector<unsigned>& localCellVector(myCrystalWindowMap_->getCrystalWindow(cell,status));
00271   if(status) 
00272     {
00273       unsigned size=localCellVector.size();      
00274       //      std::cout << " Starting from " << EBDetId(regionOfInterest_[cell].getDetId()) << std::endl;
00275       //      const std::vector<DetId>& neighbours=myCalorimeter->getNeighbours(regionOfInterest_[cell].getDetId());
00276 //      std::cout << " The neighbours are " << std::endl;
00277 //      for(unsigned ic=0;ic<neighbours.size(); ++ic)
00278 //      {
00279 //        std::cout << EBDetId(neighbours[ic]) << std::endl;
00280 //      }
00281 //      std::cout << " Done " << std::endl;
00282       for(unsigned ic=0;ic<8&&ic<size;++ic)
00283         {
00284           unsigned iq=localCellVector[ic];
00285 //        std::cout << " Testing " << EBDetId(regionOfInterest_[iq].getDetId()) << std::endl; ; 
00286 //        std::cout << " " << iq << std::endl;
00287 //        std::cout << padsatdepth_[iq] ;
00288           if(validPads_[iq]&&padsatdepth_[iq].inside(point))
00289             {
00290               //              std::cout << " Yes " << std::endl;
00291               //              myHistos->fill("h1000",niter);
00292               sp = padsatdepth_[iq].survivalProbability();
00293               //              std::cout << "Finished the loop " << niter << std::endl;
00294               //              std::cout << "Inside " << std::endl;
00295               return iq;
00296             }
00297           //      std::cout << " Not inside " << std::endl;
00298           //      std::cout << "No " << std::endl;
00299           ++niter;
00300         }
00301           }
00302   if(debug) std::cout << " not found in a quad, let's check the " << ncrackpadsatdepth_ << " cracks " << std::endl;
00303   //  std::cout << "Finished the loop " << niter << std::endl;
00304   // Let's check the cracks 
00305   //  std::cout << " Let's check the cracks " << ncrackpadsatdepth_ << " " << crackpadsatdepth_.size() << std::endl;
00306   unsigned iquad=0;
00307   unsigned iquadinside=999;
00308   while(iquad<ncrackpadsatdepth_&&!found)
00309     {
00310       //      std::cout << " Inside the while " << std::endl;
00311       if(crackpadsatdepth_[iquad].inside(point)) 
00312         {
00313           iquadinside=iquad;
00314           found=true;
00315           sp = crackpadsatdepth_[iquad].survivalProbability();
00316         }
00317       ++iquad;
00318       ++niter;
00319     }
00320   //  myHistos->fill("h1002",niter);
00321   if(!found&&debug) std::cout << " Not found in the cracks " << std::endl;
00322   return (found) ? crackpadsatdepth_[iquadinside].getNumber(): 9999;
00323 }
00324 
00325 
00326 void 
00327 EcalHitMaker::setTrackParameters(const XYZNormal& normal,
00328                                  double X0depthoffset,
00329                                  const FSimTrack& theTrack)
00330 {
00331   //  myHistos->debug("setTrackParameters");
00332   //  std::cout << " Track " << theTrack << std::endl;
00333   intersections_.clear();
00334   // This is certainly enough
00335   intersections_.reserve(50);
00336   myTrack_=&theTrack;
00337   normal_=normal.Unit();
00338   X0depthoffset_=X0depthoffset;
00339   cellLine(intersections_);
00340   buildSegments(intersections_);
00341 //  std::cout << " Segments " << segments_.size() << std::endl;
00342 //  for(unsigned ii=0; ii<segments_.size() ; ++ii)
00343 //    {
00344 //      std::cout << segments_[ii] << std::endl;
00345 //    }
00346 
00347 
00348   // This is only needed in case of electromagnetic showers 
00349   if(EMSHOWER&&onEcal_&&ecalTotalX0()>0.)
00350     {
00351       //      std::cout << "Total X0  " << ecalTotalX0() << std::endl;
00352       for(unsigned ic=0;ic<ncrystals_;++ic)
00353         {
00354           for(unsigned idir=0;idir<4;++idir)
00355             {
00356               XYZVector norm=regionOfInterest_[ic].exitingNormal(CaloDirectionOperations::Side(idir));
00357               regionOfInterest_[ic].crystalNeighbour(idir).setToBeGlued((norm.Dot(normal_)<0.));
00358             }
00359           // Now calculate the distance in X0 of the back sides of the crystals
00360           // (only for EM showers)
00361           if(EMSHOWER)
00362             {
00363               XYZVector dir=regionOfInterest_[ic].getBackCenter()-segments_[ecalFirstSegment_].entrance();
00364               double dist=dir.Dot(normal_);
00365               double absciss= dist+segments_[ecalFirstSegment_].sEntrance();
00366               std::vector<CaloSegment>::const_iterator segiterator;
00367               // First identify the correct segment
00368               //      std::cout << " Crystal " << ic  << regionOfInterest_[ic].getBackCenter() ;
00369               //      std::cout << " Entrance : " << segments_[ecalFirstSegment_].entrance()<< std::endl;
00370               //      std::cout << " Looking for the segment " << dist << std::endl;
00371               segiterator = find_if(segments_.begin(),segments_.end(),CaloSegment::inSegment(absciss));
00372               //      std::cout << " Done " << std::endl;
00373               if(segiterator==segments_.end() )
00374                 {
00375                   // in this case, we won't have any problem. No need to 
00376                   // calculate the real depth. 
00377                   regionOfInterest_[ic].setX0Back(9999);
00378                 }
00379               else
00380                 {
00381                   DetId::Detector det(segiterator->whichDetector());
00382                   if(det!=DetId::Ecal)
00383                     {
00384                       regionOfInterest_[ic].setX0Back(9999);
00385                     }
00386                   else
00387                     {
00388                       double x0=segiterator->x0FromCm(dist);
00389                       if(x0<maxX0_) maxX0_=x0;
00390                       regionOfInterest_[ic].setX0Back(x0);
00391                     }
00392                 }
00393               //  myHistos->fill("h4000",ecalentrance_.eta(), regionOfInterest_[ic].getX0Back());
00394               //}
00395               //      else
00396               //{
00397               //   // in this case, we won't have any problem. No need to 
00398               //  // calculate the real depth. 
00399               //  regionOfInterest_[ic].setX0Back(9999);
00400               //}
00401             }//EMSHOWER  
00402         } // ndir
00403       //      myHistos->fill("h6000",segments_[ecalFirstSegment_].entrance().eta(),maxX0_);
00404     }
00405   //  std::cout << "Leaving setTrackParameters" << std::endl
00406 }
00407 
00408 
00409 void 
00410 EcalHitMaker::cellLine(std::vector<CaloPoint>& cp) 
00411 {
00412   cp.clear();
00413   //  if(myTrack->onVFcal()!=2)
00414   //    {
00415   if(!central_&&onEcal_&&simulatePreshower_) preshowerCellLine(cp);
00416   if(onEcal_)ecalCellLine(EcalEntrance_,EcalEntrance_+normal_,cp);
00417   //    }
00418   
00419   XYZPoint vertex(myTrack_->vertex().position().Vect());
00420 
00421   //sort the points by distance (in the ECAL they are not necessarily ordered)
00422   XYZVector dir(0.,0.,0.);
00423   if(myTrack_->onLayer1())
00424     {
00425       vertex=(myTrack_->layer1Entrance().vertex()).Vect();
00426       dir=myTrack_->layer1Entrance().Vect().Unit();
00427     }
00428   else if(myTrack_->onLayer2())
00429     {
00430       vertex=(myTrack_->layer2Entrance().vertex()).Vect();
00431       dir=myTrack_->layer2Entrance().Vect().Unit();
00432     }  
00433   else if(myTrack_->onEcal())
00434     {
00435       vertex=(myTrack_->ecalEntrance().vertex()).Vect();
00436       dir=myTrack_->ecalEntrance().Vect().Unit();
00437     }
00438   else if(myTrack_->onHcal())
00439     {
00440       vertex=(myTrack_->hcalEntrance().vertex()).Vect();
00441       dir=myTrack_->hcalEntrance().Vect().Unit();
00442     }
00443   else if(myTrack_->onVFcal()==2)
00444     {
00445       vertex=(myTrack_->vfcalEntrance().vertex()).Vect();
00446       dir=myTrack_->vfcalEntrance().Vect().Unit();
00447     }
00448   else
00449     {
00450       std::cout << " Problem with the grid " << std::endl;
00451     }
00452  
00453   // Move the vertex for distance comparison (5cm)
00454   vertex -= 5.*dir;
00455   CaloPoint::DistanceToVertex myDistance(vertex);
00456   sort(cp.begin(),cp.end(),myDistance);
00457   
00458   // The intersections with the HCAL shouldn't need to be sorted
00459   // with the N.I it is actually a source of problems
00460   hcalCellLine(cp);
00461 
00462 //  std::cout << " Intersections ordered by distance to " << vertex << std::endl;
00463 //
00464 //  for (unsigned ic=0;ic<cp.size();++ic)
00465 //    {
00466 //      XYZVector t=cp[ic]-vertex;
00467 //      std::cout << cp[ic] << " " << t.mag() << std::endl;
00468 //    }
00469 }
00470 
00471 
00472 void 
00473 EcalHitMaker::preshowerCellLine(std::vector<CaloPoint>& cp) const
00474 {
00475   //  FSimEvent& mySimEvent = myEventMgr->simSignal();
00476   //  FSimTrack myTrack = mySimEvent.track(fsimtrack_);
00477   //  std::cout << "FsimTrack " << fsimtrack_<< std::endl;
00478   //  std::cout << " On layer 1 " << myTrack.onLayer1() << std::endl;
00479   // std::cout << " preshowerCellLine " << std::endl; 
00480   if(myTrack_->onLayer1())
00481     {
00482       XYZPoint point1=(myTrack_->layer1Entrance().vertex()).Vect();
00483       double phys_eta=myTrack_->layer1Entrance().eta();
00484       double cmthickness = 
00485         myCalorimeter->layer1Properties(1)->thickness(phys_eta);
00486 
00487       if(cmthickness>0)
00488         {
00489           XYZVector dir=myTrack_->layer1Entrance().Vect().Unit();
00490           XYZPoint point2=point1+dir*cmthickness;
00491           
00492           CaloPoint cp1(DetId::Ecal,EcalPreshower,1,point1);
00493           CaloPoint cp2(DetId::Ecal,EcalPreshower,1,point2);
00494           cp.push_back(cp1);
00495           cp.push_back(cp2);
00496         }
00497       else
00498         {
00499           //      std::cout << "Track on ECAL " << myTrack.EcalEntrance_().vertex()*0.1<< std::endl;
00500         }
00501     }
00502 
00503   //  std::cout << " On layer 2 " << myTrack.onLayer2() << std::endl;
00504   if(myTrack_->onLayer2())
00505     {
00506       XYZPoint point1=(myTrack_->layer2Entrance().vertex()).Vect();
00507       double phys_eta=myTrack_->layer2Entrance().eta();
00508       double cmthickness = 
00509         myCalorimeter->layer2Properties(1)->thickness(phys_eta);
00510       if(cmthickness>0)
00511         {
00512           XYZVector dir=myTrack_->layer2Entrance().Vect().Unit();
00513           XYZPoint point2=point1+dir*cmthickness;
00514           
00515 
00516           CaloPoint cp1(DetId::Ecal,EcalPreshower,2,point1);
00517           CaloPoint cp2(DetId::Ecal,EcalPreshower,2,point2);
00518 
00519           cp.push_back(cp1);
00520           cp.push_back(cp2);
00521         }
00522       else
00523         {
00524           //      std::cout << "Track on ECAL " << myTrack.EcalEntrance_().vertex()*0.1 << std::endl;
00525         }
00526     }
00527   //  std::cout << " Exit preshower CellLine " << std::endl;
00528 }
00529 
00530 void 
00531 EcalHitMaker::hcalCellLine(std::vector<CaloPoint>& cp) const
00532 {
00533   //  FSimEvent& mySimEvent = myEventMgr->simSignal();
00534   //  FSimTrack myTrack = mySimEvent.track(fsimtrack_);
00535   int onHcal=myTrack_->onHcal();
00536 
00537   if(onHcal<=2&&onHcal>0)
00538     {
00539       XYZPoint point1=(myTrack_->hcalEntrance().vertex()).Vect();
00540 
00541       double eta=point1.eta();
00542       // HCAL thickness in cm (assuming that the particle is coming from 000)
00543       double thickness= myCalorimeter->hcalProperties(onHcal)->thickness(eta);
00544       cp.push_back(CaloPoint(DetId::Hcal,point1));
00545       XYZVector dir=myTrack_->hcalEntrance().Vect().Unit();
00546       XYZPoint point2=point1+dir*thickness;
00547 
00548       cp.push_back(CaloPoint(DetId::Hcal,point2));
00549       
00550     }
00551   int onVFcal=myTrack_->onVFcal();
00552   if(onVFcal==2)
00553     {
00554       XYZPoint point1=(myTrack_->vfcalEntrance().vertex()).Vect();
00555       double eta=point1.eta();
00556       // HCAL thickness in cm (assuming that the particle is coming from 000)
00557       double thickness= myCalorimeter->hcalProperties(3)->thickness(eta);
00558       cp.push_back(CaloPoint(DetId::Hcal,point1));
00559       XYZVector dir=myTrack_->vfcalEntrance().Vect().Unit();
00560       if(thickness>0)
00561         {
00562           XYZPoint point2=point1+dir*thickness;
00563           cp.push_back(CaloPoint(DetId::Hcal,point2));
00564         }
00565 
00566     }
00567 }
00568 
00569 void 
00570 EcalHitMaker::ecalCellLine(const XYZPoint& a,const XYZPoint& b,std::vector<CaloPoint>& cp) 
00571 {
00572   //  std::vector<XYZPoint> corners;
00573   //  corners.resize(4);
00574   unsigned ic=0;
00575   double t;
00576   XYZPoint xp;
00577   DetId c_entrance,c_exit;
00578   bool entrancefound(false),exitfound(false);
00579   //  std::cout << " Look for intersections " << ncrystals_ << std::endl;
00580   //  std::cout << " regionOfInterest_ " << truncatedGrid_ << " " << regionOfInterest_.size() << std::endl;
00581   // try to determine the number of crystals to test
00582   // First determine the incident angle
00583   double angle=std::acos(normal_.Dot(regionOfInterest_[0].getAxis().Unit()));
00584 
00585   //  std::cout << " Normal " << normal_<< " Axis " << regionOfInterest_[0].getAxis().Unit() << std::endl;
00586   double backdistance=std::sqrt(regionOfInterest_[0].getAxis().mag2())*std::tan(angle);
00587   // 1/2.2cm = 0.45
00588 //   std::cout << " Angle " << angle << std::endl;
00589 //   std::cout << " Back distance " << backdistance << std::endl;
00590   unsigned ncrystals=(unsigned)(backdistance*0.45);
00591   unsigned highlim=(ncrystals+4);
00592   highlim*=highlim;
00593   if(highlim>ncrystals_) highlim=ncrystals_;
00594   //  unsigned lowlim=(ncrystals>2)? (ncrystals-2):0;
00595   //  std::cout << " Ncrys " << ncrystals << std::endl;
00596   
00597   
00598   while(ic<ncrystals_&&(ic<highlim||!exitfound))
00599     {
00600       // Check front side
00601       //      if(!entrancefound)
00602         {
00603           const Plane3D& plan=regionOfInterest_[ic].getFrontPlane();
00604 //        XYZVector axis1=(plan.Normal());
00605 //        XYZVector axis2=regionOfInterest_[ic].getFirstEdge();
00606           xp=intersect(plan,a,b,t,false);
00607           regionOfInterest_[ic].getFrontSide(corners);
00608           //      CrystalPad pad(9999,onEcal_,corners,regionOfInterest_[ic].getCorner(0),axis1,axis2);
00609           //      if(pad.globalinside(xp)) 
00610           if(inside3D(corners,xp))
00611             {
00612               cp.push_back(CaloPoint(regionOfInterest_[ic].getDetId(),UP,xp));
00613               entrancefound=true;
00614               c_entrance=regionOfInterest_[ic].getDetId();
00615               //      myHistos->fill("j12",highlim,ic);
00616             }
00617         }
00618       
00619       // check rear side
00620         //      if(!exitfound)
00621         {
00622           const Plane3D& plan=regionOfInterest_[ic].getBackPlane();
00623 //        XYZVector axis1=(plan.Normal());
00624 //        XYZVector axis2=regionOfInterest_[ic].getFifthEdge();
00625           xp=intersect(plan,a,b,t,false);
00626           regionOfInterest_[ic].getBackSide(corners);
00627           //      CrystalPad pad(9999,onEcal_,corners,regionOfInterest_[ic].getCorner(4),axis1,axis2);
00628           //      if(pad.globalinside(xp)) 
00629           if(inside3D(corners,xp))
00630             {
00631               cp.push_back(CaloPoint(regionOfInterest_[ic].getDetId(),DOWN,xp));
00632               exitfound=true;
00633               c_exit=regionOfInterest_[ic].getDetId();
00634               //              std::cout << " Crystal : " << ic << std::endl;
00635               //              myHistos->fill("j10",highlim,ic);
00636             }
00637         }
00638       
00639       if(entrancefound&&exitfound&&c_entrance==c_exit) return;
00640       // check lateral sides 
00641       for(unsigned iside=0;iside<4;++iside)
00642         {
00643           const Plane3D& plan=regionOfInterest_[ic].getLateralPlane(iside);
00644           xp=intersect(plan,a,b,t,false);
00645 //        XYZVector axis1=(plan.Normal());
00646 //        XYZVector axis2=regionOfInterest_[ic].getLateralEdge(iside);
00647           regionOfInterest_[ic].getLateralSide(iside,corners);
00648           //      CrystalPad pad(9999,onEcal_,corners,regionOfInterest_[ic].getCorner(iside),axis1,axis2);
00649           //      if(pad.globalinside(xp)) 
00650           if(inside3D(corners,xp))
00651             {
00652               cp.push_back(CaloPoint(regionOfInterest_[ic].getDetId(),CaloDirectionOperations::Side(iside),xp)); 
00653               //              std::cout << cp[cp.size()-1] << std::endl;
00654             }     
00655         }
00656       // Go to next crystal 
00657       ++ic;
00658     }    
00659 }
00660 
00661 
00662 void 
00663 EcalHitMaker::buildSegments(const std::vector<CaloPoint>& cp)
00664 {
00665   //  myHistos->debug();
00666   //  TimeMe theT("FamosGrid::buildSegments");
00667   unsigned size=cp.size();
00668   if(size%2!=0) 
00669     {
00670       //      std::cout << " There is a problem " << std::endl;
00671       return;
00672     }
00673   //  myHistos->debug();
00674   unsigned nsegments=size/2;
00675   segments_.reserve(nsegments);
00676   if (size==0) return;
00677   // curv abs
00678   double s=0.;
00679   double sX0=0.;
00680   double sL0=0.;
00681   
00682   unsigned ncrossedxtals = 0;
00683   unsigned is=0;
00684   while(is<nsegments)
00685     {
00686 
00687       if(cp[2*is].getDetId()!=cp[2*is+1].getDetId()&&
00688          cp[2*is].whichDetector()!=DetId::Hcal&&
00689          cp[2*is+1].whichDetector()!=DetId::Hcal) 
00690         {
00691 //        std::cout << " Problem with the segments " << std::endl;
00692 //        std::cout << cp[2*is].whichDetector() << " " << cp[2*is+1].whichDetector() << std::endl;
00693 //        std::cout << is << " " <<cp[2*is].getDetId().rawId() << std::endl; 
00694 //        std::cout << (2*is+1) << " " <<cp[2*is+1].getDetId().rawId() << std::endl; 
00695           ++is;
00696           continue;
00697         }
00698 
00699       // Check if it is a Preshower segment - Layer 1 
00700       // One segment per layer, nothing between
00701       //      myHistos->debug("Just avant Preshower"); 
00702       if(cp[2*is].whichDetector()==DetId::Ecal && cp[2*is].whichSubDetector()==EcalPreshower && cp[2*is].whichLayer()==1)
00703         {
00704           if(cp[2*is+1].whichDetector()==DetId::Ecal && cp[2*is+1].whichSubDetector()==EcalPreshower && cp[2*is+1].whichLayer()==1)
00705             {
00706               CaloSegment preshsegment(cp[2*is],cp[2*is+1],s,sX0,sL0,CaloSegment::PS,myCalorimeter);
00707               segments_.push_back(preshsegment);
00708               //              std::cout << " Added (1-1)" << preshsegment << std::endl;
00709               s+=preshsegment.length();
00710               sX0+=preshsegment.X0length();
00711               sL0+=preshsegment.L0length();
00712               X0PS1_+=preshsegment.X0length();
00713               L0PS1_+=preshsegment.L0length();
00714             }
00715           else
00716             {
00717               std::cout << " Strange segment between Preshower1 and " << cp[2*is+1].whichDetector();
00718               std::cout << std::endl;
00719             }
00720           ++is;
00721           continue;
00722         }
00723 
00724       // Check if it is a Preshower segment - Layer 2 
00725       // One segment per layer, nothing between
00726       if(cp[2*is].whichDetector()==DetId::Ecal && cp[2*is].whichSubDetector()==EcalPreshower && cp[2*is].whichLayer()==2)
00727         {
00728           if(cp[2*is+1].whichDetector()==DetId::Ecal && cp[2*is+1].whichSubDetector()==EcalPreshower && cp[2*is+1].whichLayer()==2)
00729             {
00730               CaloSegment preshsegment(cp[2*is],cp[2*is+1],s,sX0,sL0,CaloSegment::PS,myCalorimeter);
00731               segments_.push_back(preshsegment);
00732               //              std::cout << " Added (1-2)" << preshsegment << std::endl;
00733               s+=preshsegment.length();
00734               sX0+=preshsegment.X0length();
00735               sL0+=preshsegment.L0length();
00736               X0PS2_+=preshsegment.X0length();
00737               L0PS2_+=preshsegment.L0length();
00738 
00739               // material between preshower and EE
00740               if(is<(nsegments-1) && cp[2*is+2].whichDetector()==DetId::Ecal && cp[2*is+2].whichSubDetector()==EcalEndcap)
00741                                 {
00742                   CaloSegment gapsef(cp[2*is+1],cp[2*is+2],s,sX0,sL0,CaloSegment::PSEEGAP,myCalorimeter);
00743                   segments_.push_back(gapsef);
00744                   s+=gapsef.length();     
00745                   sX0+=gapsef.X0length();                                                                                          
00746                   sL0+=gapsef.L0length();                                                                                          
00747                   X0PS2EE_+=gapsef.X0length();                        
00748                   L0PS2EE_+=gapsef.L0length();      
00749                   //              std::cout << " Created  a segment " << gapsef.length()<< " " << gapsef.X0length()<< std::endl;
00750                 }             
00751             }
00752           else
00753             {
00754               std::cout << " Strange segment between Preshower2 and " << cp[2*is+1].whichDetector();
00755               std::cout << std::endl;
00756             }
00757           ++is;
00758           continue;
00759         }
00760       // Now deal with the ECAL
00761       // One segment in each crystal. Segment corresponding to cracks/gaps are added
00762       //      myHistos->debug("Just avant ECAL"); 
00763       if(cp[2*is].whichDetector()==DetId::Ecal && (cp[2*is].whichSubDetector()==EcalBarrel || cp[2*is].whichSubDetector()==EcalEndcap))
00764         {
00765           if(cp[2*is+1].whichDetector()==DetId::Ecal && (cp[2*is+1].whichSubDetector()==EcalBarrel || cp[2*is+1].whichSubDetector()==EcalEndcap) )
00766             {
00767               DetId cell2=cp[2*is+1].getDetId();
00768               // set the real entrance
00769               if (ecalFirstSegment_<0) ecalFirstSegment_=segments_.size();
00770 
00771               // !! Approximatiom : the first segment is always in a crystal
00772               if(cp[2*is].getDetId()==cell2)
00773                 {
00774                   CaloSegment segment(cp[2*is],cp[2*is+1],s,sX0,sL0,CaloSegment::PbWO4,myCalorimeter);
00775                  segments_.push_back(segment);
00776                  // std::cout << " Added (2)" << segment << std::endl;
00777                  s+=segment.length();
00778                  sX0+=segment.X0length();
00779                  sL0+=segment.L0length(); 
00780                  X0ECAL_+=segment.X0length();
00781                  L0ECAL_+=segment.L0length();
00782                  ++ncrossedxtals;
00783                  ++is; 
00784                 }
00785               else
00786                 {
00787                   std::cout << " One more bug in the segment " <<std::endl;
00788                   ++is;
00789                 }
00790               // Now check if a gap or crack should be added
00791               if(is>0 && is<nsegments)
00792                 {                 
00793                   DetId cell3=cp[2*is].getDetId();
00794                   if(cp[2*is].whichDetector()!=DetId::Hcal) 
00795                     {
00796                       // Crack inside the ECAL
00797                       bool bordercrossing=myCalorimeter->borderCrossing(cell2,cell3);
00798                       CaloSegment cracksegment(cp[2*is-1],cp[2*is],s,sX0,sL0,(bordercrossing)?CaloSegment::CRACK:CaloSegment::GAP,myCalorimeter);
00799                       segments_.push_back(cracksegment);
00800                       s+=cracksegment.length();
00801                       sX0+=cracksegment.X0length();
00802                       sL0+=cracksegment.L0length();
00803                       X0ECAL_+=cracksegment.X0length();
00804                       L0ECAL_+=cracksegment.L0length();
00805                       //   std::cout <<" Added(3) "<< cracksegment << std::endl;
00806                     }
00807                   else
00808                     {
00809                       // a segment corresponding to ECAL/HCAL transition should be
00810                       // added here
00811                       CaloSegment cracksegment(cp[2*is-1],cp[2*is],s,sX0,sL0,CaloSegment::ECALHCALGAP,myCalorimeter);
00812                       segments_.push_back(cracksegment);
00813                       s+=cracksegment.length();
00814                       sX0+=cracksegment.X0length();
00815                       sL0+=cracksegment.L0length();
00816                       X0EHGAP_+=cracksegment.X0length();
00817                       L0EHGAP_+=cracksegment.L0length();
00818                     }
00819                 }
00820               continue;
00821             }
00822           else
00823             {
00824               std::cout << " Strange segment between " << cp[2*is].whichDetector();
00825               std::cout << " and " << cp[2*is+1].whichDetector() << std::endl;
00826               ++is;
00827               continue;
00828             }     
00829         }
00830       //      myHistos->debug("Just avant HCAL"); 
00831       // HCAL
00832       if(cp[2*is].whichDetector()==DetId::Hcal&&cp[2*is+1].whichDetector()==DetId::Hcal)
00833         {
00834           CaloSegment segment(cp[2*is],cp[2*is+1],s,sX0,sL0,CaloSegment::HCAL,myCalorimeter);
00835           segments_.push_back(segment);
00836           s+=segment.length();
00837           sX0+=segment.X0length();
00838           sL0+=segment.L0length(); 
00839           X0HCAL_+=segment.X0length();
00840           L0HCAL_+=segment.L0length();
00841           //      std::cout <<" Added(4) "<< segment << std::endl;
00842           ++is; 
00843         }
00844     }
00845 //  std::cout << " PS1 " << X0PS1_ << " " << L0PS1_ << std::endl;
00846 //  std::cout << " PS2 " << X0PS2_ << " " << L0PS2_ << std::endl;
00847 //  std::cout << " ECAL " << X0ECAL_ << " " << L0ECAL_ << std::endl;
00848 //  std::cout << " HCAL " << X0HCAL_ << " " << L0HCAL_ << std::endl;
00849   
00850   totalX0_ = X0PS1_+X0PS2_+X0PS2EE_+X0ECAL_+X0EHGAP_+X0HCAL_;
00851   totalL0_ = L0PS1_+L0PS2_+L0PS2EE_+L0ECAL_+L0EHGAP_+L0HCAL_;
00852   //  myHistos->debug("Just avant le fill"); 
00853 
00854   #ifdef DEBUGCELLLINE
00855   myHistos->fill("h200",fabs(EcalEntrance_.eta()),X0ECAL_);
00856   myHistos->fill("h210",EcalEntrance_.phi(),X0ECAL_);
00857   if(X0ECAL_<20)
00858     myHistos->fill("h212",EcalEntrance_.phi(),X0ECAL_);
00859 //  if(X0ECAL_<1.) 
00860 //    {
00861 //      for(unsigned ii=0; ii<segments_.size() ; ++ii)
00862 //      {
00863 //        std::cout << segments_[ii] << std::endl;
00864 //      }
00865 //    }
00866   myHistos->fillByNumber("h30",ncrossedxtals,EcalEntrance_.eta(),X0ECAL_);
00867   
00868   double zvertex = myTrack_->vertex().position().z();
00869 
00870   myHistos->fill("h310",EcalEntrance_.eta(),X0ECAL_);
00871   if(X0ECAL_<22)   myHistos->fill("h410",EcalEntrance_.phi());
00872   myHistos->fill("h400",zvertex,X0ECAL_);
00873   #endif
00874   //  std::cout << " Finished the segments " << std::endl;
00875     
00876 }
00877 
00878 void 
00879 EcalHitMaker::buildGeometry()
00880 {
00881   configuredGeometry_ = false;
00882   ncrystals_ = CellsWindow_.size();
00883   // create the vector with of pads with the appropriate size
00884   padsatdepth_.resize(ncrystals_);
00885   
00886   // This is fully correct in the barrel. 
00887   ny_= phisize_;
00888   nx_=ncrystals_/ny_;
00889   std::vector<unsigned> empty;
00890   empty.resize(ny_,0);
00891   myCrystalNumberArray_.reserve((unsigned)nx_);
00892   for(unsigned inx=0;inx<(unsigned)nx_;++inx)
00893     {      
00894       myCrystalNumberArray_.push_back(empty);
00895     }
00896 
00897   hits_.resize(ncrystals_,0.);
00898   regionOfInterest_.clear();
00899   regionOfInterest_.resize(ncrystals_);
00900   validPads_.resize(ncrystals_);
00901   for(unsigned ic=0;ic<ncrystals_;++ic)
00902     {
00903       myCalorimeter->buildCrystal(CellsWindow_[ic],regionOfInterest_[ic]);
00904       regionOfInterest_[ic].setNumber(ic);
00905       DetIdMap_.insert(std::pair<DetId,unsigned>(CellsWindow_[ic],ic));     
00906     }
00907   
00908   // Computes the map of the neighbours
00909   myCrystalWindowMap_ = new CrystalWindowMap(myCalorimeter,regionOfInterest_);
00910 }
00911 
00912 
00913 
00914 // depth is in X0 , L0 (depending on EMSHOWER/HADSHOWER) or in CM if inCM
00915 bool 
00916 EcalHitMaker::getPads(double depth,bool inCm) 
00917 {
00918   //std::cout << " New depth " << depth << std::endl;
00919   // The first time, the relationship between crystals must be calculated
00920   // but only in the case of EM showers
00921 
00922   if(EMSHOWER && !configuredGeometry_) configureGeometry();
00923 
00924 
00925   radiusFactor_ = (EMSHOWER) ? moliereRadius*radiusCorrectionFactor_:interactionLength;
00926   detailedShowerTail_ = false;
00927   if(EMSHOWER)
00928     currentdepth_ = depth+X0depthoffset_;
00929   else
00930     currentdepth_ = depth;
00931   
00932 //  if(currentdepth_>maxX0_+ps1TotalX0()+ps2TotalX0()) 
00933 //    {
00934 //      currentdepth_=maxX0_+ps1TotalX0()+ps2TotalX0()-1.; // the -1 is for safety
00935 //      detailedShowerTail_=true;
00936 //    }
00937     
00938 //  std::cout << " FamosGrid::getQuads " << currentdepth_ << " " << maxX0_ << std::endl;
00939 
00940   ncrackpadsatdepth_=0;
00941   
00942   xmin_=ymin_=999;
00943   xmax_=ymax_=-999;
00944   double locxmin,locxmax,locymin,locymax;
00945 
00946   // Get the depth of the pivot
00947   std::vector<CaloSegment>::const_iterator segiterator;
00948   // First identify the correct segment
00949 
00950   if(inCm) // centimeter
00951     {
00952       segiterator = find_if(segments_.begin(),segments_.end(),CaloSegment::inSegment(currentdepth_));
00953     }
00954   else
00955     {
00956       // EM shower 
00957       if(EMSHOWER)
00958         segiterator = find_if(segments_.begin(),segments_.end(),CaloSegment::inX0Segment(currentdepth_));
00959       
00960       //Hadron shower 
00961       if(HADSHOWER)
00962         segiterator = find_if(segments_.begin(),segments_.end(),CaloSegment::inL0Segment(currentdepth_));
00963     }
00964   if(segiterator==segments_.end()) 
00965     {
00966       std::cout << " FamosGrid: Could not go at such depth " << depth << std::endl;
00967       std::cout << " EMSHOWER " << EMSHOWER << std::endl;
00968       std::cout << " Track " << *myTrack_ << std::endl;
00969       std::cout << " Segments " << segments_.size() << std::endl;
00970       for(unsigned ii=0; ii<segments_.size() ; ++ii)
00971         {
00972           std::cout << segments_[ii] << std::endl;
00973         }
00974       
00975       return false;
00976     }
00977   //  std::cout << *segiterator << std::endl;
00978   
00979   if(segiterator->whichDetector()!=DetId::Ecal)
00980     {
00981       std::cout << " In  " << segiterator->whichDetector() << std::endl;
00982       //      buildPreshower();
00983       return false;
00984     }
00985   
00986   //  std::cout << *segiterator << std::endl;
00987   // get the position of the origin
00988   
00989   XYZPoint origin;
00990   if(inCm)
00991     {
00992       origin=segiterator->positionAtDepthincm(currentdepth_);
00993     }
00994   else
00995     {
00996       if(EMSHOWER)
00997         origin=segiterator->positionAtDepthinX0(currentdepth_);
00998       if(HADSHOWER)
00999         origin=segiterator->positionAtDepthinL0(currentdepth_);
01000     }
01001   //  std::cout << " currentdepth_ " << currentdepth_ << " " << origin << std::endl;
01002   XYZVector newaxis=pivot_.getFirstEdge().Cross(normal_);
01003 
01004 //  std::cout  << "Normal " << normal_ << std::endl;
01005 //  std::cout << " New axis " << newaxis << std::endl;
01006 
01007 //  std::cout << " ncrystals  " << ncrystals << std::endl;
01008   plan_ = Plane3D((Vector)normal_,(Point)origin);
01009 
01010   unsigned nquads=0;
01011   double sign=(central_) ? -1.: 1.;
01012   Transform3DR trans((Point)origin,(Point)(origin+normal_),(Point)(origin+newaxis),
01013                      Point(0,0,0), Point(0.,0.,sign),      Point(0.,1.,0.));
01014   for(unsigned ic=0;ic<ncrystals_;++ic)
01015     {
01016       //      std::cout << " Building geometry for " << regionOfInterest_[ic].getCellID() << std::endl;
01017       XYZPoint a,b;
01018 
01019       //      std::cout << " Origin " << origin << std::endl;
01020 
01021       //      std::vector<XYZPoint> corners;
01022       //      corners.reserve(4);
01023       double dummyt;
01024       bool hasbeenpulled=false;
01025       bool behindback=false;
01026       for(unsigned il=0;il<4;++il)
01027         {
01028           // a is the il-th front corner of the crystal. b is the corresponding rear corner
01029           regionOfInterest_[ic].getLateralEdges(il,a,b);
01030           
01031           // pull the surface if necessary (only in the front of the crystals) 
01032           XYZPoint aprime=a;
01033           if(pulled(origin,normal_,a)) 
01034             {
01035               b=aprime;
01036               hasbeenpulled=true;
01037             }
01038           
01039           // compute the intersection. 
01040           // Check that the intersection is in the [a,b] segment  if HADSHOWER
01041           // if EMSHOWER the intersection is calculated as if the crystals were infinite
01042           XYZPoint xx=(EMSHOWER)?intersect(plan_,a,b,dummyt,false):intersect(plan_,a,b,dummyt,true);
01043           
01044           if(dummyt>1) behindback=true;
01045           //      std::cout << " Intersect " << il << " " << a << " " << b << " " << plan_ << " " << xx << std::endl;
01046           // check that the intersection actually exists 
01047           if(xx.mag2()!=0)
01048             {
01049               corners[il] = xx;
01050             }     
01051         }    
01052       //      std::cout << " ncorners " << corners.size() << std::endl;
01053       if(behindback&&EMSHOWER) detailedShowerTail_=true;
01054       // If the quad is completly defined. Store it ! 
01055       if(corners.size()==4)
01056         {
01057           padsatdepth_[ic]=CrystalPad(ic,corners,trans,bfactor_);
01058           // Parameter to be tuned
01059           if(hasbeenpulled) padsatdepth_[ic].setSurvivalProbability(pulledPadProbability_);
01060           validPads_[ic]=true;
01061           ++nquads;     
01062           // In principle, this should be done after the quads reorganization. But it would cost one more loop
01063           //  quadsatdepth_[ic].extrems(locxmin,locxmax,locymin,locymax);
01064           //  if(locxmin<xmin_) xmin_=locxmin;
01065           //  if(locymin<ymin_) ymin_=locymin;
01066           //  if(locxmax>xmax_) xmax_=locxmax;
01067           //  if(locymax>ymax_) ymax_=locymax;
01068         }
01069       else
01070         {
01071           padsatdepth_[ic]=CrystalPad();
01072           validPads_[ic]=false;
01073         }
01074     }
01075   //  std::cout << " Number of quads " << quadsatdepth_.size() << std::endl; 
01076   if(doreorg_)reorganizePads();
01077   //  std::cout << "Finished to reorganize " << std::endl;
01078   npadsatdepth_=nquads;
01079   //  std::cout << " prepareCellIDMap " << std::endl;
01080   
01081 
01082   // Resize the Quads to allow for some numerical inaccuracy 
01083   // in the "inside" function
01084   for(unsigned ic=0;ic<ncrystals_;++ic) 
01085     {
01086     
01087       if (!validPads_[ic]) continue;
01088 
01089       if(EMSHOWER) padsatdepth_[ic].resetCorners();
01090 
01091       padsatdepth_[ic].extrems(locxmin,locxmax,locymin,locymax);
01092       if(locxmin<xmin_) xmin_=locxmin;
01093       if(locymin<ymin_) ymin_=locymin;
01094       if(locxmax>xmax_) xmax_=locxmax;
01095       if(locymax>ymax_) ymax_=locymax;
01096     }
01097   
01098   sizex_=(xmax_-xmin_)/nx_;
01099   sizey_=(ymax_-ymin_)/ny_;
01100 
01101 
01102   // Make sure that sizex_ and sizey_ are set before running prepareCellIDMap
01103   prepareCrystalNumberArray();
01104 
01105   //  std::cout << " Finished  prepareCellIDMap " << std::endl;
01106   ncrackpadsatdepth_=crackpadsatdepth_.size();
01107 
01108   return true;
01109 }
01110 
01111 
01112 
01113 
01114 void 
01115 EcalHitMaker::configureGeometry()
01116 { 
01117   configuredGeometry_=true;
01118   for(unsigned ic=0;ic<ncrystals_;++ic)
01119     {
01120       //      std::cout << " Building " << cellids_[ic] << std::endl;
01121       for(unsigned idir=0;idir<8;++idir)
01122         {
01123 
01124           unsigned oppdir=CaloDirectionOperations::oppositeDirection(idir);
01125           // Is there something else to do ? 
01126           // The relationship with the neighbour may have been set previously.
01127           if(regionOfInterest_[ic].crystalNeighbour(idir).status()>=0) 
01128             {
01129               //              std::cout << " Nothing to do " << std::endl;
01130               continue ;
01131             }
01132 
01133           const DetId & oldcell(regionOfInterest_[ic].getDetId());
01134           CaloDirection dir=CaloDirectionOperations::neighbourDirection(idir);
01135           DetId newcell(oldcell);
01136           if(!myCalorimeter->move(newcell,dir))
01137             {
01138               // no neighbour in this direction
01139               regionOfInterest_[ic].crystalNeighbour(idir).setStatus(-1);
01140               continue;
01141             }
01142           // Determine the number of this neighbour
01143           //      std::cout << " The neighbour is " << newcell << std::endl;
01144           std::map<DetId,unsigned>::const_iterator niter(DetIdMap_.find(newcell));
01145           if(niter==DetIdMap_.end())
01146             {
01147               //              std::cout << " The neighbour is not in the map " << std::endl;
01148               regionOfInterest_[ic].crystalNeighbour(idir).setStatus(-1);
01149               continue;
01150             }
01151           // Now there is a neighbour     
01152           //      std::cout << " The neighbour is " << niter->second << " " << cellids_[niter->second] << std::endl;
01153           regionOfInterest_[ic].crystalNeighbour(idir).setNumber(niter->second);
01154           //      std::cout << " Managed to set crystalNeighbour " << ic << " " << idir << std::endl;
01155           //      std::cout << " Trying  " << niter->second << " " << oppdir << std::endl;
01156           regionOfInterest_[niter->second].crystalNeighbour(oppdir).setNumber(ic);
01157           //      std::cout << " Crack/gap " << std::endl;
01158           if(myCalorimeter->borderCrossing(oldcell,newcell))
01159             {
01160               regionOfInterest_[ic].crystalNeighbour(idir).setStatus(1);
01161               regionOfInterest_[niter->second].crystalNeighbour(oppdir).setStatus(1);
01162               //              std::cout << " Crack ! " << std::endl;
01163             }
01164           else
01165             {
01166               regionOfInterest_[ic].crystalNeighbour(idir).setStatus(0);
01167               regionOfInterest_[niter->second].crystalNeighbour(oppdir).setStatus(0);
01168               //              std::cout << " Gap" << std::endl;
01169             }       
01170         }
01171     }
01172     // Magnetic field a la Charlot
01173   double theta=EcalEntrance_.theta();
01174   if(theta>M_PI_2) theta=M_PI-theta;
01175   bfactor_=1./(1.+0.133*theta);
01176   // the effect of the magnetic field in the EC is currently ignored 
01177   if(myCalorimeter->magneticField()==0. || !central_) bfactor_=1.;
01178 }
01179 
01180 // project fPoint on the plane (original,normal)
01181 bool 
01182 EcalHitMaker::pulled(const XYZPoint & origin, 
01183                      const XYZNormal& normal, 
01184                      XYZPoint & fPoint) const 
01185 {
01186   // check if fPoint is behind the origin
01187   double dotproduct=normal.Dot(fPoint-origin);
01188   if(dotproduct<=0.) return false;
01189 
01190   //norm of normal is 1 
01191   fPoint -= (1+dotproduct)*normal;
01192   return true;
01193 }
01194 
01195 
01196 void 
01197 EcalHitMaker::prepareCrystalNumberArray()
01198 {
01199   for(unsigned iq=0;iq<npadsatdepth_;++iq)
01200     {
01201       if(!validPads_[iq]) continue;
01202       unsigned d1,d2;
01203       convertIntegerCoordinates(padsatdepth_[iq].center().x(),padsatdepth_[iq].center().y(),d1,d2);
01204       myCrystalNumberArray_[d1][d2]=iq;
01205     }
01206 }
01207 
01208 void EcalHitMaker::convertIntegerCoordinates(double x, double y,unsigned &ix,unsigned &iy) const 
01209 {
01210   int tix=(int)((x-xmin_)/sizex_);
01211   int tiy=(int)((y-ymin_)/sizey_);
01212   ix=iy=9999;
01213   if(tix>=0) ix=(unsigned)tix;
01214   if(tiy>=0) iy=(unsigned)tiy;
01215 }
01216 
01217 const std::map<uint32_t,float>& EcalHitMaker::getHits() 
01218 {
01219   if (hitmaphasbeencalculated_) return hitMap_;
01220   for(unsigned ic=0;ic<ncrystals_;++ic)
01221     {
01222       hitMap_.insert(std::pair<uint32_t,double>(regionOfInterest_[ic].getDetId().rawId(),hits_[ic]));
01223     }
01224   hitmaphasbeencalculated_=true;
01225   return hitMap_;
01226 }
01227 
01228 
01230 
01231 void 
01232 EcalHitMaker::reorganizePads()
01233 {
01234 
01235   // Some cleaning first 
01236   //  std::cout << " Starting reorganize " << std::endl;
01237   crackpadsatdepth_.clear();
01238   crackpadsatdepth_.reserve(etasize_*phisize_);
01239   ncrackpadsatdepth_=0;
01240   std::vector<neighbour> gaps;
01241   std::vector<std::vector<neighbour> > cracks;  
01242   //  cracks.clear();
01243   cracks.resize(ncrystals_);
01244 
01245 
01246   for(unsigned iq=0;iq<ncrystals_;++iq)
01247     {
01248       if(!validPads_[iq]) continue;
01249 
01250       gaps.clear();
01251       //   std::cout << padsatdepth_[iq] << std::endl;
01252       //check all the directions
01253       for(unsigned iside=0;iside<4;++iside)
01254         {
01255           //      std::cout << " To be glued " << iside << " " << regionOfInterest_[iq].crystalNeighbour(iside).toBeGlued() << std::endl;
01256           CaloDirection thisside=CaloDirectionOperations::Side(iside);
01257           if(regionOfInterest_[iq].crystalNeighbour(iside).toBeGlued())
01258             {
01259               // look for the neighbour and check that it exists
01260               int neighbourstatus=regionOfInterest_[iq].crystalNeighbour(iside).status();
01261               if(neighbourstatus<0)
01262                 continue;
01263 
01264               unsigned neighbourNumber=regionOfInterest_[iq].crystalNeighbour(iside).number();
01265               if(!validPads_[neighbourNumber]) continue;
01266               // there is a crack between 
01267               if(neighbourstatus==1)
01268                 {
01269                   //              std::cout << " 1 Crack : " << thisside << " " << cellids_[iq]<< " " << cellids_[neighbourNumber] << std::endl;
01270                   cracks[iq].push_back(neighbour(thisside,neighbourNumber));
01271                 } // else it is a gap 
01272               else
01273                 {
01274                   gaps.push_back(neighbour(thisside,neighbourNumber));
01275                 }
01276             }
01277         }
01278       // Now lift the gaps
01279       gapsLifting(gaps,iq);
01280     }
01281 
01282   unsigned ncracks=cracks.size();
01283   //  std::cout << " Cracks treatment : " << cracks.size() << std::endl;
01284   for(unsigned icrack=0;icrack<ncracks;++icrack)
01285     {
01286       //      std::cout << " Crack number " << crackiter->first << std::endl;
01287       cracksPads(cracks[icrack],icrack);
01288     }
01289 
01290 }
01291 
01292 //dir 2 = N,E,W,S
01293 CLHEP::Hep2Vector & 
01294 EcalHitMaker::correspondingEdge(neighbour& myneighbour,CaloDirection dir2 ) 
01295 {
01296   CaloDirection dir=CaloDirectionOperations::oppositeSide(myneighbour.first);
01297   CaloDirection corner=CaloDirectionOperations::add2d(dir,dir2);
01298   //  std::cout << "Corresponding Edge " << dir<< " " << dir2 << " " << corner << std::endl;
01299   return padsatdepth_[myneighbour.second].edge(corner);  
01300 }
01301 
01302 bool EcalHitMaker::diagonalEdge(unsigned myPad, CaloDirection dir,CLHEP::Hep2Vector & point)
01303 {
01304   unsigned idir=CaloDirectionOperations::neighbourDirection(dir);
01305   if(regionOfInterest_[myPad].crystalNeighbour(idir).status()<0)
01306     return false;
01307   unsigned nneighbour=regionOfInterest_[myPad].crystalNeighbour(idir).number();
01308   if(!validPads_[nneighbour])
01309     {
01310       //      std::cout << " Wasn't able to move " << std::endl;
01311       return false;
01312     }
01313   point = padsatdepth_[nneighbour].edge(CaloDirectionOperations::oppositeSide(dir));
01314   return true;
01315 }
01316 
01317 bool EcalHitMaker::unbalancedDirection(const std::vector<neighbour>& dirs,unsigned & unb,unsigned & dir1, unsigned & dir2)
01318 {
01319   if(dirs.size()==1) return false;
01320   if(dirs.size()%2==0) return false;
01321   CaloDirection tmp;
01322   tmp=CaloDirectionOperations::add2d(dirs[0].first,dirs[1].first);
01323   if(tmp==NONE) 
01324     {
01325       unb=2;
01326       dir1=0;
01327       dir2=1;
01328       return true;
01329     }
01330   tmp=CaloDirectionOperations::add2d(dirs[0].first,dirs[2].first);
01331   if(tmp==NONE)
01332     { 
01333       unb=1;
01334       dir1=0;
01335       dir2=2;
01336       return true;
01337     }
01338   unb=0;
01339   dir1=1;
01340   dir2=2;
01341   return true;
01342 }
01343 
01344 
01345 void 
01346 EcalHitMaker::gapsLifting(std::vector<neighbour>& gaps,unsigned iq)
01347 {
01348   //  std::cout << " Entering gapsLifting "  << std::endl;
01349   CrystalPad & myPad = padsatdepth_[iq];
01350   unsigned ngaps=gaps.size();
01351   static bool debug=false;
01352   if(ngaps==1)
01353     {
01354       if(debug)
01355         {
01356           std::cout << " Avant " << ngaps << " " <<gaps[0].first<< std::endl;
01357           std::cout << myPad << std::endl;
01358         }
01359       if(gaps[0].first==NORTH||gaps[0].first==SOUTH)
01360         {
01361           CaloDirection dir1=CaloDirectionOperations::add2d(gaps[0].first,EAST);
01362           CaloDirection dir2=CaloDirectionOperations::add2d(gaps[0].first,WEST);
01363           myPad.edge(dir1)=correspondingEdge(gaps[0],EAST);
01364           myPad.edge(dir2)=correspondingEdge(gaps[0],WEST);
01365         }
01366       else
01367         {
01368           CaloDirection dir1=CaloDirectionOperations::add2d(gaps[0].first,NORTH);
01369           CaloDirection dir2=CaloDirectionOperations::add2d(gaps[0].first,SOUTH);
01370           myPad.edge(dir1)=correspondingEdge(gaps[0],NORTH);
01371           myPad.edge(dir2)=correspondingEdge(gaps[0],SOUTH);
01372         }
01373       if(debug)
01374         {
01375           std::cout << " Apres " << std::endl;
01376           std::cout << myPad << std::endl;
01377         }
01378     }
01379   else
01380     if(ngaps==2)
01381       {
01382         if(debug)
01383           {
01384             std::cout << " Avant " << ngaps << " " <<gaps[0].first<< " " <<gaps[1].first << std::endl;
01385             std::cout << myPad << std::endl;
01386             std::cout << " Voisin 1 " << (gaps[0].second) << std::endl;
01387             std::cout << " Voisin 2 " << (gaps[1].second) << std::endl;
01388           }
01389         CaloDirection corner0=CaloDirectionOperations::add2d(gaps[0].first,gaps[1].first);
01390 
01391         CLHEP::Hep2Vector point(0.,0.);
01392         if(corner0!=NONE&&diagonalEdge(iq,corner0,point))
01393           {
01394             CaloDirection corner1=CaloDirectionOperations::add2d(CaloDirectionOperations::oppositeSide(gaps[0].first),gaps[1].first);
01395             CaloDirection corner2=CaloDirectionOperations::add2d(gaps[0].first,CaloDirectionOperations::oppositeSide(gaps[1].first));
01396             myPad.edge(corner0) = point;
01397             myPad.edge(corner1) = correspondingEdge(gaps[1],CaloDirectionOperations::oppositeSide(gaps[0].first));
01398             myPad.edge(corner2) = correspondingEdge(gaps[0],CaloDirectionOperations::oppositeSide(gaps[1].first));
01399           }
01400         else
01401           if(corner0==NONE)
01402             {
01403               if(gaps[0].first==EAST||gaps[0].first==WEST)
01404                 {
01405                   CaloDirection corner1=CaloDirectionOperations::add2d(gaps[0].first,NORTH);
01406                   CaloDirection corner2=CaloDirectionOperations::add2d(gaps[0].first,SOUTH);
01407                   myPad.edge(corner1)=correspondingEdge(gaps[0],NORTH);
01408                   myPad.edge(corner2)=correspondingEdge(gaps[0],SOUTH);
01409                   
01410                   corner1=CaloDirectionOperations::add2d(gaps[1].first,NORTH);
01411                   corner2=CaloDirectionOperations::add2d(gaps[1].first,SOUTH);
01412                   myPad.edge(corner1)=correspondingEdge(gaps[1],NORTH);
01413                   myPad.edge(corner2)=correspondingEdge(gaps[1],SOUTH);
01414                 }
01415               else
01416                 {
01417                   CaloDirection corner1=CaloDirectionOperations::add2d(gaps[0].first,EAST);
01418                   CaloDirection corner2=CaloDirectionOperations::add2d(gaps[0].first,WEST);
01419                   myPad.edge(corner1)=correspondingEdge(gaps[0],EAST);
01420                   myPad.edge(corner2)=correspondingEdge(gaps[0],WEST);
01421                   
01422                   corner1=CaloDirectionOperations::add2d(gaps[1].first,EAST);
01423                   corner2=CaloDirectionOperations::add2d(gaps[1].first,WEST);
01424                   myPad.edge(corner1)=correspondingEdge(gaps[1],EAST);
01425                   myPad.edge(corner2)=correspondingEdge(gaps[1],WEST);
01426                 }
01427             }
01428         if(debug)
01429           {
01430             std::cout << " Apres " << std::endl;
01431             std::cout << myPad << std::endl;
01432           }
01433       }
01434     else
01435       if(ngaps==3)
01436         {
01437           // in this case the four corners have to be changed 
01438           unsigned iubd,idir1,idir2;
01439           CaloDirection diag;
01440           CLHEP::Hep2Vector point(0.,0.);
01441           //      std::cout << " Yes : 3 gaps" << std::endl;
01442           if(unbalancedDirection(gaps,iubd,idir1,idir2))                
01443             {
01444               CaloDirection ubd(gaps[iubd].first),dir1(gaps[idir1].first);
01445               CaloDirection dir2(gaps[idir2].first);
01446               
01447               //              std::cout << " Avant " << std::endl << myPad << std::endl;
01448               //              std::cout << ubd << " " << dir1 << " " << dir2 << std::endl;
01449               diag=CaloDirectionOperations::add2d(ubd,dir1);
01450               if(diagonalEdge(iq,diag,point))
01451                 myPad.edge(diag)=point;
01452               diag=CaloDirectionOperations::add2d(ubd,dir2);
01453               if(diagonalEdge(iq,diag,point))
01454                 myPad.edge(diag)=point;
01455               CaloDirection oppside=CaloDirectionOperations::oppositeSide(ubd);
01456               myPad.edge(CaloDirectionOperations::add2d(oppside,dir1))=correspondingEdge(gaps[idir1],oppside);
01457               myPad.edge(CaloDirectionOperations::add2d(oppside,dir2))=correspondingEdge(gaps[idir2],oppside);
01458               //              std::cout << " Apres " << std::endl << myPad << std::endl;
01459             }                         
01460         }
01461       else
01462         if(ngaps==4)
01463           {
01464             //      std::cout << " Waouh :4 gaps" << std::endl;
01465             //      std::cout << " Avant " << std::endl;
01466             //      std::cout << myPad<< std::endl;
01467             CLHEP::Hep2Vector point(0.,0.);
01468             if(diagonalEdge(iq,NORTHEAST,point)) 
01469               myPad.edge(NORTHEAST)=point;
01470             if(diagonalEdge(iq,NORTHWEST,point)) 
01471               myPad.edge(NORTHWEST)=point;
01472             if(diagonalEdge(iq,SOUTHWEST,point)) 
01473               myPad.edge(SOUTHWEST)=point;
01474             if(diagonalEdge(iq,SOUTHEAST,point)) 
01475               myPad.edge(SOUTHEAST)=point;      
01476             //      std::cout << " Apres " << std::endl;
01477             //      std::cout << myPad<< std::endl;
01478           }   
01479 }
01480 
01481 void 
01482 EcalHitMaker::cracksPads(std::vector<neighbour> & cracks, unsigned iq)
01483 {
01484   //  std::cout << " myPad " << &myPad << std::endl;
01485   unsigned ncracks=cracks.size();
01486   CrystalPad & myPad = padsatdepth_[iq];
01487   for(unsigned ic=0;ic<ncracks;++ic)
01488     {
01489       //      std::vector<CLHEP::Hep2Vector> mycorners;
01490       //      mycorners.reserve(4);
01491       switch(cracks[ic].first)
01492         {
01493         case NORTH:
01494           {
01495             mycorners[0] = (padsatdepth_[cracks[ic].second].edge(SOUTHWEST));
01496             mycorners[1] = (padsatdepth_[cracks[ic].second].edge(SOUTHEAST));
01497             mycorners[2] = (myPad.edge(NORTHEAST));
01498             mycorners[3] = (myPad.edge(NORTHWEST));
01499           }
01500           break;
01501         case SOUTH:
01502           {
01503             mycorners[0] = (myPad.edge(SOUTHWEST));
01504             mycorners[1] = (myPad.edge(SOUTHEAST));
01505             mycorners[2] = (padsatdepth_[cracks[ic].second].edge(NORTHEAST));
01506             mycorners[3] = (padsatdepth_[cracks[ic].second].edge(NORTHWEST));
01507           }
01508           break;
01509         case EAST:
01510           {
01511             mycorners[0] = (myPad.edge(NORTHEAST));
01512             mycorners[1] = (padsatdepth_[cracks[ic].second].edge(NORTHWEST));
01513             mycorners[2] = (padsatdepth_[cracks[ic].second].edge(SOUTHWEST));
01514             mycorners[3] = (myPad.edge(SOUTHEAST));
01515           }
01516           break;
01517         case WEST:
01518           {
01519             mycorners[0] = (padsatdepth_[cracks[ic].second].edge(NORTHEAST));
01520             mycorners[1] = (myPad.edge(NORTHWEST));
01521             mycorners[2] = (myPad.edge(SOUTHWEST));
01522             mycorners[3] = (padsatdepth_[cracks[ic].second].edge(SOUTHEAST));
01523           }
01524           break;
01525         default:
01526           {
01527           }
01528         }  
01529       CrystalPad crackpad(ic,mycorners);
01530       // to be tuned. A simpleconfigurable should be used
01531       crackpad.setSurvivalProbability(crackPadProbability_);   
01532       crackpadsatdepth_.push_back(crackpad);
01533     }
01534   //  std::cout << " Finished cracksPads " << std::endl;
01535 }
01536 
01537 
01538 bool EcalHitMaker::inside3D(const std::vector<XYZPoint>& corners, 
01539                             const XYZPoint& p) const
01540 {
01541   // corners and p are in the same plane
01542   // p is inside "corners" if the four crossproducts (corners[i]xcorners[i+1]) 
01543   // are in the same direction
01544 
01545   XYZVector crossproduct(0.,0.,0.),previouscrossproduct(0.,0.,0.);
01546 
01547   for(unsigned ip=0;ip<4 ; ++ip) {
01548     crossproduct = (corners[ip]-p).Cross(corners[(ip+1)%4]-p);
01549     if(ip==0)
01550         previouscrossproduct=crossproduct;
01551     else
01552       if (crossproduct.Dot(previouscrossproduct)<0.) return false;       
01553   }
01554 
01555   return true;
01556 }