CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_4_5_patch3/src/SimG4CMS/Calo/src/HCalSD.cc

Go to the documentation of this file.
00001 
00002 // File: HCalSD.cc
00003 // Description: Sensitive Detector class for calorimeters
00005 
00006 #include "SimG4CMS/Calo/interface/HCalSD.h"
00007 #include "SimG4CMS/Calo/interface/HcalTestNumberingScheme.h"
00008 #include "SimG4Core/Notification/interface/TrackInformation.h"
00009 #include "DataFormats/HcalDetId/interface/HcalDetId.h"
00010 #include "DetectorDescription/Core/interface/DDFilter.h"
00011 #include "DetectorDescription/Core/interface/DDFilteredView.h"
00012 #include "DetectorDescription/Core/interface/DDLogicalPart.h"
00013 #include "DetectorDescription/Core/interface/DDMaterial.h"
00014 #include "DetectorDescription/Core/interface/DDValue.h"
00015 #include "FWCore/Utilities/interface/Exception.h"
00016 
00017 #include "FWCore/ServiceRegistry/interface/Service.h"
00018 #include "CommonTools/UtilAlgos/interface/TFileService.h"
00019 
00020 #include "G4LogicalVolumeStore.hh"
00021 #include "G4LogicalVolume.hh"
00022 #include "G4Step.hh"
00023 #include "G4Track.hh"
00024 #include "G4ParticleTable.hh"
00025 
00026 #include <iostream>
00027 #include <fstream>
00028 #include <iomanip>
00029 
00030 //#define DebugLog
00031 //#define mkdebug
00032 //#define mkodebug
00033 
00034 HCalSD::HCalSD(G4String name, const DDCompactView & cpv,
00035                SensitiveDetectorCatalog & clg, 
00036                edm::ParameterSet const & p, const SimTrackManager* manager) : 
00037          CaloSD(name, cpv, clg, p, manager,
00038          p.getParameter<edm::ParameterSet>("HCalSD").getParameter<int>("TimeSliceUnit"),
00039          p.getParameter<edm::ParameterSet>("HCalSD").getParameter<bool>("IgnoreTrackID")), 
00040          numberingFromDDD(0), numberingScheme(0), showerLibrary(0), hfshower(0), 
00041          showerParam(0), showerPMT(0), showerBundle(0)
00042 {
00043 
00044   //static SimpleConfigurable<bool>   on1(false, "HCalSD:UseBirkLaw");
00045   //static SimpleConfigurable<double> bk1(0.013, "HCalSD:BirkC1");
00046   //static SimpleConfigurable<double> bk2(0.0568,"HCalSD:BirkC2");
00047   //static SimpleConfigurable<double> bk3(1.75,  "HCalSD:BirkC3");
00048   // Values from NIM 80 (1970) 239-244: as implemented in Geant3
00049   //static SimpleConfigurable<bool> on2(true,"HCalSD:UseShowerLibrary");
00050   edm::ParameterSet m_HF  = p.getParameter<edm::ParameterSet>("HFShower");
00051   applyFidCut      = m_HF.getParameter<bool>("ApplyFiducialCut");
00052   edm::ParameterSet m_HC = p.getParameter<edm::ParameterSet>("HCalSD");
00053   useBirk          = m_HC.getParameter<bool>("UseBirkLaw");
00054   birk1            = m_HC.getParameter<double>("BirkC1")*(g/(MeV*cm2));
00055   birk2            = m_HC.getParameter<double>("BirkC2");
00056   birk3            = m_HC.getParameter<double>("BirkC3");
00057   useShowerLibrary = m_HC.getParameter<bool>("UseShowerLibrary");
00058   useParam         = m_HC.getParameter<bool>("UseParametrize");
00059   bool testNumber  = m_HC.getParameter<bool>("TestNumberingScheme");
00060   usePMTHit        = m_HC.getParameter<bool>("UsePMTHits");
00061   betaThr          = m_HC.getParameter<double>("BetaThreshold");
00062   eminHitHB        = m_HC.getParameter<double>("EminHitHB")*MeV;
00063   eminHitHE        = m_HC.getParameter<double>("EminHitHE")*MeV;
00064   eminHitHO        = m_HC.getParameter<double>("EminHitHO")*MeV;
00065   eminHitHF        = m_HC.getParameter<double>("EminHitHF")*MeV;
00066   useFibreBundle   = m_HC.getParameter<bool>("UseFibreBundleHits");
00067   useHF            = m_HC.getUntrackedParameter<bool>("UseHF",true);
00068   bool forTBH2     = m_HC.getUntrackedParameter<bool>("ForTBH2",false);
00069   useLayerWt       = m_HC.getUntrackedParameter<bool>("UseLayerWt",false);
00070   std::string file = m_HC.getUntrackedParameter<std::string>("WtFile","None");
00071 #ifdef mkdebug
00072   std::cout<<"HCalSD::Con: * Called* applyFidCut="<<applyFidCut<<std::endl;
00073 #endif
00074 
00075 #ifdef DebugLog
00076   LogDebug("HcalSim") << "***************************************************" 
00077                       << "\n"
00078                       << "*                                                 *"
00079                       << "\n"
00080                       << "* Constructing a HCalSD  with name " << name << "\n"
00081                       << "*                                                 *"
00082                       << "\n"
00083                       << "***************************************************";
00084 #endif
00085   edm::LogInfo("HcalSim") << "HCalSD:: Use of HF code is set to " << useHF
00086                           << "\nUse of shower parametrization set to "
00087                           << useParam << "\nUse of shower library is set to " 
00088                           << useShowerLibrary << "\nUse PMT Hit is set to "
00089                           << usePMTHit << " with beta Threshold "<< betaThr
00090                           << "\nUSe of FibreBundle Hit set to "<<useFibreBundle
00091                           << "\n         Use of Birks law is set to      " 
00092                           << useBirk << "  with three constants kB = "
00093                           << birk1 << ", C1 = " << birk2 << ", C2 = " << birk3;
00094   edm::LogInfo("HcalSim") << "HCalSD:: Suppression Flag " << suppressHeavy
00095                           << " protons below " << kmaxProton << " MeV,"
00096                           << " neutrons below " << kmaxNeutron << " MeV and"
00097                           << " ions below " << kmaxIon << " MeV\n"
00098                           << "         Threshold for storing hits in HB: "
00099                           << eminHitHB << " HE: " << eminHitHE << " HO: "
00100                           << eminHitHO << " HF: " << eminHitHF;
00101 
00102   numberingFromDDD = new HcalNumberingFromDDD(name, cpv);
00103   HcalNumberingScheme* scheme;
00104   if (testNumber || forTBH2) 
00105        scheme = dynamic_cast<HcalNumberingScheme*>(new HcalTestNumberingScheme(forTBH2));
00106   else scheme = new HcalNumberingScheme();
00107   setNumberingScheme(scheme);
00108 
00109   const G4LogicalVolumeStore * lvs = G4LogicalVolumeStore::GetInstance();
00110   std::vector<G4LogicalVolume *>::const_iterator lvcite;
00111   G4LogicalVolume* lv;
00112   std::string attribute, value;
00113 #ifdef mkdebug
00114   std::cout<<"HCalSD::Con:HF="<<useHF<<",P="<<useParam<<",L="<<useShowerLibrary<<std::endl;
00115 #endif
00116   if (useHF)
00117   {
00118     if (useParam) showerParam = new HFShowerParam(name, cpv, p);
00119     else
00120     {
00121       if (useShowerLibrary) showerLibrary = new HFShowerLibrary(name, cpv, p);
00122       hfshower  = new HFShower(name, cpv, p, 0);
00123     }
00124 
00125     // HF volume names
00126     attribute = "Volume";
00127     value     = "HF";
00128     DDSpecificsFilter filter0;
00129     DDValue           ddv0(attribute, value, 0);
00130     filter0.setCriteria(ddv0, DDSpecificsFilter::equals);
00131     DDFilteredView fv0(cpv);
00132     fv0.addFilter(filter0);
00133     hfNames = getNames(fv0);
00134     fv0.firstChild();
00135     DDsvalues_type sv0(fv0.mergedSpecifics());
00136     std::vector<double> temp =  getDDDArray("Levels",sv0);
00137     edm::LogInfo("HcalSim") << "HCalSD: Names to be tested for " << attribute 
00138                             << " = " << value << " has " << hfNames.size() << " elements";
00139     for (unsigned int i=0; i < hfNames.size(); ++i)
00140     {
00141       G4String namv = hfNames[i];
00142       lv            = 0;
00143       for(lvcite=lvs->begin(); lvcite!=lvs->end(); lvcite++) if((*lvcite)->GetName()==namv)
00144       {
00145         lv = (*lvcite);
00146         break;
00147       }
00148       hfLV.push_back(lv);
00149       int level = static_cast<int>(temp[i]);
00150       hfLevels.push_back(level);
00151       edm::LogInfo("HcalSim") << "HCalSD:  HF[" << i << "] = " << hfNames[i]
00152                               << " LV " << hfLV[i] << " at level " << hfLevels[i];
00153     }
00154   
00155     // HF Fibre volume names
00156     value     = "HFFibre";
00157     DDSpecificsFilter filter1;
00158     DDValue           ddv1(attribute,value,0);
00159     filter1.setCriteria(ddv1, DDSpecificsFilter::equals);
00160     DDFilteredView fv1(cpv);
00161     fv1.addFilter(filter1);
00162     fibreNames = getNames(fv1);
00163 #ifdef mkdebug
00164     std::cout<<"HCalSD::Con: ===> fibNamesSize = "<<fibreNames.size()<<std::endl;
00165 #endif
00166     edm::LogInfo("HcalSim") << "HCalSD: Names to be tested for " << attribute 
00167                             << " = " << value << ":";
00168 #ifdef mkdebug
00169     std::cout<<"HCalSD::Con: fibNamesSize = "<<fibreNames.size()<<std::endl;
00170 #endif
00171     for (unsigned int i=0; i<fibreNames.size(); ++i)
00172     {
00173       G4String namv = fibreNames[i];
00174 #ifdef mkdebug
00175       std::cout<<"HCalSD::Con: FibreName["<<i<<"]="<<namv<<std::endl;
00176 #endif 
00177       lv            = 0;
00178       for (lvcite = lvs->begin(); lvcite != lvs->end(); ++lvcite) 
00179       {
00180 #ifdef mkdebug
00181         std::cout<<"HCalSD::Con: NameStored="<<(*lvcite)->GetName()<<std::endl;
00182 #endif
00183         if ((*lvcite)->GetName() == namv)
00184         {
00185           lv = (*lvcite);
00186           break;
00187         }
00188 #ifdef mkdebug
00189         std::cout<<"HCalSD::Con: New FibreLV_Name="<<lv->GetName()<<std::endl;
00190 #endif
00191       }
00192       fibreLV.push_back(lv);
00193       edm::LogInfo("HcalSim") << "HCalSD:  (" << i << ") " << fibreNames[i]
00194                               << " LV " << fibreLV[i];
00195     }
00196   
00197     // HF PMT volume names
00198     value     = "HFPMT";
00199     DDSpecificsFilter filter3;
00200     DDValue           ddv3(attribute,value,0);
00201     filter3.setCriteria(ddv3,DDSpecificsFilter::equals);
00202     DDFilteredView fv3(cpv);
00203     fv3.addFilter(filter3);
00204     std::vector<G4String> pmtNames = getNames(fv3);
00205     edm::LogInfo("HcalSim") << "HCalSD: Names to be tested for " << attribute << " = "
00206                             << value << " have " << pmtNames.size() << " entries";
00207     for (unsigned int i=0; i<pmtNames.size(); ++i)
00208     {
00209       G4String namv = pmtNames[i];
00210       lv            = 0;
00211       for (lvcite = lvs->begin(); lvcite != lvs->end(); ++lvcite) 
00212         if ((*lvcite)->GetName() == namv)
00213       {
00214         lv = (*lvcite);
00215         break;
00216       }
00217       pmtLV.push_back(lv);
00218       edm::LogInfo("HcalSim") << "HCalSD:  (" << i << ") " << pmtNames[i]
00219                               << " LV " << pmtLV[i];
00220     }
00221     if (pmtNames.size() > 0) showerPMT = new HFShowerPMT (name, cpv, p);
00222   
00223     // HF Fibre bundles
00224     value     = "HFFibreBundleStraight";
00225     DDSpecificsFilter filter4;
00226     DDValue           ddv4(attribute,value,0);
00227     filter4.setCriteria(ddv4,DDSpecificsFilter::equals);
00228     DDFilteredView fv4(cpv);
00229     fv4.addFilter(filter4);
00230     std::vector<G4String> fibreNames = getNames(fv4);
00231     edm::LogInfo("HcalSim") << "HCalSD: Names to be tested for " << attribute
00232                             << " = " << value << " have " << fibreNames.size()
00233                             << " entries";
00234     for (unsigned int i=0; i<fibreNames.size(); ++i)
00235     {
00236       G4String namv = fibreNames[i];
00237       lv            = 0;
00238       for (lvcite = lvs->begin(); lvcite != lvs->end(); lvcite++) 
00239         if ((*lvcite)->GetName() == namv)
00240       {
00241         lv = (*lvcite);
00242         break;
00243       }
00244       fibre1LV.push_back(lv);
00245       edm::LogInfo("HcalSim") << "HCalSD:  (" << i << ") " << fibreNames[i]
00246                               << " LV " << fibre1LV[i];
00247     }
00248 
00249     value     = "HFFibreBundleConical";
00250     DDSpecificsFilter filter5;
00251     DDValue           ddv5(attribute,value,0);
00252     filter5.setCriteria(ddv5,DDSpecificsFilter::equals);
00253     DDFilteredView fv5(cpv);
00254     fv5.addFilter(filter5);
00255     fibreNames = getNames(fv5);
00256     edm::LogInfo("HcalSim") << "HCalSD: Names to be tested for " << attribute << " = "
00257                             << value << " have " << fibreNames.size() << " entries";
00258     for (unsigned int i=0; i<fibreNames.size(); ++i)
00259     {
00260       G4String namv = fibreNames[i];
00261       lv            = 0;
00262       for (lvcite = lvs->begin(); lvcite != lvs->end(); ++lvcite) 
00263         if ((*lvcite)->GetName() == namv)
00264       {
00265         lv = (*lvcite);
00266         break;
00267       }
00268       fibre2LV.push_back(lv);
00269       edm::LogInfo("HcalSim") << "HCalSD:  (" << i << ") " << fibreNames[i]
00270                               << " LV " << fibre2LV[i];
00271     }
00272     if (fibre1LV.size() > 0 || fibre2LV.size() > 0) 
00273       showerBundle = new HFShowerFibreBundle (name, cpv, p);
00274   }
00275 
00276   //Material list for HB/HE/HO sensitive detectors
00277   attribute = "ReadOutName";
00278   DDSpecificsFilter filter2;
00279   DDValue           ddv2(attribute,name,0);
00280   filter2.setCriteria(ddv2,DDSpecificsFilter::equals);
00281   DDFilteredView fv2(cpv);
00282   fv2.addFilter(filter2);
00283   bool dodet = fv2.firstChild();
00284 
00285   DDsvalues_type sv(fv2.mergedSpecifics());
00286   //Layer0 Weight
00287   layer0wt = getDDDArray("Layer0Wt",sv);
00288   edm::LogInfo("HcalSim") << "HCalSD: " << layer0wt.size() << " Layer0Wt";
00289   for (unsigned int it=0; it<layer0wt.size(); ++it)
00290     edm::LogInfo("HcalSim") << "HCalSD: [" << it << "] " << layer0wt[it];
00291 
00292   const G4MaterialTable * matTab = G4Material::GetMaterialTable();
00293   std::vector<G4Material*>::const_iterator matite;
00294   while (dodet)
00295   {
00296     const DDLogicalPart & log = fv2.logicalPart();
00297     G4String namx = log.name().name();
00298     if (!isItHF(namx) && !isItFibre(namx))
00299     {
00300       bool notIn = true;
00301       for (unsigned int i=0; i<matNames.size(); ++i)
00302       {
00303         if (!strcmp(matNames[i].c_str(),log.material().name().name().c_str()))
00304         {
00305           notIn = false;
00306           break;
00307         }
00308       }
00309       if (notIn)
00310       {
00311         namx = log.material().name().name();
00312         matNames.push_back(namx);
00313         G4Material* mat;
00314         for (matite = matTab->begin(); matite != matTab->end(); ++matite)
00315         {
00316           if ((*matite)->GetName() == namx)
00317           {
00318             mat = (*matite);
00319             break;
00320           }
00321         }
00322         materials.push_back(mat);
00323       }
00324     }
00325     dodet = fv2.next();
00326   }
00327 
00328   edm::LogInfo("HcalSim") << "HCalSD: Material names for " << attribute 
00329                           << " = " << name << ":";
00330   for (unsigned int i=0; i<matNames.size(); ++i)
00331     edm::LogInfo("HcalSim") << "HCalSD: (" << i << ") " << matNames[i]
00332                             << " pointer " << materials[i];
00333 
00334   mumPDG = mupPDG = 0;
00335   
00336   if (useLayerWt) readWeightFromFile(file);
00337 
00338   for (int i=0;  i<9; ++i) hit_[i] = time_[i]= dist_[i] = 0;
00339 
00340 #ifdef DebugLog
00341   edm::Service<TFileService> tfile;
00342 
00343   if ( tfile.isAvailable() )
00344   {
00345     static std::string labels[9] = {"HB", "HE", "HO", "HF Absorber", "HF PMT",
00346                                     "HF Absorber Long", "HF Absorber Short",
00347                                     "HF PMT Long", "HF PMT Short"};
00348     char name[20], title[60];
00349     for (int i=0; i<9; ++i)
00350     {
00351       sprintf (title, "Hit energy in %s", labels[i].c_str());
00352       sprintf (name, "HCalSDHit%d", i);
00353       hit_[i] = tfile->make<TH1F>(name, title, 2000, 0., 2000.);
00354       sprintf (title, "Energy (MeV)");
00355       hit_[i]->GetXaxis()->SetTitle(title);
00356       hit_[i]->GetYaxis()->SetTitle("Hits");
00357       sprintf (title, "Time of the hit in %s", labels[i].c_str());
00358       sprintf (name, "HCalSDTime%d", i);
00359       time_[i] = tfile->make<TH1F>(name, title, 2000, 0., 2000.);
00360       sprintf (title, "Time (ns)");
00361       time_[i]->GetXaxis()->SetTitle(title);
00362       time_[i]->GetYaxis()->SetTitle("Hits");
00363       sprintf (title, "Longitudinal profile in %s", labels[i].c_str());
00364       sprintf (name, "HCalSDDist%d", i);
00365       dist_[i] = tfile->make<TH1F>(name, title, 2000, 0., 2000.);
00366       sprintf (title, "Distance (mm)");
00367       dist_[i]->GetXaxis()->SetTitle(title);
00368       dist_[i]->GetYaxis()->SetTitle("Hits");
00369     }
00370   }
00371 #endif
00372 }
00373 
00374 HCalSD::~HCalSD()
00375 { 
00376   if (numberingFromDDD) delete numberingFromDDD;
00377   if (numberingScheme)  delete numberingScheme;
00378   if (showerLibrary)    delete showerLibrary;
00379   if (hfshower)         delete hfshower;
00380   if (showerParam)      delete showerParam;
00381   if (showerPMT)        delete showerPMT;
00382   if (showerBundle)     delete showerBundle;
00383 }
00384 
00385 bool HCalSD::ProcessHits(G4Step * aStep, G4TouchableHistory * )
00386 { //  TimeMe t1( theHitTimer, false);
00387 #ifdef mkodebug
00388   static int nofH=0;
00389   int afH=0;
00390   int bfH=0;
00391 #endif 
00392   NaNTrap( aStep ) ;
00393   
00394   if (aStep == NULL) return true;
00395   else
00396   {
00397     G4LogicalVolume* lv= aStep->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume();
00398     G4String nameVolume = lv->GetName();
00399     if (isItHF(aStep))
00400     {
00401 #ifdef mkdebug
00402       if(nameVolume != "HVQX")
00403              std::cout<<"HCalSD::ProcessHits: P="<<useParam<<",LV="<<nameVolume<<std::endl;
00404 #endif 
00405       G4int parCode =aStep->GetTrack()->GetDefinition()->GetPDGEncoding();
00406       if (useParam)
00407       {
00408 #ifdef mkodebug
00409         bfH=getNumberOfHits();
00410 #endif 
00411 #ifdef DebugLog
00412         LogDebug("HcalSim") << "HCalSD: Hit from parametrization in " << nameVolume
00413                             << " for Track " << aStep->GetTrack()->GetTrackID() <<" ("
00414                             << aStep->GetTrack()->GetDefinition()->GetParticleName() <<")";
00415 #endif
00416 #ifdef mkdebug
00417         std::cout<<"HCalSD::ProcessHits: isItFibre = "<<isItFibre(lv)<<", applyFidCut="
00418                  <<applyFidCut<<std::endl;
00419 #endif 
00420         getFromParam(aStep);
00421 #ifdef mkodebug
00422         afH=getNumberOfHits();
00423         if(afH>bfH) std::cout<<"HCalSD::ProcessHits:*AfterParamS*, NHits="<<afH<<std::endl;
00424 #endif 
00425       }
00426       else
00427       {
00428         bool notaMuon = true;
00429         if (parCode == mupPDG || parCode == mumPDG ) notaMuon = false;
00430 #ifdef mkdebug
00431         if (isItFibre(lv)) std::cout<<"HCalSD::ProcessHits:L="<<useShowerLibrary
00432                                     << ", notMuon=" << notaMuon << std::endl;
00433 #endif 
00434         if (useShowerLibrary && notaMuon)
00435         {
00436 #ifdef mkodebug
00437           bfH=getNumberOfHits();
00438 #endif 
00439 #ifdef DebugLog
00440           LogDebug("HcalSim") << "HCalSD: Starts shower library from " 
00441                               << nameVolume << " for Track " 
00442                               << aStep->GetTrack()->GetTrackID() << " ("
00443                               <<aStep->GetTrack()->GetDefinition()->GetParticleName()<<")";
00444 #endif
00445           getFromLibrary(aStep);
00446 #ifdef mkodebug
00447           afH=getNumberOfHits();
00448           if(afH>bfH) std::cout<<"HCalSD::ProcessHits:*AfterLibr*, NHits="<<afH<<std::endl;
00449 #endif 
00450         }
00451         else if (isItFibre(lv) || applyFidCut)
00452         {
00453 #ifdef mkodebug
00454           bfH=getNumberOfHits();
00455 #endif 
00456 #ifdef DebugLog
00457           LogDebug("HcalSim") << "HCalSD: Hit at Fibre in " << nameVolume 
00458                               << " for Track " 
00459                               << aStep->GetTrack()->GetTrackID() <<" ("
00460                               <<aStep->GetTrack()->GetDefinition()->GetParticleName()<<")";
00461 #endif
00462           hitForFibre(aStep);
00463 #ifdef mkodebug
00464           afH=getNumberOfHits();
00465           if(afH>bfH) std::cout<<"HCalSD::ProcessHits:*AfterFibr*, NHits="<<afH<<std::endl;
00466 #endif 
00467         }
00468       }
00469     }
00470     else if (isItPMT(lv))
00471     {
00472 #ifdef mkodebug
00473       bfH=getNumberOfHits();
00474 #endif 
00475 #ifdef DebugLog
00476       LogDebug("HcalSim") << "HCalSD: Hit from PMT parametrization from " 
00477                           <<  nameVolume << " for Track " 
00478                           << aStep->GetTrack()->GetTrackID() << " ("
00479                           << aStep->GetTrack()->GetDefinition()->GetParticleName() << ")";
00480 #endif
00481       if (usePMTHit && showerPMT) getHitPMT(aStep);
00482 #ifdef mkodebug
00483       afH=getNumberOfHits();
00484       if(afH>bfH) std::cout<<"***>>***>>HCalSD::ProcessHits:*PMT*, NHits="<<afH<<std::endl;
00485 #endif 
00486     }
00487     else if (isItStraightBundle(lv) || isItConicalBundle(lv))
00488     {
00489 #ifdef mkodebug
00490       bfH=getNumberOfHits();
00491 #endif 
00492 #ifdef DebugLog
00493       LogDebug("HcalSim") << "HCalSD: Hit from FibreBundle from "
00494                           << nameVolume << " for Track " 
00495                           << aStep->GetTrack()->GetTrackID() << " ("
00496                           << aStep->GetTrack()->GetDefinition()->GetParticleName() << ")";
00497 #endif
00498       if (useFibreBundle && showerBundle) getHitFibreBundle(aStep, isItConicalBundle(lv));
00499 #ifdef mkodebug
00500       afH=getNumberOfHits();
00501       if(afH>bfH) std::cout<<"*>>*>>HCalSD::ProcessHits:*FBundle*, NHits="<<afH<<std::endl;
00502 #endif 
00503     }
00504     else
00505     {
00506 #ifdef mkodebug
00507       bfH=getNumberOfHits();
00508 #endif 
00509 #ifdef DebugLog
00510       LogDebug("HcalSim") << "HCalSD: Hit from standard path from " 
00511                           <<  nameVolume << " for Track " 
00512                           << aStep->GetTrack()->GetTrackID() << " ("
00513                           << aStep->GetTrack()->GetDefinition()->GetParticleName() << ")";
00514 #endif
00515       if (getStepInfo(aStep) && !applyFidCut)
00516       {
00517 #ifdef DebugLog
00518         if (edepositEM+edepositHAD > 0)
00519           plotProfile(aStep, aStep->GetPreStepPoint()->GetPosition(),
00520                       edepositEM+edepositHAD,aStep->GetPostStepPoint()->GetGlobalTime(),0);
00521 #endif
00522         if (hitExists() == false && edepositEM+edepositHAD>0.) currentHit = createNewHit();
00523       }
00524 #ifdef mkodebug
00525       afH=getNumberOfHits();
00526       if(afH>bfH) std::cout<< "**>>**>>HCalSD::ProcessHits:*Else*, NHits=" << afH
00527                            << ", eEM=" << edepositEM<< ", eHD=" << edepositHAD <<std::endl;
00528 #endif 
00529     }
00530 #ifdef mkodebug
00531     bfH=getNumberOfHits();
00532     if(bfH>nofH) std::cout<<"HCalSD::ProcessHits:*Done* increased#OfHits="<<bfH<<std::endl;
00533     nofH=bfH;
00534 #endif 
00535     return true;
00536   }
00537 } 
00538 
00539 double HCalSD::getEnergyDeposit(G4Step* aStep)
00540 {
00541   double destep = aStep->GetTotalEnergyDeposit();
00542   double weight = 1;
00543   const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
00544   int depth = (touch->GetReplicaNumber(0))%10;
00545   int det   = ((touch->GetReplicaNumber(1))/1000)-3;
00546   if (depth==0 && (det==0 || det==1)) weight = layer0wt[det];
00547   if (useLayerWt)
00548   {
00549     int lay   = (touch->GetReplicaNumber(0)/10)%100 + 1;
00550     G4ThreeVector hitPoint = aStep->GetPreStepPoint()->GetPosition();
00551     weight = layerWeight(det+3, hitPoint, depth, lay);
00552   }
00553 
00554   if (suppressHeavy)
00555   {
00556     G4Track* theTrack = aStep->GetTrack();
00557     TrackInformation * trkInfo = (TrackInformation *)(theTrack->GetUserInformation());
00558     if (trkInfo)
00559     {
00560       int pdg = theTrack->GetDefinition()->GetPDGEncoding();
00561       if (!(trkInfo->isPrimary()))          // Only secondary particles
00562       {
00563         double ke = theTrack->GetKineticEnergy()/MeV;
00564         if ( pdg/1000000000 == 1 && (pdg/10000)%100 > 0 &&
00565                  (pdg/10)%100 > 0    && ke <kmaxIon         ) weight = 0;
00566         if ((pdg == 2212) && (ke < kmaxProton))     weight = 0;
00567         if ((pdg == 2112) && (ke < kmaxNeutron))    weight = 0;
00568 #ifdef DebugLog
00569         if (weight == 0) 
00570           edm::LogInfo("HcalSim") << "HCalSD:Ignore Track " 
00571                                   << theTrack->GetTrackID() << " Type " 
00572                                   << theTrack->GetDefinition()->GetParticleName()
00573                                   << " Kinetic Energy " << ke << " MeV";
00574 #endif
00575       }
00576     }
00577   }
00578 #ifdef DebugLog
00579   double weight0 = weight;
00580 #endif
00581   if (useBirk)
00582   {
00583     G4Material* mat = aStep->GetPreStepPoint()->GetMaterial();
00584     if (isItScintillator(mat)) weight *= getAttenuation(aStep, birk1, birk2, birk3);
00585   }
00586   double wt1    = getResponseWt(theTrack);
00587 #ifdef DebugLog
00588   edm::LogInfo("HcalSim") << "HCalSD: Detector " << det+3 << " Depth " << depth
00589                           << " weight " << weight0 << " " << weight << " " << wt1; 
00590 #endif
00591   return weight*wt1*destep;
00592 }
00593 
00594 uint32_t HCalSD::setDetUnitId(G4Step * aStep)
00595 { 
00596   G4StepPoint* preStepPoint = aStep->GetPreStepPoint(); 
00597   const G4VTouchable* touch = preStepPoint->GetTouchable();
00598   G4ThreeVector hitPoint    = preStepPoint->GetPosition();
00599 
00600   int depth = (touch->GetReplicaNumber(0))%10 + 1;
00601   int lay   = (touch->GetReplicaNumber(0)/10)%100 + 1;
00602   int det   = (touch->GetReplicaNumber(1))/1000;
00603 
00604   return setDetUnitId (det, hitPoint, depth, lay);
00605 }
00606 
00607 void HCalSD::setNumberingScheme(HcalNumberingScheme * scheme)
00608 {
00609   if (scheme != 0)
00610   {
00611     edm::LogInfo("HcalSim") << "HCalSD: updates numbering scheme for " << GetName() <<"\n";
00612     if (numberingScheme) delete numberingScheme;
00613     numberingScheme = scheme;
00614   }
00615 }
00616 
00617 void HCalSD::initRun()
00618 {
00619   G4ParticleTable * theParticleTable = G4ParticleTable::GetParticleTable();
00620   G4String          particleName;
00621   mumPDG = theParticleTable->FindParticle(particleName="mu-")->GetPDGEncoding();
00622   mupPDG = theParticleTable->FindParticle(particleName="mu+")->GetPDGEncoding();
00623 #ifdef DebugLog
00624   edm::LogInfo("HcalSim") << "HCalSD: Particle code for mu- = " << mumPDG
00625                           << " for mu+ = " << mupPDG;
00626 #endif
00627   if (showerLibrary) showerLibrary->initRun(theParticleTable);
00628   if (showerParam)   showerParam->initRun(theParticleTable);
00629   if (hfshower)      hfshower->initRun(theParticleTable);
00630 }
00631 
00632 bool HCalSD::filterHit(CaloG4Hit* aHit, double time)
00633 {
00634   double threshold=0;
00635   DetId theId(aHit->getUnitID());
00636   switch (theId.subdetId())
00637   {
00638   case HcalBarrel:
00639     threshold = eminHitHB; break;
00640   case HcalEndcap:
00641     threshold = eminHitHE; break;
00642   case HcalOuter:
00643     threshold = eminHitHO; break;
00644   case HcalForward:
00645     threshold = eminHitHF; break;
00646   default:
00647     break;
00648   }
00649   return ((time <= tmaxHit) && (aHit->getEnergyDeposit() > threshold));
00650 }
00651 
00652 
00653 uint32_t HCalSD::setDetUnitId (int det, G4ThreeVector pos, int depth, int lay=1)
00654 { 
00655   uint32_t id = 0;
00656   if (numberingFromDDD)
00657   {
00658     //get the ID's as eta, phi, depth, ... indices
00659     HcalNumberingFromDDD::HcalID tmp = numberingFromDDD->unitID(det, pos, depth, lay);
00660     //get the ID
00661     if (numberingScheme) id = numberingScheme->getUnitID(tmp);
00662   }
00663   return id;
00664 }
00665 
00666 std::vector<double> HCalSD::getDDDArray(const std::string & str,
00667                                         const DDsvalues_type & sv)
00668 {
00669 #ifdef DebugLog
00670   LogDebug("HcalSim") << "HCalSD:getDDDArray called for " << str;
00671 #endif
00672   DDValue value(str);
00673   if (DDfetch(&sv,value))
00674   {
00675 #ifdef DebugLog
00676     LogDebug("HcalSim") << value;
00677 #endif
00678     const std::vector<double> & fvec = value.doubles();
00679     int nval = fvec.size();
00680     if (nval < 1)
00681     {
00682       edm::LogError("HcalSim") << "HCalSD : # of " << str << " bins " << nval
00683           << " < 2 ==> illegal";
00684       throw cms::Exception("Unknown", "HCalSD") << "nval < 2 for array " << str << "\n";
00685     }
00686     
00687     return fvec;
00688   }
00689   else
00690   {
00691     edm::LogError("HcalSim") << "HCalSD :  cannot get array " << str;
00692     throw cms::Exception("Unknown", "HCalSD") << "cannot get array " << str << "\n";
00693   }
00694 }
00695 
00696 std::vector<G4String> HCalSD::getNames(DDFilteredView& fv) // Get Fibre volumes names
00697 {
00698   std::vector<G4String> tmp;
00699   bool dodet = fv.firstChild();
00700 #ifdef mkdebug
00701   std::cout<<"HCalSD::getNames: dodet = "<<dodet<<std::endl;
00702 #endif
00703   while (dodet)
00704   {
00705     const DDLogicalPart & log = fv.logicalPart();
00706     bool ok = true;
00707 
00708     for (unsigned int i=0; i<tmp.size(); ++i)
00709     {
00710       if (!strcmp(tmp[i].c_str(), log.name().name().c_str()))
00711       {
00712         ok = false;
00713         break;
00714       }
00715     }
00716 #ifdef mkdebug
00717     if (ok)std::cout<<"HCalSD::getNames: Size= "<<tmp.size()<<", Name= "<<log.name().name()
00718                     <<std::endl;
00719 #endif
00720     if (ok) tmp.push_back(log.name().name());
00721     dodet = fv.next();
00722   }
00723   return tmp;
00724 }
00725 
00726 bool HCalSD::isItHF(G4Step * aStep)
00727 {
00728   const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
00729   int levels = (touch->GetHistoryDepth()) + 1;
00730   for (unsigned int it=0; it < hfNames.size(); ++it)
00731   {
00732     if (levels >= hfLevels[it])
00733     {
00734       G4LogicalVolume* lv = touch->GetVolume(levels-hfLevels[it])->GetLogicalVolume();
00735       if (lv == hfLV[it]) return true;
00736     }
00737   }
00738   return false;
00739 }
00740 
00741 bool HCalSD::isItHF (G4String name)
00742 {
00743   std::vector<G4String>::const_iterator it = hfNames.begin();
00744   for (; it != hfNames.end(); ++it) if (name == *it) return true;
00745   return false;
00746 }
00747 
00748 bool HCalSD::isItFibre (G4LogicalVolume* lv)
00749 {
00750   std::vector<G4LogicalVolume*>::const_iterator ite = fibreLV.begin();
00751   for (; ite != fibreLV.end(); ++ite) if (lv == *ite) return true;
00752   return false;
00753 }
00754 
00755 bool HCalSD::isItFibre (G4String name)
00756 {
00757   std::vector<G4String>::const_iterator it = fibreNames.begin();
00758   for (; it != fibreNames.end(); ++it) if (name == *it) return true;
00759   return false;
00760 }
00761 
00762 bool HCalSD::isItPMT (G4LogicalVolume* lv)
00763 {
00764   std::vector<G4LogicalVolume*>::const_iterator ite = pmtLV.begin();
00765   for (; ite != pmtLV.end(); ++ite) if (lv == *ite) return true;
00766   return false;
00767 }
00768 
00769 bool HCalSD::isItStraightBundle (G4LogicalVolume* lv)
00770 {
00771   std::vector<G4LogicalVolume*>::const_iterator ite = fibre1LV.begin();
00772   for (; ite != fibre1LV.end(); ++ite) if (lv == *ite) return true;
00773   return false;
00774 }
00775 
00776 bool HCalSD::isItConicalBundle (G4LogicalVolume* lv)
00777 {
00778   std::vector<G4LogicalVolume*>::const_iterator ite = fibre2LV.begin();
00779   for (; ite != fibre2LV.end(); ++ite) if (lv == *ite) return true;
00780   return false;
00781 }
00782 
00783 bool HCalSD::isItScintillator (G4Material* mat)
00784 {
00785   std::vector<G4Material*>::const_iterator ite = materials.begin();
00786   for (; ite != materials.end(); ++ite) if (mat == *ite) return true;
00787   return false;
00788 }
00789 
00790 void HCalSD::getFromLibrary (G4Step* aStep)
00791 {
00792   preStepPoint  = aStep->GetPreStepPoint(); 
00793   theTrack      = aStep->GetTrack();   
00794   int det       = 5;
00795   bool ok;
00796 
00797   std::vector<HFShowerLibrary::Hit> hits = showerLibrary->getHits(aStep, ok);
00798 
00799   double etrack    = preStepPoint->GetKineticEnergy();
00800   int    primaryID = setTrackID(aStep);
00801   /*
00802   int    primaryID = 0;
00803   if (etrack >= energyCut) {
00804     primaryID    = theTrack->GetTrackID();
00805   } else {
00806     primaryID    = theTrack->GetParentID();
00807     if (primaryID == 0) primaryID = theTrack->GetTrackID();
00808   }
00809   */
00810 
00811   // Reset entry point for new primary
00812   posGlobal = preStepPoint->GetPosition();
00813   resetForNewPrimary(posGlobal, etrack);
00814 
00815   G4int particleCode = theTrack->GetDefinition()->GetPDGEncoding();
00816   if (particleCode==emPDG || particleCode==epPDG || particleCode==gammaPDG)
00817   {
00818     edepositEM  = 1.*GeV;
00819     edepositHAD = 0.;
00820   }
00821   else
00822   {
00823     edepositEM  = 0.;
00824     edepositHAD = 1.*GeV;
00825   }
00826 #ifdef DebugLog
00827   edm::LogInfo("HcalSim") << "HCalSD::getFromLibrary " <<hits.size() 
00828                           << " hits for " << GetName() << " of " << primaryID 
00829                           << " with " << theTrack->GetDefinition()->GetParticleName() 
00830                           << " of " << preStepPoint->GetKineticEnergy()/GeV << " GeV";
00831 #endif
00832   for (unsigned int i=0; i<hits.size(); ++i)
00833   {
00834     G4ThreeVector hitPoint = hits[i].position;
00835     int depth              = hits[i].depth;
00836     double time            = hits[i].time;
00837     unsigned int unitID    = setDetUnitId(det, hitPoint, depth);
00838     currentID.setID(unitID, time, primaryID, 0);
00839 #ifdef DebugLog
00840     plotProfile(aStep, hitPoint, 1.0*GeV, time, depth);
00841 #endif
00842    
00843     // check if it is in the same unit and timeslice as the previous one
00844     if (currentID == previousID) updateHit(currentHit);
00845     else if (!checkHit()) currentHit = createNewHit();
00846   }
00847 
00848   //Now kill the current track
00849   if (ok)
00850   {
00851     theTrack->SetTrackStatus(fStopAndKill);
00852     G4TrackVector tv = *(aStep->GetSecondary());
00853     for (unsigned int kk=0; kk<tv.size(); ++kk)
00854       if (tv[kk]->GetVolume() == preStepPoint->GetPhysicalVolume())
00855         tv[kk]->SetTrackStatus(fStopAndKill);
00856   }
00857 }
00858 
00859 void HCalSD::hitForFibre (G4Step* aStep) // if not ParamShower
00860 {
00861   preStepPoint  = aStep->GetPreStepPoint();
00862   theTrack      = aStep->GetTrack();
00863   int primaryID = setTrackID(aStep);
00864 
00865   int det   = 5;
00866 #ifdef mkdebug
00867   std::cout<<"HCalSD::HitForFibre: *Before getHits*, det="<<det<<std::endl;
00868 #endif 
00869   std::vector<HFShower::Hit> hits = hfshower->getHits(aStep);
00870 
00871   G4int particleCode = theTrack->GetDefinition()->GetPDGEncoding();
00872   if (particleCode==emPDG || particleCode==epPDG || particleCode==gammaPDG)
00873   {
00874     edepositEM  = 1.*GeV;
00875     edepositHAD = 0.;
00876   }
00877   else
00878   {
00879     edepositEM  = 0.;
00880     edepositHAD = 1.*GeV;
00881   }
00882  
00883 #ifdef DebugLog
00884   edm::LogInfo("HcalSim") << "HCalSD::hitForFibre " << hits.size() 
00885      << " hits for " << GetName() << " of " << primaryID 
00886      << " with " << theTrack->GetDefinition()->GetParticleName() 
00887      << " of " << preStepPoint->GetKineticEnergy()/GeV 
00888      << " GeV in detector type " << det;
00889 #endif
00890   if (hits.size() > 0)
00891   {
00892     G4ThreeVector hitPoint = preStepPoint->GetPosition();
00893     int depth = (preStepPoint->GetTouchable()->GetReplicaNumber(0))%10; // Old/Default
00894 
00895     for (unsigned int i=0; i<hits.size(); ++i)
00896     {
00897       G4ThreeVector hitPoint = hits[i].position;
00898       if (applyFidCut) depth = hits[i].depth;
00899       double time            = hits[i].time;
00900       unsigned int unitID = setDetUnitId(det, hitPoint, depth);
00901       currentID.setID(unitID, time, primaryID, 0);
00902 #ifdef DebugLog
00903       plotProfile(aStep, hitPoint, edepositEM, time, depth);
00904 #endif
00905       // check if it is in the same unit and timeslice as the previous one
00906       if (currentID == previousID) updateHit(currentHit);
00907       else
00908       {
00909         posGlobal = preStepPoint->GetPosition();
00910         if (!checkHit()) currentHit = createNewHit();
00911       }
00912     }
00913   }
00914 }
00915 
00916 void HCalSD::getFromParam (G4Step* aStep)
00917 {
00918   std::vector<HFShowerParam::Hit> hits = showerParam->getHits(aStep);
00919   int nHit = static_cast<int>(hits.size());
00920 
00921   if (nHit > 0)
00922   {
00923     preStepPoint  = aStep->GetPreStepPoint();
00924     int primaryID = setTrackID(aStep);
00925    
00926     int det   = 5;
00927 #ifdef DebugLog
00928     edm::LogInfo("HcalSim") << "HCalSD::getFromParam " << nHit << " hits for " 
00929                             << GetName() << " of " << primaryID << " with " 
00930                             <<  aStep->GetTrack()->GetDefinition()->GetParticleName()
00931                             << " of " << preStepPoint->GetKineticEnergy()/GeV 
00932                             << " GeV in detector type " << det;
00933 #endif
00934     for (int i = 0; i<nHit; ++i)
00935     {
00936       G4ThreeVector hitPoint = hits[i].position;
00937       int depth              = hits[i].depth;
00938       double time            = hits[i].time;
00939       unsigned int unitID    = setDetUnitId(det, hitPoint, depth);
00940       currentID.setID(unitID, time, primaryID, 0);
00941       edepositEM             = hits[i].edep*GeV; 
00942       edepositHAD            = 0.;
00943 #ifdef DebugLog
00944       plotProfile(aStep, hitPoint, edepositEM, time, depth);
00945 #endif
00946 
00947       // check if it is in the same unit and timeslice as the previous one
00948       if (currentID == previousID) updateHit(currentHit);
00949       else
00950       {
00951         posGlobal = preStepPoint->GetPosition();
00952         if (!checkHit()) currentHit = createNewHit();
00953       }
00954     }
00955   }
00956 }
00957 
00958 void HCalSD::getHitPMT (G4Step * aStep)
00959 {
00960   preStepPoint = aStep->GetPreStepPoint();
00961   theTrack     = aStep->GetTrack();
00962   double edep  = showerPMT->getHits(aStep);
00963 
00964   if (edep >= 0)
00965   {
00966     double etrack    = preStepPoint->GetKineticEnergy();
00967     int    primaryID = 0;
00968     if (etrack >= energyCut) primaryID    = theTrack->GetTrackID();
00969     else
00970     {
00971       primaryID    = theTrack->GetParentID();
00972       if (primaryID == 0) primaryID = theTrack->GetTrackID();
00973     }
00974     // Reset entry point for new primary
00975     posGlobal = preStepPoint->GetPosition();
00976     resetForNewPrimary(posGlobal, etrack);
00977 
00978     //
00979     int    det      = static_cast<int>(HcalForward);
00980     G4ThreeVector hitPoint = preStepPoint->GetPosition();   
00981     double rr       = (hitPoint.x()*hitPoint.x() + hitPoint.y()*hitPoint.y());
00982     double phi      = (rr == 0. ? 0. :atan2(hitPoint.y(),hitPoint.x()));
00983     double etaR     = showerPMT->getRadius();
00984     int depth       = 1;
00985     if (etaR < 0)
00986     {
00987       depth         = 2;
00988       etaR          =-etaR;
00989     }
00990     if (hitPoint.z() < 0) etaR =-etaR;
00991 #ifdef DebugLog
00992     edm::LogInfo("HcalSim") << "HCalSD::Hit for Detector " << det << " etaR "
00993                             << etaR << " phi " << phi/deg << " depth " <<depth;
00994 #endif
00995     double time = (aStep->GetPostStepPoint()->GetGlobalTime());
00996     uint32_t unitID = 0;
00997     if (numberingFromDDD)
00998     {
00999       HcalNumberingFromDDD::HcalID tmp = numberingFromDDD->unitID(det,etaR,phi,
01000           depth,1);
01001       if (numberingScheme) unitID = numberingScheme->getUnitID(tmp);
01002     }
01003     currentID.setID(unitID, time, primaryID, 1);
01004 
01005     edepositHAD = aStep->GetTotalEnergyDeposit();
01006     edepositEM  =-edepositHAD + (edep*GeV);
01007 #ifdef DebugLog
01008     plotProfile(aStep, hitPoint, edep*GeV, time, depth);
01009     double beta = preStepPoint->GetBeta();
01010     LogDebug("HcalSim") << "HCalSD::getHitPMT 1 hit for " << GetName() 
01011                         << " of " << primaryID << " with " 
01012                         << theTrack->GetDefinition()->GetParticleName()
01013                         << " of " << preStepPoint->GetKineticEnergy()/GeV 
01014                         << " GeV with velocity " << beta << " UnitID "
01015                         << std::hex << unitID << std::dec;
01016 #endif
01017     // check if it is in the same unit and timeslice as the previous one
01018     if      (currentID == previousID) updateHit(currentHit);
01019     else if (!checkHit())             currentHit = createNewHit();
01020   }
01021 }
01022 
01023 void HCalSD::getHitFibreBundle (G4Step* aStep, bool type)
01024 {
01025   preStepPoint = aStep->GetPreStepPoint();
01026   theTrack     = aStep->GetTrack();
01027   double edep  = showerBundle->getHits(aStep, type);
01028 
01029   if (edep >= 0)
01030   {
01031     double etrack    = preStepPoint->GetKineticEnergy();
01032     int    primaryID = 0;
01033     if (etrack >= energyCut) primaryID    = theTrack->GetTrackID();
01034     else
01035     {
01036       primaryID    = theTrack->GetParentID();
01037       if (primaryID == 0) primaryID = theTrack->GetTrackID();
01038     }
01039     // Reset entry point for new primary
01040     posGlobal = preStepPoint->GetPosition();
01041     resetForNewPrimary(posGlobal, etrack);
01042 
01043     //
01044     int    det      = static_cast<int>(HcalForward);
01045     G4ThreeVector hitPoint = preStepPoint->GetPosition();   
01046     double rr       = (hitPoint.x()*hitPoint.x() + hitPoint.y()*hitPoint.y());
01047     double phi      = (rr == 0. ? 0. :atan2(hitPoint.y(),hitPoint.x()));
01048     double etaR     = showerBundle->getRadius();
01049     int depth       = 1;
01050     if (etaR < 0)
01051     {
01052       depth         = 2;
01053       etaR          =-etaR;
01054     }
01055     if (hitPoint.z() < 0) etaR =-etaR;
01056 #ifdef DebugLog
01057     edm::LogInfo("HcalSim") << "HCalSD::Hit for Detector " << det << " etaR "
01058                             << etaR << " phi " << phi/deg << " depth " <<depth;
01059 #endif
01060     double time = (aStep->GetPostStepPoint()->GetGlobalTime());
01061     uint32_t unitID = 0;
01062     if (numberingFromDDD)
01063     {
01064       HcalNumberingFromDDD::HcalID tmp = numberingFromDDD->unitID(det,etaR,phi,depth,1);
01065       if (numberingScheme) unitID = numberingScheme->getUnitID(tmp);
01066     }
01067     if (type) currentID.setID(unitID, time, primaryID, 3);
01068     else      currentID.setID(unitID, time, primaryID, 2);
01069 
01070     edepositHAD = aStep->GetTotalEnergyDeposit();
01071     edepositEM  =-edepositHAD + (edep*GeV);
01072 #ifdef DebugLog
01073     plotProfile(aStep, hitPoint, edep*GeV, time, depth);
01074     double beta = preStepPoint->GetBeta();
01075     LogDebug("HcalSim") << "HCalSD::getHitFibreBundle 1 hit for " << GetName() 
01076                         << " of " << primaryID << " with " 
01077                         << theTrack->GetDefinition()->GetParticleName()
01078                         << " of " << preStepPoint->GetKineticEnergy()/GeV 
01079                         << " GeV with velocity " << beta << " UnitID "
01080                         << std::hex << unitID << std::dec;
01081 #endif
01082     // check if it is in the same unit and timeslice as the previous one
01083     if (currentID == previousID) updateHit(currentHit);
01084     else if (!checkHit()) currentHit = createNewHit();
01085   } // non-zero energy deposit
01086 }
01087 
01088 int HCalSD::setTrackID (G4Step* aStep)
01089 {
01090   theTrack     = aStep->GetTrack();
01091 
01092   double etrack = preStepPoint->GetKineticEnergy();
01093   TrackInformation * trkInfo = (TrackInformation *)(theTrack->GetUserInformation());
01094   int      primaryID = trkInfo->getIDonCaloSurface();
01095   if (primaryID == 0)
01096   {
01097 #ifdef DebugLog
01098     edm::LogInfo("HcalSim") << "HCalSD: Problem with primaryID **** set by "
01099        << "force to TkID **** " <<theTrack->GetTrackID();
01100 #endif
01101     primaryID = theTrack->GetTrackID();
01102   }
01103 
01104   if (primaryID != previousID.trackID())
01105     resetForNewPrimary(preStepPoint->GetPosition(), etrack);
01106 
01107   return primaryID;
01108 }
01109 
01110 void HCalSD::readWeightFromFile(std::string fName)
01111 {
01112   std::ifstream infile;
01113   int entry=0;
01114   infile.open(fName.c_str(), std::ios::in);
01115   if (infile)
01116   {
01117     int    det, zside, etaR, phi, lay;
01118     double wt;
01119     while (infile >> det >> zside >> etaR >> phi >> lay >> wt)
01120     {
01121       uint32_t id = HcalTestNumbering::packHcalIndex(det,zside,1,etaR,phi,lay);
01122       layerWeights.insert(std::pair<uint32_t,double>(id,wt));
01123       ++entry;
01124 #ifdef DebugLog
01125       edm::LogInfo("HcalSim") << "HCalSD::readWeightFromFile:Entry " << entry
01126                               << " ID " << std::hex << id << std::dec << " ("
01127                               << det << "/" << zside << "/1/" << etaR << "/"
01128                               << phi << "/" << lay << ") Weight " << wt;
01129 #endif
01130     }
01131     infile.close();
01132   }
01133   edm::LogInfo("HcalSim") << "HCalSD::readWeightFromFile: reads " << entry
01134                           << " weights from " << fName;
01135   if (entry <= 0) useLayerWt = false;
01136 }
01137 
01138 double HCalSD::layerWeight(int det, G4ThreeVector pos, int depth, int lay)
01139 { 
01140   double wt = 1.;
01141   if (numberingFromDDD)
01142   {
01143     //get the ID's as eta, phi, depth, ... indices
01144     HcalNumberingFromDDD::HcalID tmp = numberingFromDDD->unitID(det, pos, depth, lay);
01145     uint32_t id = HcalTestNumbering::packHcalIndex(tmp.subdet, tmp.zside, 1,
01146                                                    tmp.etaR, tmp.phis, tmp.lay);
01147     std::map<uint32_t,double>::const_iterator ite = layerWeights.find(id);
01148     if (ite != layerWeights.end()) wt = ite->second;
01149 #ifdef DebugLog
01150     edm::LogInfo("HcalSim") << "HCalSD::layerWeight: ID " << std::hex << id 
01151                             << std::dec << " (" << tmp.subdet << "/"  
01152                             << tmp.zside << "/1/" << tmp.etaR << "/" 
01153                             << tmp.phis << "/"  << tmp.lay << ") Weight " <<wt;
01154 #endif
01155   }
01156   return wt;
01157 }
01158 
01159 void HCalSD::plotProfile(G4Step* aStep, G4ThreeVector global, double edep,
01160                          double time, int id)
01161 { 
01162 
01163   const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
01164   static G4String modName[8] = {"HEModule", "HVQF" , "HBModule", "MBAT",
01165                                 "MBBT"    , "MBBTC", "MBBT_R1P", "MBBT_R1M"};
01166   G4ThreeVector local;
01167   bool found=false;
01168   double depth=-2000;
01169   int idx = 4;
01170   for (int n=0; n<touch->GetHistoryDepth(); ++n)
01171   {
01172     G4String name = touch->GetVolume(n)->GetName();
01173 #ifdef DebugLog
01174     LogDebug("HcalSim") << "plotProfile Depth " << n << " Name " << name;
01175 #endif
01176     for (unsigned int ii=0; ii<8; ++ii)
01177       if (name == modName[ii])
01178     {
01179       found = true;
01180       int dn = touch->GetHistoryDepth() - n;
01181       local = touch->GetHistory()->GetTransform(dn).TransformPoint(global);
01182       if      (ii == 0) {depth = local.z() - 4006.5; idx = 1;}
01183       else if (ii == 1) {depth = local.z() + 825.0;  idx = 3;}
01184       else if (ii == 2) {depth = local.x() - 1775.;  idx = 0;}
01185       else              {depth = local.y() + 15.;    idx = 2;}
01186       break;
01187     }
01188     if (found) break;
01189   }
01190   if (!found) depth = std::abs(global.z()) - 11500;
01191 #ifdef DebugLog
01192   LogDebug("HcalSim") << "plotProfile Found " << found << " Global " << global
01193                       << " Local " << local << " depth " << depth << " ID " << id
01194                       << " EDEP " << edep << " Time " << time;
01195 #endif
01196   if (hit_[idx]  != 0) hit_[idx]->Fill(edep);
01197   if (time_[idx] != 0) time_[idx]->Fill(time,edep);
01198   if (dist_[idx] != 0) dist_[idx]->Fill(depth,edep);
01199   int jd = 2*idx + id - 7;
01200   if (jd >= 0 && jd < 4)
01201   {
01202     jd += 5;
01203     if (hit_[jd]  != 0) hit_[jd]->Fill(edep);
01204     if (time_[jd] != 0) time_[jd]->Fill(time,edep);
01205     if (dist_[jd] != 0) dist_[jd]->Fill(depth,edep);
01206   }
01207 }