CMS 3D CMS Logo

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