CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch9/src/SimG4CMS/HcalTestBeam/src/HcalTB06Analysis.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 // Package:     HcalTestBeam
00004 // Class  :     HcalTB06Analysis
00005 //
00006 // Implementation:
00007 //     Main analysis class for Hcal Test Beam 2006 Analysis
00008 //
00009 // Original Author:
00010 //         Created:  Tue May 16 10:14:34 CEST 2006
00011 // $Id: HcalTB06Analysis.cc,v 1.7 2009/05/24 15:31:15 fabiocos Exp $
00012 //
00013   
00014 // system include files
00015 #include <cmath>
00016 #include <iostream>
00017 #include <iomanip>
00018 
00019 // user include files
00020 #include "SimG4CMS/HcalTestBeam/interface/HcalTB06Analysis.h"
00021 
00022 #include "SimG4Core/Notification/interface/BeginOfRun.h"
00023 #include "SimG4Core/Notification/interface/BeginOfEvent.h"
00024 #include "SimG4Core/Notification/interface/EndOfEvent.h"
00025 
00026 // to retreive hits
00027 #include "DataFormats/HcalDetId/interface/HcalSubdetector.h"
00028 #include "SimG4CMS/Calo/interface/CaloG4Hit.h"
00029 #include "SimG4CMS/Calo/interface/CaloG4HitCollection.h"
00030 #include "DataFormats/Math/interface/Point3D.h"
00031 
00032 #include "FWCore/Framework/interface/Event.h"
00033 #include "FWCore/Framework/interface/EventSetup.h"
00034 #include "FWCore/Framework/interface/ESHandle.h"
00035 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00036 
00037 #include "G4SDManager.hh"
00038 #include "G4VProcess.hh"
00039 #include "G4HCofThisEvent.hh"
00040 #include "CLHEP/Units/GlobalSystemOfUnits.h"
00041 #include "CLHEP/Units/GlobalPhysicalConstants.h"
00042 
00043 //
00044 // constructors and destructor
00045 //
00046 
00047 HcalTB06Analysis::HcalTB06Analysis(const edm::ParameterSet &p): histo(0) {
00048 
00049   edm::ParameterSet m_Anal = p.getParameter<edm::ParameterSet>("HcalTB06Analysis");
00050   names          = m_Anal.getParameter<std::vector<std::string> >("Names");
00051   beamOffset     =-m_Anal.getParameter<double>("BeamPosition")*cm;
00052   double fMinEta = m_Anal.getParameter<double>("MinEta");
00053   double fMaxEta = m_Anal.getParameter<double>("MaxEta");
00054   double fMinPhi = m_Anal.getParameter<double>("MinPhi");
00055   double fMaxPhi = m_Anal.getParameter<double>("MaxPhi");
00056   double beamEta = (fMaxEta+fMinEta)/2.;
00057   double beamPhi = (fMaxPhi+fMinPhi)/2.;
00058   double beamThet= 2*atan(exp(-beamEta));
00059   if (beamPhi < 0) beamPhi += twopi;
00060   iceta          = (int)(beamEta/0.087) + 1;
00061   icphi          = (int)(fabs(beamPhi)/0.087) + 5;
00062   if (icphi > 72) icphi -= 73;
00063 
00064   produces<PHcalTB06Info>();
00065 
00066   beamline_RM = new G4RotationMatrix;
00067   beamline_RM->rotateZ(-beamPhi);
00068   beamline_RM->rotateY(-beamThet);
00069  
00070   edm::LogInfo("HcalTBSim") << "HcalTB06:: Initialised as observer of BeginOf"
00071                             << "Job/BeginOfRun/BeginOfEvent/G4Step/EndOfEvent"
00072                             << " with Parameter values:\n \tbeamOffset = " 
00073                             << beamOffset << "\ticeta = " << iceta 
00074                             << "\ticphi = " << icphi << "\n\tbeamline_RM = "
00075                             << *beamline_RM;
00076 
00077   init();
00078 
00079   histo  = new HcalTB06Histo(m_Anal);
00080 } 
00081    
00082 HcalTB06Analysis::~HcalTB06Analysis() {
00083 
00084   edm::LogInfo("HcalTBSim") << "\n -------->  Total number of selected entries"
00085                             << " : " << count << "\nPointers:: Histo " <<histo;
00086   if (histo)   {
00087     delete histo;
00088     histo  = 0;
00089   }
00090 }
00091 
00092 //
00093 // member functions
00094 //
00095 
00096 void HcalTB06Analysis::produce(edm::Event& e, const edm::EventSetup&) {
00097 
00098   std::auto_ptr<PHcalTB06Info> product(new PHcalTB06Info);
00099   fillEvent(*product);
00100   e.put(product);
00101 }
00102 
00103 void HcalTB06Analysis::init() {
00104 
00105   // counter 
00106   count = 0;
00107   evNum = 0;
00108   clear();
00109 }
00110 
00111 void HcalTB06Analysis::update(const BeginOfRun * run) {
00112 
00113   int irun = (*run)()->GetRunID();
00114   edm::LogInfo("HcalTBSim") <<" =====> Begin of Run = " << irun;
00115  
00116 }
00117 
00118 void HcalTB06Analysis::update(const BeginOfEvent * evt) {
00119  
00120   evNum = (*evt) ()->GetEventID ();
00121   clear();
00122   edm::LogInfo("HcalTBSim") << "HcalTB06Analysis: =====> Begin of event = "
00123                             << evNum;
00124 }
00125 
00126 void HcalTB06Analysis::update(const G4Step * aStep) {
00127 
00128   if (aStep != NULL) {
00129     //Get Step properties
00130     G4ThreeVector thePreStepPoint  = aStep->GetPreStepPoint()->GetPosition();
00131     G4ThreeVector thePostStepPoint;
00132 
00133     // Get Tracks properties
00134     G4Track*      aTrack   = aStep->GetTrack();
00135     int           trackID  = aTrack->GetTrackID();
00136     int           parentID = aTrack->GetParentID();
00137     G4ThreeVector position = aTrack->GetPosition();
00138     G4ThreeVector momentum = aTrack->GetMomentum();
00139     G4String      partType = aTrack->GetDefinition()->GetParticleType();
00140     G4String      partSubType = aTrack->GetDefinition()->GetParticleSubType();
00141     int    partPDGEncoding = aTrack->GetDefinition()->GetPDGEncoding();
00142 #ifdef ddebug
00143     bool   isPDGStable = aTrack->GetDefinition()->GetPDGStable();
00144 #endif
00145     double pDGlifetime = aTrack->GetDefinition()->GetPDGLifeTime();
00146     double gammaFactor = aStep->GetPreStepPoint()->GetGamma();
00147 
00148     if (!pvFound) { //search for v1
00149       double stepDeltaEnergy = aStep->GetDeltaEnergy ();
00150       double kinEnergy = aTrack->GetKineticEnergy ();
00151       
00152       // look for DeltaE > 10% kinEnergy of particle, or particle death - Ek=0
00153       if (trackID == 1 && parentID == 0 && 
00154           ((kinEnergy == 0.) || (fabs (stepDeltaEnergy / kinEnergy) > 0.1))) {
00155         int pvType = -1;
00156         if (kinEnergy == 0.) {
00157           pvType = 0;
00158         } else {
00159           if (fabs (stepDeltaEnergy / kinEnergy) > 0.1) pvType = 1;
00160         }
00161         pvFound    = true;
00162         pvPosition = position;
00163         pvMomentum = momentum;
00164         // Rotated coord.system:
00165         pvUVW      = (*beamline_RM)*(pvPosition);
00166 
00167         //Volume name requires some checks:
00168         G4String thePostPVname = "NoName";
00169         G4StepPoint * thePostPoint = aStep->GetPostStepPoint ();
00170         if (thePostPoint) {
00171           thePostStepPoint = thePostPoint->GetPosition();
00172           G4VPhysicalVolume * thePostPV = thePostPoint->GetPhysicalVolume ();
00173           if (thePostPV) thePostPVname = thePostPV->GetName ();
00174         }
00175 #ifdef ddebug
00176         LogDebug("HcalTBSim") << "HcalTB06Analysis:: V1 found at: " 
00177                               << thePostStepPoint << " G4VPhysicalVolume: " 
00178                               << thePostPVname;
00179 #endif      
00180         LogDebug("HcalTBSim") << "HcalTB06Analysis::fill_v1Pos: Primary Track "
00181                               << "momentum: " << pvMomentum << " psoition " 
00182                               << pvPosition << " u/v/w " << pvUVW;
00183       }
00184     } else { 
00185       // watch for secondaries originating @v1, including the surviving primary
00186       if ((trackID != 1 && parentID == 1 &&
00187            (aTrack->GetCurrentStepNumber () == 1) && 
00188            (thePreStepPoint == pvPosition)) || 
00189           (trackID == 1 && thePreStepPoint == pvPosition)) {
00190 #ifdef ddebug
00191         LogDebug("HcalTBSim") << "HcalTB06Analysis::A secondary...  PDG:" 
00192                               << partPDGEncoding << " TrackID:" << trackID
00193                               << " ParentID:" << parentID << " stable: "  
00194                               << isPDGStable << " Tau: " << pDGlifetime 
00195                               << " cTauGamma=" 
00196                               << c_light*pDGlifetime*gammaFactor*1000.
00197                               << "um" << " GammaFactor: " << gammaFactor;
00198 #endif      
00199         secTrackID.push_back(trackID);
00200         secPartID.push_back(partPDGEncoding);
00201         secMomentum.push_back(momentum);
00202         secEkin.push_back(aTrack->GetKineticEnergy());
00203 
00204         // Check for short-lived secondaries: cTauGamma<100um
00205         double ctaugamma_um = c_light * pDGlifetime * gammaFactor * 1000.;
00206         if ((ctaugamma_um>0.) && (ctaugamma_um<100.)) {//short-lived secondary
00207           shortLivedSecondaries.push_back(trackID);
00208       } else {//normal secondary - enter into the V1-calorimetric tree
00209         //          histos->fill_v1cSec (aTrack);
00210       }
00211       }
00212       // Also watch for tertiary particles coming from 
00213       // short-lived secondaries from V1
00214       if (aTrack->GetCurrentStepNumber() == 1) {
00215         if (shortLivedSecondaries.size() > 0) {
00216           int pid = parentID;
00217           std::vector<int>::iterator pos1= shortLivedSecondaries.begin();
00218           std::vector<int>::iterator pos2 = shortLivedSecondaries.end();
00219           std::vector<int>::iterator pos;
00220           for (pos = pos1; pos != pos2; pos++) {
00221             if (*pos == pid) {//ParentID is on the list of short-lived 
00222               // secondary 
00223 #ifdef ddebug
00224               LogDebug("HcalTBSim") << "HcalTB06Analysis:: A tertiary...  PDG:"
00225                                     << partPDGEncoding << " TrackID:" <<trackID
00226                                     << " ParentID:" << parentID << " stable: "
00227                                     << isPDGStable << " Tau: " << pDGlifetime
00228                                     << " cTauGamma=" 
00229                                     << c_light*pDGlifetime*gammaFactor*1000. 
00230                                     << "um GammaFactor: " << gammaFactor;
00231 #endif
00232             }
00233           }
00234         }
00235       }
00236     }
00237   }
00238 }
00239 
00240 void HcalTB06Analysis::update(const EndOfEvent * evt) {
00241 
00242   count++;
00243 
00244   //fill the buffer
00245   LogDebug("HcalTBSim") << "HcalTB06Analysis::Fill event " 
00246                         << (*evt)()->GetEventID();
00247   fillBuffer (evt);
00248   
00249   //Final Analysis
00250   LogDebug("HcalTBSim") << "HcalTB06Analysis::Final analysis";  
00251   finalAnalysis();
00252 
00253   int iEvt = (*evt)()->GetEventID();
00254   if (iEvt < 10) 
00255     edm::LogInfo("HcalTBSim") << "HcalTB06Analysis:: Event " << iEvt;
00256   else if ((iEvt < 100) && (iEvt%10 == 0)) 
00257     edm::LogInfo("HcalTBSim") << "HcalTB06Analysis:: Event " << iEvt;
00258   else if ((iEvt < 1000) && (iEvt%100 == 0)) 
00259     edm::LogInfo("HcalTBSim") << "HcalTB06Analysis:: Event " << iEvt;
00260   else if ((iEvt < 10000) && (iEvt%1000 == 0)) 
00261     edm::LogInfo("HcalTBSim") << "HcalTB06Analysis:: Event " << iEvt;
00262 }
00263 
00264 void HcalTB06Analysis::fillBuffer(const EndOfEvent * evt) {
00265 
00266   std::vector<CaloHit> hhits;
00267   int                  idHC, j;
00268   CaloG4HitCollection* theHC;
00269   std::map<int,float,std::less<int> > primaries;
00270   double               etot1=0, etot2=0;
00271 
00272   // Look for the Hit Collection of HCal
00273   G4HCofThisEvent* allHC = (*evt)()->GetHCofThisEvent();
00274   std::string sdName = names[0];
00275   idHC  = G4SDManager::GetSDMpointer()->GetCollectionID(sdName);
00276   theHC = (CaloG4HitCollection*) allHC->GetHC(idHC);
00277   LogDebug("HcalTBSim") << "HcalTB06Analysis:: Hit Collection for " << sdName
00278                         << " of ID " << idHC << " is obtained at " << theHC;
00279 
00280   if (idHC >= 0 && theHC > 0) {
00281     hhits.reserve(theHC->entries());
00282     for (j = 0; j < theHC->entries(); j++) {
00283       CaloG4Hit* aHit = (*theHC)[j]; 
00284       double e        = aHit->getEnergyDeposit()/GeV;
00285       double time     = aHit->getTimeSlice();
00286       math::XYZPoint pos  = aHit->getEntry();
00287       unsigned int id = aHit->getUnitID();
00288       double theta    = pos.theta();
00289       double eta      = -log(tan(theta * 0.5));
00290       double phi      = pos.phi();
00291       CaloHit hit(2,1,e,eta,phi,time,id);
00292       hhits.push_back(hit);
00293       primaries[aHit->getTrackID()]+= e;
00294       etot1 += e;
00295 #ifdef ddebug
00296       LogDebug("HcalTBSim") << "HcalTB06Analysis:: Hcal Hit i/p " << j 
00297                             << "  ID 0x" << std::hex << id << std::dec 
00298                             << " time " << std::setw(6) << time << " theta "
00299                             << std::setw(8) << theta << " eta " << std::setw(8)
00300                             << eta << " phi " << std::setw(8) << phi << " e " 
00301                             << std::setw(8) << e;
00302 #endif
00303     }
00304   }
00305 
00306   // Add hits in the same channel within same time slice
00307   std::vector<CaloHit>::iterator itr;
00308   int nHit = hhits.size();
00309   std::vector<CaloHit*> hits(nHit);
00310   for (j = 0, itr = hhits.begin(); itr != hhits.end(); j++, itr++) {
00311     hits[j] = &hhits[j];
00312   }
00313   sort(hits.begin(),hits.end(),CaloHitIdMore());
00314   std::vector<CaloHit*>::iterator k1, k2;
00315   int nhit = 0;
00316   for (k1 = hits.begin(); k1 != hits.end(); k1++) {
00317     int      det    = (**k1).det();
00318     int      layer  = (**k1).layer();
00319     double   ehit   = (**k1).e();
00320     double   eta    = (**k1).eta();
00321     double   phi    = (**k1).phi();
00322     double   jitter = (**k1).t();
00323     uint32_t unitID = (**k1).id();
00324     int      jump  = 0;
00325     for (k2 = k1+1; k2 != hits.end() && fabs(jitter-(**k2).t())<1 &&
00326            unitID==(**k2).id(); k2++) {
00327       ehit += (**k2).e();
00328       jump++;
00329     }
00330     nhit++;
00331     CaloHit hit(det, layer, ehit, eta, phi, jitter, unitID);
00332     hcalHitCache.push_back(hit);
00333     etot2 += ehit;
00334     k1    += jump;
00335 #ifdef ddebug
00336     LogDebug("HcalTBSim") << "HcalTB06Analysis:: Hcal Hit store " << nhit 
00337                           << "  ID 0x" << std::hex  << unitID  << std::dec 
00338                           << " time " << std::setw(6) << jitter << " eta "
00339                           << std::setw(8) << eta << " phi " << std::setw(8) 
00340                           << phi  << " e " << std::setw(8) << ehit;
00341 #endif
00342   }
00343   LogDebug("HcalTBSim") << "HcalTB06Analysis:: Stores " << nhit << " HCal hits"
00344                         << " from " << nHit << " input hits E(Hcal) " << etot1 
00345                         << " " << etot2;
00346   
00347   // Look for the Hit Collection of ECal
00348   std::vector<CaloHit> ehits;
00349   sdName= names[1];
00350   idHC  = G4SDManager::GetSDMpointer()->GetCollectionID(sdName);
00351   theHC = (CaloG4HitCollection*) allHC->GetHC(idHC);
00352   etot1 = etot2 = 0;
00353   LogDebug("HcalTBSim") << "HcalTB06Analysis:: Hit Collection for " << sdName
00354                         << " of ID " << idHC << " is obtained at " << theHC;
00355   if (idHC >= 0 && theHC > 0) {
00356     ehits.reserve(theHC->entries());
00357     for (j = 0; j < theHC->entries(); j++) {
00358       CaloG4Hit* aHit = (*theHC)[j]; 
00359       double e        = aHit->getEnergyDeposit()/GeV;
00360       double time     = aHit->getTimeSlice();
00361       math::XYZPoint pos  = aHit->getEntry();
00362       unsigned int id = aHit->getUnitID();
00363       double theta    = pos.theta();
00364       double eta      = -log(tan(theta * 0.5));
00365       double phi      = pos.phi();
00366       if (e < 0 || e > 100000.) e = 0;
00367       CaloHit hit(1,0,e,eta,phi,time,id);
00368       ehits.push_back(hit);
00369       primaries[aHit->getTrackID()]+= e;
00370       etot1 += e;
00371 #ifdef ddebug
00372       LogDebug("HcalTBSim") << "HcalTB06Analysis:: Ecal Hit i/p " << j 
00373                             << "  ID 0x" << std::hex << id  << std::dec 
00374                             << " time " << std::setw(6) << time << " theta " 
00375                             << std::setw(8) << theta  << " eta " <<std::setw(8)
00376                             << eta  << " phi " << std::setw(8) << phi << " e "
00377                             << std::setw(8) << e;
00378 #endif
00379     }
00380   }
00381 
00382   // Add hits in the same channel within same time slice
00383   nHit = ehits.size();
00384   std::vector<CaloHit*> hite(nHit);
00385   for (j = 0, itr = ehits.begin(); itr != ehits.end(); j++, itr++) {
00386     hite[j] = &ehits[j];
00387   }
00388   sort(hite.begin(),hite.end(),CaloHitIdMore());
00389   nhit = 0;
00390   for (k1 = hite.begin(); k1 != hite.end(); k1++) {
00391     int      det    = (**k1).det();
00392     int      layer  = (**k1).layer();
00393     double   ehit   = (**k1).e();
00394     double   eta    = (**k1).eta();
00395     double   phi    = (**k1).phi();
00396     double   jitter = (**k1).t();
00397     uint32_t unitID = (**k1).id();
00398     int      jump  = 0;
00399     for (k2 = k1+1; k2 != hite.end() && fabs(jitter-(**k2).t())<1 &&
00400            unitID==(**k2).id(); k2++) {
00401       ehit += (**k2).e();
00402       jump++;
00403     }
00404     nhit++;
00405     CaloHit hit(det, layer, ehit, eta, phi, jitter, unitID);
00406     ecalHitCache.push_back(hit);
00407     etot2 += ehit;
00408     k1    += jump;
00409 #ifdef ddebug
00410     LogDebug("HcalTBSim") << "HcalTB06Analysis:: Ecal Hit store " << nhit
00411                           << "  ID 0x" << std::hex << unitID  << std::dec 
00412                           << " time " << std::setw(6) << jitter << " eta "
00413                           << std::setw(8) << eta << " phi " << std::setw(8)
00414                           << phi << " e " << std::setw(8) << ehit;
00415 #endif
00416   }
00417   LogDebug("HcalTBSim") << "HcalTB06Analysis:: Stores " << nhit << " ECal hits"
00418                         << " from " << nHit << " input hits E(Ecal) " << etot1 
00419                         << " " << etot2;
00420 
00421   // Find Primary info:
00422   nPrimary    = (int)(primaries.size());
00423   int trackID = 0;
00424   G4PrimaryParticle* thePrim=0;
00425   int nvertex = (*evt)()->GetNumberOfPrimaryVertex();
00426   LogDebug("HcalTBSim") << "HcalTB06Analysis:: Event has " << nvertex 
00427                         << " verteices";
00428   if (nvertex<=0)
00429     edm::LogInfo("HcalTBSim") << "HcalTB06Analysis::EndOfEvent ERROR: no "
00430                               << "vertex found for event " << evNum;
00431 
00432   for (int i = 0 ; i<nvertex; i++) {
00433     G4PrimaryVertex* avertex = (*evt)()->GetPrimaryVertex(i);
00434     if (avertex == 0) {
00435       edm::LogInfo("HcalTBSim") << "HcalTB06Analysis::EndOfEvent ERR: pointer "
00436                                 << "to vertex = 0 for event " << evNum;
00437     } else {
00438       LogDebug("HcalTBSim") << "HcalTB06Analysis::Vertex number :" << i << " "
00439                             << avertex->GetPosition();
00440       int npart = avertex->GetNumberOfParticle();
00441       if (npart == 0)
00442         edm::LogWarning("HcalTBSim") << "HcalTB06Analysis::End Of Event ERR: "
00443                                      << "no primary!";
00444       if (thePrim==0) thePrim=avertex->GetPrimary(trackID);
00445     }
00446   }
00447     
00448   if (thePrim != 0) {
00449     double px = thePrim->GetPx();
00450     double py = thePrim->GetPy();
00451     double pz = thePrim->GetPz();
00452     double p  = std::sqrt(pow(px,2.)+pow(py,2.)+pow(pz,2.));
00453     pInit     = p/GeV;
00454     if (p==0) 
00455       edm::LogWarning("HcalTBSim") << "HcalTB06Analysis:: EndOfEvent ERR: "
00456                                    << "primary has p=0 ";
00457     else {
00458       double costheta = pz/p;
00459       double theta = acos(std::min(std::max(costheta,-1.),1.));
00460       etaInit = -log(tan(theta/2));
00461       if (px != 0 || py != 0) phiInit = atan2(py,px);  
00462     }
00463     particleType = thePrim->GetPDGcode();
00464   } else 
00465     edm::LogWarning("HcalTBSim") << "HcalTB06Analysis::EndOfEvent ERR: could "
00466                                  << "not find primary";
00467 
00468 }
00469 
00470 void HcalTB06Analysis::finalAnalysis() {
00471 
00472   //Beam Information
00473   histo->fillPrimary(pInit, etaInit, phiInit);
00474 
00475   // Total Energy
00476   eecals = ehcals = 0.;
00477   for (unsigned int i=0; i<hcalHitCache.size(); i++) {
00478     ehcals += hcalHitCache[i].e();
00479   }
00480   for (unsigned int i=0; i<ecalHitCache.size(); i++) {
00481     eecals += ecalHitCache[i].e();
00482   }
00483   etots = eecals + ehcals;
00484   LogDebug("HcalTBSim") << "HcalTB06Analysis:: Energy deposit at Sim Level "
00485                         << "(Total) " << etots << " (ECal) " << eecals 
00486                         << " (HCal) " << ehcals;
00487   histo->fillEdep(etots, eecals, ehcals);
00488 }
00489 
00490 
00491 void HcalTB06Analysis::fillEvent (PHcalTB06Info& product) {
00492 
00493   //Beam Information
00494   product.setPrimary(nPrimary, particleType, pInit, etaInit, phiInit);
00495 
00496   // Total Energy
00497   product.setEdep(etots, eecals, ehcals);
00498 
00499   //Energy deposits in the crystals and towers
00500   for (unsigned int i=0; i<hcalHitCache.size(); i++) 
00501     product.saveHit(hcalHitCache[i].id(), hcalHitCache[i].eta(),
00502                     hcalHitCache[i].phi(), hcalHitCache[i].e(),
00503                     hcalHitCache[i].t());
00504   for (unsigned int i=0; i<ecalHitCache.size(); i++) 
00505     product.saveHit(ecalHitCache[i].id(), ecalHitCache[i].eta(),
00506                     ecalHitCache[i].phi(), ecalHitCache[i].e(),
00507                     ecalHitCache[i].t());
00508 
00509   //Vertex associated quantities
00510   product.setVtxPrim(evNum, pvType, pvPosition.x(), pvPosition.y(), 
00511                      pvPosition.z(), pvUVW.x(), pvUVW.y(), pvUVW.z(),
00512                      pvMomentum.x(), pvMomentum.y(), pvMomentum.z());
00513   for (unsigned int i = 0; i < secTrackID.size(); i++) {
00514     product.setVtxSec(secTrackID[i], secPartID[i], secMomentum[i].x(),
00515                       secMomentum[i].y(), secMomentum[i].z(), secEkin[i]);
00516   }
00517 }
00518 
00519 void HcalTB06Analysis::clear() {
00520 
00521   pvFound = false;
00522   pvType  =-2;
00523   pvPosition = G4ThreeVector();
00524   pvMomentum = G4ThreeVector();
00525   pvUVW      = G4ThreeVector();
00526   secTrackID.clear();
00527   secPartID.clear();
00528   secMomentum.clear();
00529   secEkin.clear();
00530   shortLivedSecondaries.clear();
00531 
00532   ecalHitCache.erase(ecalHitCache.begin(), ecalHitCache.end()); 
00533   hcalHitCache.erase(hcalHitCache.begin(), hcalHitCache.end()); 
00534   nPrimary = particleType = 0;
00535   pInit = etaInit = phiInit = 0;
00536 }