CMS 3D CMS Logo

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