CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/SimG4CMS/Calo/src/HFShowerParam.cc

Go to the documentation of this file.
00001 
00002 // File: HFShowerParam.cc
00003 // Description: Parametrized version of HF hits
00005 
00006 #include "SimG4CMS/Calo/interface/HFShowerParam.h"
00007 #include "SimG4CMS/Calo/interface/HFFibreFiducial.h"
00008 #include "DetectorDescription/Core/interface/DDFilter.h"
00009 #include "DetectorDescription/Core/interface/DDFilteredView.h"
00010 
00011 #include "FWCore/Utilities/interface/Exception.h"
00012 #include "FWCore/ServiceRegistry/interface/Service.h"
00013 #include "CommonTools/UtilAlgos/interface/TFileService.h"
00014 
00015 #include "G4VPhysicalVolume.hh"
00016 #include "G4Step.hh"
00017 #include "G4Track.hh"
00018 #include "G4NavigationHistory.hh"
00019 #include "Randomize.hh"
00020 
00021 #include "CLHEP/Units/GlobalPhysicalConstants.h"
00022 #include "CLHEP/Units/GlobalSystemOfUnits.h"
00023 
00024 //#define DebugLog
00025 
00026 HFShowerParam::HFShowerParam(std::string & name, const DDCompactView & cpv,
00027                              edm::ParameterSet const & p) : showerLibrary(0),
00028                                                             fibre(0), 
00029                                                             gflash(0),
00030                                                             fillHisto(false) {
00031 
00032   edm::ParameterSet m_HF  = p.getParameter<edm::ParameterSet>("HFShower");
00033   pePerGeV                = m_HF.getParameter<double>("PEPerGeV");
00034   trackEM                 = m_HF.getParameter<bool>("TrackEM");
00035   bool useShowerLibrary   = m_HF.getParameter<bool>("UseShowerLibrary");
00036   bool useGflash          = m_HF.getParameter<bool>("UseHFGflash");
00037   edMin                   = m_HF.getParameter<double>("EminLibrary");
00038   onlyLong                = m_HF.getParameter<bool>("OnlyLong");
00039   ref_index               = m_HF.getParameter<double>("RefIndex");
00040   double lambdaMean       = m_HF.getParameter<double>("LambdaMean");
00041   applyFidCut             = m_HF.getParameter<bool>("ApplyFiducialCut");
00042   parametrizeLast         = m_HF.getUntrackedParameter<bool>("ParametrizeLast",false);
00043   edm::LogInfo("HFShower") << "HFShowerParam::Use of shower library is set to "
00044                            << useShowerLibrary << " Use of Gflash is set to "
00045                            << useGflash << " P.E. per GeV " << pePerGeV
00046                            << ", ref. index of fibre " << ref_index 
00047                            << ", Track EM Flag " << trackEM << ", edMin "
00048                            << edMin << " GeV, use of Short fibre info in"
00049                            << " shower library set to " << !(onlyLong)
00050                            << ", use of parametrization for last part set to "
00051                            << parametrizeLast << ", Mean lambda " << lambdaMean
00052                            << ", Application of Fiducial Cut " << applyFidCut;
00053 
00054 #ifdef DebugLog
00055   edm::Service<TFileService> tfile;
00056   if (tfile.isAvailable()) {
00057     fillHisto = true;
00058     LogDebug("HFShower") << "HFShowerParam::Save histos in directory "
00059                          << "ProfileFromParam";
00060     TFileDirectory showerDir = tfile->mkdir("ProfileFromParam");
00061     hzv             = showerDir.make<TH1F>("hzv","Longitudinal Profile;Radiation Length;Number of Spots",330,0.0,1650.0);
00062     em_2d_1         = showerDir.make<TH2F>("em_2d_1","Lateral Profile vs. Shower Depth;cm;Events",800,800.0,1600.0,100,50.0,150.0);
00063     em_long_1       = showerDir.make<TH1F>("em_long_1","Longitudinal Profile;Radiation Length;Number of Spots",800,800.0,1600.0);
00064     em_long_1_tuned = showerDir.make<TH1F>("em_long_1_tuned","Longitudinal Profile;Radiation Length;Number of Spots",800,800.0,1600.0);
00065     em_lateral_1    = showerDir.make<TH1F>("em_lateral_1","Lateral Profile;cm;Events",100,50.0,150.0);
00066     em_2d_2         = showerDir.make<TH2F>("em_2d_2","Lateral Profile vs. Shower Depth;cm;Events",800,800.0,1600.0,100,50.0,150.0);
00067     em_long_2       = showerDir.make<TH1F>("em_long_2","Longitudinal Profile;Radiation Length;Number of Spots",800,800.0,1600.0);
00068     em_lateral_2    = showerDir.make<TH1F>("em_lateral_2","Lateral Profile;cm;Events",100,50.0,150.0);
00069     em_long_gflash  = showerDir.make<TH1F>("em_long_gflash","Longitudinal Profile From GFlash;cm;Number of Spots",800,800.0,1600.0);
00070     em_long_sl      = showerDir.make<TH1F>("em_long_sl","Longitudinal Profile From Shower Library;cm;Number of Spots",800,800.0,1600.0);
00071   } else {
00072     fillHisto = false;
00073     LogDebug("HFShower") << "HFShowerParam::No file is available for "
00074                              << "saving histos so the flag is set to false";
00075   }
00076 #endif
00077   
00078   G4String attribute = "ReadOutName";
00079   G4String value     = name;
00080   DDSpecificsFilter filter;
00081   DDValue           ddv(attribute,value,0);
00082   filter.setCriteria(ddv,DDSpecificsFilter::equals);
00083   DDFilteredView fv(cpv);
00084   fv.addFilter(filter);
00085   bool dodet = fv.firstChild();
00086   if (dodet) {
00087     DDsvalues_type sv(fv.mergedSpecifics());
00088     //Special Geometry parameters
00089     gpar      = getDDDArray("gparHF",sv);
00090     edm::LogInfo("HFShower") << "HFShowerParam: " <<gpar.size() <<" gpar (cm)";
00091     for (unsigned int ig=0; ig<gpar.size(); ig++)
00092       edm::LogInfo("HFShower") << "HFShowerParam: gpar[" << ig << "] = "
00093                                << gpar[ig]/cm << " cm";
00094   } else {
00095     edm::LogError("HFShower") << "HFShowerParam: cannot get filtered "
00096                               << " view for " << attribute << " matching "
00097                               << name;
00098     throw cms::Exception("Unknown", "HFShowerParam")
00099       << "cannot match " << attribute << " to " << name <<"\n";
00100   }
00101   
00102   if (useShowerLibrary) showerLibrary = new HFShowerLibrary(name, cpv, p);
00103   if (useGflash)        gflash        = new HFGflash(p);
00104   fibre = new HFFibre(name, cpv, p);
00105   attLMeanInv = fibre->attLength(lambdaMean);
00106   edm::LogInfo("HFShower") << "att. length used for (lambda=" << lambdaMean
00107                            << ") = " << 1/(attLMeanInv*cm) << " cm";
00108 }
00109 
00110 HFShowerParam::~HFShowerParam() {
00111   if (fibre)         delete fibre;
00112   if (gflash)        delete gflash;
00113   if (showerLibrary) delete showerLibrary;
00114 }
00115 
00116 void HFShowerParam::initRun(G4ParticleTable * theParticleTable) {
00117 
00118   emPDG = theParticleTable->FindParticle("e-")->GetPDGEncoding();
00119   epPDG = theParticleTable->FindParticle("e+")->GetPDGEncoding();
00120   gammaPDG = theParticleTable->FindParticle("gamma")->GetPDGEncoding();
00121 #ifdef DebugLog
00122   edm::LogInfo("HFShower") << "HFShowerParam: Particle code for e- = " << emPDG
00123                            << " for e+ = " << epPDG << " for gamma = " 
00124                            << gammaPDG;
00125 #endif
00126   if (showerLibrary) showerLibrary->initRun(theParticleTable);
00127 }
00128 
00129 std::vector<HFShowerParam::Hit> HFShowerParam::getHits(G4Step * aStep) {
00130 
00131   G4StepPoint * preStepPoint  = aStep->GetPreStepPoint(); 
00132   G4Track *     track    = aStep->GetTrack();   
00133   G4ThreeVector hitPoint = preStepPoint->GetPosition();   
00134   G4int         particleCode = track->GetDefinition()->GetPDGEncoding();
00135   double        zv = std::abs(hitPoint.z()) - gpar[4] - 0.5*gpar[1];
00136   G4ThreeVector localPoint = G4ThreeVector(hitPoint.x(),hitPoint.y(),zv);
00137 
00138   double pin    = (preStepPoint->GetTotalEnergy())/GeV;
00139   double zint   = hitPoint.z(); 
00140   double zz     = std::abs(zint) - gpar[4];
00141 
00142 #ifdef DebugLog
00143   edm::LogInfo("HFShower") << "HFShowerParam: getHits " 
00144                            << track->GetDefinition()->GetParticleName()
00145                            << " of energy " << pin << " GeV" 
00146                            << " Pos x,y,z = " << hitPoint.x() << "," 
00147                            << hitPoint.y() << "," << zint << " (" << zz << ","
00148                            << localPoint.z() << ", " 
00149                            << (localPoint.z()+0.5*gpar[1]) << ") Local " 
00150                            << localPoint;
00151 #endif
00152   std::vector<HFShowerParam::Hit> hits;
00153   HFShowerParam::Hit hit;
00154   hit.position = hitPoint;
00155 
00156   // look for other charged particles
00157   bool   other = false;
00158   double pBeta = track->GetDynamicParticle()->GetTotalMomentum() / track->GetDynamicParticle()->GetTotalEnergy();
00159   if (particleCode != emPDG && particleCode != epPDG && particleCode != gammaPDG ) {
00160     if (track->GetDefinition()->GetPDGCharge() != 0 && pBeta > (1/ref_index) &&
00161         aStep->GetTotalEnergyDeposit() > 0) other = true;
00162   }
00163 
00164   // take only e+-/gamma/or special particles
00165   if (particleCode == emPDG || particleCode == epPDG ||
00166       particleCode == gammaPDG || other) {
00167     // Leave out the last part
00168     double edep = 0.;
00169     bool   kill = false;
00170     if ((!trackEM) && ((zz<(gpar[1]-gpar[2])) || parametrizeLast) && (!other)){
00171       edep = pin;
00172       kill = true;
00173     } else if (track->GetDefinition()->GetPDGCharge() != 0 && 
00174                pBeta > (1/ref_index)) {
00175       edep = (aStep->GetTotalEnergyDeposit())/GeV;
00176     }
00177     std::string path = "ShowerLibrary";
00178 #ifdef DebugLog
00179     edm::LogInfo("HFShower") << "HFShowerParam: getHits edep = " << edep
00180                              << ", Kill = " << kill << ", pin = " << pin 
00181                              << ", edMin = " << edMin << " Other " << other;
00182 #endif
00183     if (edep > 0) {
00184       if ((showerLibrary || gflash) && kill && pin > edMin && (!other)) {
00185         if (showerLibrary) {
00186           std::vector<HFShowerLibrary::Hit> hitSL = showerLibrary->getHits(aStep,kill, onlyLong);
00187           for (unsigned int i=0; i<hitSL.size(); i++) {
00188             bool ok = true;
00189             if (applyFidCut) {
00190               int npmt = HFFibreFiducial:: PMTNumber(hitSL[i].position);
00191               if (npmt <= 0) ok = false;
00192             }
00193             if (ok) {
00194               hit.position = hitSL[i].position;
00195               hit.depth    = hitSL[i].depth;
00196               hit.time     = hitSL[i].time;
00197               hit.edep     = 1;
00198               hits.push_back(hit);
00199 #ifdef DebugLog
00200               if (fillHisto) {
00201                 hzv->Fill(zv);
00202                 em_long_sl->Fill(hit.position.z()/cm);
00203                 if (hit.depth == 1) {
00204                   em_2d_1->Fill(hit.position.z()/cm, sqrt(pow(hit.position.x()/cm,2)+pow(hit.position.y()/cm,2)));
00205                   em_lateral_1->Fill(sqrt(pow(hit.position.x()/cm,2)+pow(hit.position.y()/cm,2)));
00206                   em_long_1->Fill(hit.position.z()/cm);
00207                 } else if (hit.depth == 2) {
00208                   em_2d_2->Fill(hit.position.z()/cm, sqrt(pow(hit.position.x()/cm,2)+pow(hit.position.y()/cm,2)));
00209                   em_lateral_2->Fill(sqrt(pow(hit.position.x()/cm,2)+pow(hit.position.y()/cm,2)));
00210                   em_long_2->Fill(hit.position.z()/cm);
00211                 }
00212               }
00213 
00214               edm::LogInfo("HFShower") << "HFShowerParam: Hit at depth " 
00215                                        << hit.depth << " with edep " << hit.edep
00216                                        << " Time " << hit.time;
00217 #endif
00218             }
00219           }
00220         } else {
00221           std::vector<HFGflash::Hit>  hitSL = gflash->gfParameterization(aStep,kill, onlyLong);
00222           for (unsigned int i=0; i<hitSL.size(); i++) {
00223             bool ok = true;
00224             G4ThreeVector pe_effect(hitSL[i].position.x(), hitSL[i].position.y(), hitSL[i].position.z());
00225             double zv  = std::abs(pe_effect.z()) - gpar[4];
00226             //depth
00227             int depth    = 1;
00228             int npmt     = 0;
00229             if (zv < 0 )      ok = false;
00230             if (zv > gpar[1]) ok = false;
00231             if (applyFidCut) {
00232               npmt = HFFibreFiducial:: PMTNumber(pe_effect);
00233               if (npmt <= 0) ok = false;
00234               else if (npmt > 24) { // a short fibre
00235                 if (zv > gpar[0]) depth = 2; 
00236                 else              ok = false;
00237               }
00238 #ifdef DebugLog
00239               edm::LogInfo("HFShower") << "HFShowerParam: npmt " << npmt 
00240                                        << " zv " << std::abs(pe_effect.z()) 
00241                                        << ":" << gpar[4] << ":" << zv << ":" 
00242                                        << gpar[0] << " ok " << ok << " depth "
00243                                        << depth;
00244 #endif
00245             } else {
00246               if (G4UniformRand() > 0.5) depth = 2;
00247               if (depth == 2 && zv < gpar[0]) ok = false;
00248             }
00249             //attenuation
00250             double dist = fibre->zShift(localPoint,depth,0); // distance to PMT
00251             double r1   = G4UniformRand();
00252 #ifdef DebugLog
00253             edm::LogInfo("HFShower") << "HFShowerParam:Distance to PMT (" <<npmt
00254                                      << ") " << dist << ", exclusion flag " 
00255                                      << (r1 > exp(-attLMeanInv*zv));
00256 #endif
00257             if (r1 > exp(-attLMeanInv*dist)) ok = false;
00258             if (ok) {
00259               double time = fibre->tShift(localPoint,depth,0); //remaining part
00260 
00261               hit.position = hitSL[i].position;
00262               hit.depth    = depth;
00263               hit.time     = time + hitSL[i].time;
00264               hit.edep     = 1;
00265               hits.push_back(hit);
00266 #ifdef DebugLog
00267               if (fillHisto) {
00268                 em_long_gflash->Fill(pe_effect.z()/cm, hitSL[i].edep);
00269                 hzv->Fill(zv);
00270                 if (hit.depth == 1) {
00271                   em_2d_1->Fill(hit.position.z()/cm, sqrt(pow(hit.position.x()/cm,2)+pow(hit.position.y()/cm,2)));
00272                   em_lateral_1->Fill(sqrt(pow(hit.position.x()/cm,2)+pow(hit.position.y()/cm,2)));
00273                   em_long_1->Fill(hit.position.z()/cm);
00274                 } else if (hit.depth == 2) {
00275                   em_2d_2->Fill(hit.position.z()/cm, sqrt(pow(hit.position.x()/cm,2)+pow(hit.position.y()/cm,2)));
00276                   em_lateral_2->Fill(sqrt(pow(hit.position.x()/cm,2)+pow(hit.position.y()/cm,2)));
00277                   em_long_2->Fill(hit.position.z()/cm);
00278                 }
00279               }
00280 
00281               edm::LogInfo("HFShower") << "HFShowerParam: Hit at depth " 
00282                                        << hit.depth << " with edep " 
00283                                        << hit.edep << " Time "  << hit.time;
00284 #endif
00285             }
00286           }
00287         }
00288       } else {
00289         path          = "Rest";
00290         edep         *= pePerGeV;
00291         double tSlice = (aStep->GetPostStepPoint()->GetGlobalTime());
00292         double time = fibre->tShift(localPoint,1,0); // remaining part
00293         hit.depth   = 1;
00294         hit.time    = tSlice+time;
00295         hit.edep    = edep;
00296         hits.push_back(hit);
00297 #ifdef DebugLog
00298         edm::LogInfo("HFShower") << "HFShowerParam: Hit at depth 1 with edep " 
00299                                  << edep << " Time " << tSlice << ":" << time
00300                                  << ":" << hit.time;
00301 #endif
00302         if (zz >= gpar[0]) {
00303           time      = fibre->tShift(localPoint,2,0);
00304           hit.depth = 2;
00305           hit.time  = tSlice+time;
00306           hits.push_back(hit);
00307 #ifdef DebugLog
00308           edm::LogInfo("HFShower") <<"HFShowerParam: Hit at depth 2 with edep "
00309                                    << edep << " Time " << tSlice << ":" << time
00310                                    << hit.time;
00311 #endif
00312         }
00313       }
00314 #ifdef DebugLog
00315       for (unsigned int ii=0; ii<hits.size(); ii++) {
00316         double zv = std::abs(hits[ii].position.z());
00317         if (zv > 12790)
00318           edm::LogInfo("HFShower") << "HFShowerParam: Abnormal hit along " 
00319                                    << path << " in " 
00320                                    << preStepPoint->GetPhysicalVolume()->GetLogicalVolume()->GetName()
00321                                    << " at " << hits[ii].position << " zz " 
00322                                    << zv << " Edep " << edep << " due to " 
00323                                    << track->GetDefinition()->GetParticleName()
00324                                    << " time " << hit.time;
00325       }
00326 #endif
00327       if (kill) {
00328         track->SetTrackStatus(fStopAndKill);
00329         G4TrackVector tv = *(aStep->GetSecondary());
00330         for (unsigned int kk=0; kk<tv.size(); kk++) {
00331           if (tv[kk]->GetVolume() == preStepPoint->GetPhysicalVolume())
00332             tv[kk]->SetTrackStatus(fStopAndKill);
00333         }
00334       }
00335 #ifdef DebugLog
00336       edm::LogInfo("HFShower") << "HFShowerParam: getHits kill (" << kill
00337                                << ") track " << track->GetTrackID() 
00338                                << " at " << hitPoint
00339                                << " and deposit " << edep << " " << hits.size()
00340                                << " times" << " ZZ " << zz << " " << gpar[0];
00341 #endif
00342     }
00343   }
00344     
00345   return hits;
00346 }
00347 
00348 std::vector<double> HFShowerParam::getDDDArray(const std::string & str, 
00349                                                const DDsvalues_type & sv) {
00350 
00351 #ifdef DebugLog
00352   LogDebug("HFShower") << "HFShowerParam:getDDDArray called for " << str;
00353 #endif
00354   DDValue value(str);
00355   if (DDfetch(&sv,value)) {
00356 #ifdef DebugLog
00357     LogDebug("HFShower") << value;
00358 #endif
00359     const std::vector<double> & fvec = value.doubles();
00360     int nval = fvec.size();
00361     if (nval < 2) {
00362       edm::LogError("HFShower") << "HFShowerParam : # of " << str 
00363                                 << " bins " << nval << " < 2 ==> illegal";
00364       throw cms::Exception("Unknown", "HFShowerParam")
00365         << "nval < 2 for array " << str << "\n";
00366     }
00367 
00368     return fvec;
00369   } else {
00370     edm::LogError("HFShower") << "HFShowerParam : cannot get array " << str;
00371     throw cms::Exception("Unknown", "HFShowerParam") 
00372       << "cannot get array " << str << "\n";
00373   }
00374 }