00001
00002 #include "SimG4CMS/Calo/interface/CaloG4Hit.h"
00003 #include "SimG4CMS/Calo/interface/CaloG4HitCollection.h"
00004 #include "DataFormats/Math/interface/Point3D.h"
00005
00006 #include "SimG4CMS/Forward/interface/ZdcTestAnalysis.h"
00007 #include "SimG4CMS/Forward/interface/ZdcNumberingScheme.h"
00008
00009 #include "TFile.h"
00010 #include <cmath>
00011 #include <iostream>
00012 #include <iomanip>
00013
00014 enum ntzdcs_elements {
00015 ntzdcs_evt, ntzdcs_trackid, ntzdcs_charge, ntzdcs_pdgcode, ntzdcs_x, ntzdcs_y, ntzdcs_z, ntzdcs_stepl,
00016 ntzdcs_stepe, ntzdcs_eta, ntzdcs_phi, ntzdcs_vpx, ntzdcs_vpy, ntzdcs_vpz, ntzdcs_idx, ntzdcs_idl,
00017 ntzdcs_pvtype, ntzdcs_ncherphot
00018 };
00019
00020 enum ntzdce_elements {
00021 ntzdce_evt,ntzdce_ihit,ntzdce_fiberid,ntzdce_zside,ntzdce_subdet,ntzdce_layer,ntzdce_fiber,ntzdce_channel,
00022 ntzdce_enem,ntzdce_enhad,ntzdce_hitenergy,ntzdce_x,ntzdce_y,ntzdce_z,ntzdce_time,ntzdce_etot
00023 };
00024
00025 ZdcTestAnalysis::ZdcTestAnalysis(const edm::ParameterSet &p){
00026
00027 edm::ParameterSet m_Anal = p.getParameter<edm::ParameterSet>("ZdcTestAnalysis");
00028 verbosity = m_Anal.getParameter<int>("Verbosity");
00029 doNTzdcstep = m_Anal.getParameter<int>("StepNtupleFlag");
00030 doNTzdcevent = m_Anal.getParameter<int>("EventNtupleFlag");
00031 stepNtFileName = m_Anal.getParameter<std::string>("StepNtupleFileName");
00032 eventNtFileName = m_Anal.getParameter<std::string>("EventNtupleFileName");
00033
00034 if (verbosity > 0)
00035 std::cout<<std::endl;
00036 std::cout<<"============================================================================"<<std::endl;
00037 std::cout << "ZdcTestAnalysis:: Initialized as observer"<< std::endl;
00038 if (doNTzdcstep > 0){
00039 std::cout <<" Step Ntuple will be created"<< std::endl;
00040 std::cout <<" Step Ntuple file: "<<stepNtFileName<<std::endl;
00041 }
00042 if (doNTzdcevent > 0){
00043 std::cout <<" Event Ntuple will be created"<< std::endl;
00044 std::cout <<" Step Ntuple file: "<<stepNtFileName<<std::endl;
00045 }
00046 std::cout<<"============================================================================"<<std::endl;
00047 std::cout<<std::endl;
00048
00049 if (doNTzdcstep > 0)
00050 zdcstepntuple =
00051 new TNtuple("NTzdcstep","NTzdcstep",
00052 "evt:trackid:charge:pdgcode:x:y:z:stepl:stepe:eta:phi:vpx:vpy:vpz:idx:idl:pvtype:ncherphot");
00053
00054 if (doNTzdcevent > 0)
00055 zdceventntuple =
00056 new TNtuple("NTzdcevent","NTzdcevent",
00057 "evt:ihit:fiberid:zside:subdet:layer:fiber:channel:enem:enhad:hitenergy:x:y:z:time:etot");
00058
00059
00060 }
00061
00062 ZdcTestAnalysis::~ZdcTestAnalysis() {
00063
00064 finish();
00065 if (verbosity > 0) {
00066 std::cout << std::endl << "ZdcTestAnalysis Dextructor --------> End of ZdcTestAnalysis : "
00067 << std::endl << std::endl;
00068 }
00069
00070
00071
00072
00073 std::cout<<"ZdcTestAnalysis: End of process"<<std::endl;
00074 }
00075
00076
00077 void ZdcTestAnalysis::update(const BeginOfJob * job) {
00078
00079 std::cout<<"beggining of job"<<std::endl;;
00080 }
00081
00082
00083
00084 void ZdcTestAnalysis::update(const BeginOfRun * run) {
00085
00086
00087 std::cout << std::endl << "ZdcTestAnalysis: Begining of Run"<< std::endl;
00088 if (doNTzdcstep) {
00089 std::cout << "ZDCTestAnalysis: output step file created"<< std::endl;
00090 TString stepfilename = stepNtFileName;
00091 zdcOutputStepFile = new TFile(stepfilename,"RECREATE");
00092
00093 }
00094
00095 if (doNTzdcevent) {
00096 std::cout << "ZDCTestAnalysis: output event file created"<< std::endl;
00097 TString stepfilename = eventNtFileName;
00098 zdcOutputEventFile = new TFile(stepfilename,"RECREATE");
00099 }
00100
00101 eventIndex = 0;
00102 }
00103
00104
00105
00106
00107 void ZdcTestAnalysis::update(const BeginOfEvent * evt) {
00108
00109 std::cout << "ZdcTest: Processing Event Number: "<<eventIndex<< std::endl;
00110 eventIndex++;
00111 stepIndex = 0;
00112 }
00113
00114
00115 void ZdcTestAnalysis::update(const G4Step * aStep) {
00116
00117 stepIndex++;
00118
00119 if (doNTzdcstep) {
00120 G4StepPoint * preStepPoint = aStep->GetPreStepPoint();
00121
00122 G4double stepL = aStep->GetStepLength();
00123 G4double stepE = aStep->GetTotalEnergyDeposit();
00124
00125 if (verbosity >= 2)
00126 std::cout << "Step " << stepL << "," << stepE <<std::endl;
00127
00128 G4Track * theTrack = aStep->GetTrack();
00129 G4int theTrackID = theTrack->GetTrackID();
00130 G4double theCharge = theTrack->GetDynamicParticle()->GetCharge();
00131 G4String particleType = theTrack->GetDefinition()->GetParticleName();
00132 G4int pdgcode = theTrack->GetDefinition()->GetPDGEncoding();
00133
00134 G4ThreeVector vert_mom = theTrack->GetVertexMomentumDirection();
00135 G4double vpx = vert_mom.x();
00136 G4double vpy = vert_mom.y();
00137 G4double vpz = vert_mom.z();
00138 double eta = 0.5 * log( (1.+vpz) / (1.-vpz) );
00139 double phi = atan2(vpy,vpx);
00140
00141 G4ThreeVector hitPoint = preStepPoint->GetPosition();
00142 G4ThreeVector localPoint = theTrack->GetTouchable()->GetHistory()->
00143 GetTopTransform().TransformPoint(hitPoint);
00144
00145 const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
00146 int idx = touch->GetReplicaNumber(0);
00147 int idLayer = -1;
00148 int thePVtype = -1;
00149
00150 int historyDepth = touch->GetHistoryDepth();
00151
00152 if (historyDepth > 0) {
00153 std::vector<int> theReplicaNumbers(historyDepth);
00154 std::vector<G4VPhysicalVolume*> thePhysicalVolumes(historyDepth);
00155 std::vector<G4String> thePVnames(historyDepth);
00156 std::vector<G4LogicalVolume*> theLogicalVolumes(historyDepth);
00157 std::vector<G4String> theLVnames(historyDepth);
00158 std::vector<G4Material*> theMaterials(historyDepth);
00159 std::vector<G4String> theMaterialNames(historyDepth);
00160
00161 for (int jj = 0; jj < historyDepth; jj++) {
00162 theReplicaNumbers[jj] = touch->GetReplicaNumber(jj);
00163 thePhysicalVolumes[jj] = touch->GetVolume(jj);
00164 thePVnames[jj] = thePhysicalVolumes[jj]->GetName();
00165 theLogicalVolumes[jj] = thePhysicalVolumes[jj]->GetLogicalVolume();
00166 theLVnames[jj] = theLogicalVolumes[jj]->GetName();
00167 theMaterials[jj] = theLogicalVolumes[jj]->GetMaterial();
00168 theMaterialNames[jj] = theMaterials[jj]->GetName();
00169 if (verbosity >= 2)
00170 std::cout << " GHD " << jj << ": " << theReplicaNumbers[jj] << ","
00171 << thePVnames[jj] << "," << theLVnames[jj] << ","
00172 << theMaterialNames[jj] << std::endl;
00173 }
00174
00175 idLayer = theReplicaNumbers[1];
00176 if (thePVnames[0] == "ZDC_EMLayer")
00177 thePVtype = 1;
00178 else if (thePVnames[0] == "ZDC_EMAbsorber")
00179 thePVtype = 2;
00180 else if (thePVnames[0] == "ZDC_EMFiber")
00181 thePVtype = 3;
00182 else if (thePVnames[0] == "ZDC_HadLayer")
00183 thePVtype = 7;
00184 else if (thePVnames[0] == "ZDC_HadAbsorber")
00185 thePVtype = 8;
00186 else if (thePVnames[0] == "ZDC_HadFiber")
00187 thePVtype = 9;
00188 else if (thePVnames[0] == "ZDC_PhobosLayer")
00189 thePVtype = 11;
00190 else if (thePVnames[0] == "ZDC_PhobosAbsorber")
00191 thePVtype = 12;
00192 else if (thePVnames[0] == "ZDC_PhobosFiber")
00193 thePVtype = 13;
00194 else {
00195 thePVtype = 0;
00196 if (verbosity >= 2)
00197 std::cout << " pvtype=0 hd=" << historyDepth << " " << theReplicaNumbers[0] << ","
00198 << thePVnames[0] << "," << theLVnames[0] << "," << theMaterialNames[0] << std::endl;
00199 }
00200 }
00201 else if (historyDepth == 0) {
00202 int theReplicaNumber = touch->GetReplicaNumber(0);
00203 G4VPhysicalVolume* thePhysicalVolume = touch->GetVolume(0);
00204 G4String thePVname = thePhysicalVolume->GetName();
00205 G4LogicalVolume * theLogicalVolume = thePhysicalVolume->GetLogicalVolume();
00206 G4String theLVname = theLogicalVolume->GetName();
00207 G4Material * theMaterial = theLogicalVolume->GetMaterial();
00208 G4String theMaterialName = theMaterial->GetName();
00209 if (verbosity >= 2)
00210 std::cout << " hd=0 " << theReplicaNumber << ","
00211 << thePVname << "," << theLVname << ","
00212 << theMaterialName << std::endl;
00213 }
00214 else {
00215 std::cout << " hd<0: hd=" << historyDepth << std::endl;
00216 }
00217
00218
00219 double NCherPhot = -1;
00220 zdcsteparray[ntzdcs_evt] = (float)eventIndex;
00221 zdcsteparray[ntzdcs_trackid] = (float)theTrackID;
00222 zdcsteparray[ntzdcs_charge] = theCharge;
00223 zdcsteparray[ntzdcs_pdgcode] = (float)pdgcode;
00224 zdcsteparray[ntzdcs_x] = hitPoint.x();
00225 zdcsteparray[ntzdcs_y] = hitPoint.y();
00226 zdcsteparray[ntzdcs_z] = hitPoint.z();
00227 zdcsteparray[ntzdcs_stepl] = stepL;
00228 zdcsteparray[ntzdcs_stepe] = stepE;
00229 zdcsteparray[ntzdcs_eta] = eta;
00230 zdcsteparray[ntzdcs_phi] = phi;
00231 zdcsteparray[ntzdcs_vpx] = vpx;
00232 zdcsteparray[ntzdcs_vpy] = vpy;
00233 zdcsteparray[ntzdcs_vpz] = vpz;
00234 zdcsteparray[ntzdcs_idx] = (float)idx;
00235 zdcsteparray[ntzdcs_idl] = (float)idLayer;
00236 zdcsteparray[ntzdcs_pvtype] = thePVtype;
00237 zdcsteparray[ntzdcs_ncherphot] = NCherPhot;
00238 zdcstepntuple->Fill(zdcsteparray);
00239
00240 }
00241 }
00242
00243 void ZdcTestAnalysis::update(const EndOfEvent * evt) {
00244
00245
00246
00247 std::cout << std::endl << "ZdcTest::upDate(const EndOfEvent * evt) - event #" << (*evt)()->GetEventID()
00248 << std::endl << " # of aSteps followed in event = " << stepIndex << std::endl;
00249
00250
00251 G4HCofThisEvent* allHC = (*evt)()->GetHCofThisEvent();
00252 std::cout << " accessed all HC";
00253
00254 int theZDCHCid = G4SDManager::GetSDMpointer()->GetCollectionID("ZDCHITS");
00255 std::cout << " - theZDCHCid = " << theZDCHCid;
00256
00257 CaloG4HitCollection* theZDCHC = (CaloG4HitCollection*) allHC->GetHC(theZDCHCid);
00258 std::cout << " - theZDCHC = " << theZDCHC << std::endl;
00259
00260 ZdcNumberingScheme * theZdcNumScheme = new ZdcNumberingScheme(1);
00261
00262 float ETot=0., SEnergy=0.;
00263 int maxTime=0;
00264 int fiberID=0;
00265 unsigned int unsignedfiberID=0;
00266 std::map<int,float,std::less<int> > energyInFibers;
00267 std::map<int,float,std::less<int> > primaries;
00268 float totalEnergy = 0;
00269 int nentries = theZDCHC->entries();
00270 std::cout << " theZDCHC has " << nentries << " entries" << std::endl;
00271
00272 if (doNTzdcevent) {
00273 if (nentries > 0) {
00274 for (int ihit = 0; ihit < nentries; ihit++) {
00275 CaloG4Hit* caloHit = (*theZDCHC)[ihit];
00276 totalEnergy += caloHit->getEnergyDeposit();
00277 }
00278
00279 for (int ihit = 0; ihit < nentries; ihit++) {
00280 CaloG4Hit* aHit = (*theZDCHC)[ihit];
00281 fiberID = aHit->getUnitID();
00282 unsignedfiberID = aHit->getUnitID();
00283 double enEm = aHit->getEM();
00284 double enHad = aHit->getHadr();
00285 math::XYZPoint hitPoint = aHit->getPosition();
00286 double hitEnergy = aHit->getEnergyDeposit();
00287 if (verbosity >= 1)
00288 std::cout << " entry #" << ihit << ": fiberID=0x" << std::hex
00289 << fiberID << std::dec << "; enEm=" << enEm
00290 << "; enHad=" << enHad << "; hitEnergy=" << hitEnergy
00291 << "z=" << hitPoint.z() << std::endl;
00292 energyInFibers[fiberID]+= enEm + enHad;
00293 primaries[aHit->getTrackID()]+= enEm + enHad;
00294 float time = aHit->getTimeSliceID();
00295 if (time > maxTime) maxTime=(int)time;
00296
00297 int thesubdet, thelayer, thefiber, thechannel, thez;
00298 theZdcNumScheme->unpackZdcIndex(fiberID, thesubdet, thelayer, thefiber, thechannel, thez);
00299 int unsignedsubdet, unsignedlayer, unsignedfiber, unsignedchannel, unsignedz;
00300 theZdcNumScheme->unpackZdcIndex(unsignedfiberID, unsignedsubdet,
00301 unsignedlayer, unsignedfiber, unsignedchannel, unsignedz);
00302
00303
00304
00305
00306
00307
00308 zdceventarray[ntzdce_evt] = (float)eventIndex;
00309 zdceventarray[ntzdce_ihit] = (float)ihit;
00310 zdceventarray[ntzdce_fiberid] = (float)fiberID;
00311 zdceventarray[ntzdce_zside] = (float)thez;
00312 zdceventarray[ntzdce_subdet] = (float)thesubdet;
00313 zdceventarray[ntzdce_layer] = (float)thelayer;
00314 zdceventarray[ntzdce_fiber] = (float)thefiber;
00315 zdceventarray[ntzdce_channel] = (float)thechannel;
00316 zdceventarray[ntzdce_enem] = enEm;
00317 zdceventarray[ntzdce_enhad] = enHad;
00318 zdceventarray[ntzdce_hitenergy] = hitEnergy;
00319 zdceventarray[ntzdce_x] = hitPoint.x();
00320 zdceventarray[ntzdce_y] = hitPoint.y();
00321 zdceventarray[ntzdce_z] = hitPoint.z();
00322 zdceventarray[ntzdce_time] = time;
00323 zdceventarray[ntzdce_etot] = totalEnergy;
00324 zdceventntuple->Fill(zdceventarray);
00325 }
00326
00327 for (std::map<int,float,std::less<int> >::iterator is = energyInFibers.begin();
00328 is!= energyInFibers.end(); is++)
00329 {
00330 ETot = (*is).second;
00331 SEnergy += ETot;
00332 }
00333
00334
00335 int trackID = 0;
00336 G4PrimaryParticle* thePrim=0;
00337 G4int nvertex = (*evt)()->GetNumberOfPrimaryVertex();
00338 std::cout << "Event has " << nvertex << " vertex" << std::endl;
00339 if (nvertex==0)
00340 std::cout << "ZdcTest End Of Event ERROR: no vertex" << std::endl;
00341
00342 for (int i = 0 ; i<nvertex; i++) {
00343
00344 G4PrimaryVertex* avertex = (*evt)()->GetPrimaryVertex(i);
00345 if (avertex == 0)
00346 std::cout << "ZdcTest End Of Event ERR: pointer to vertex = 0"
00347 << std::endl;
00348 std::cout << "Vertex number :" <<i << std::endl;
00349 int npart = avertex->GetNumberOfParticle();
00350 if (npart ==0)
00351 std::cout << "ZdcTest End Of Event ERR: no primary!" << std::endl;
00352 if (thePrim==0) thePrim=avertex->GetPrimary(trackID);
00353 }
00354
00355 double px=0.,py=0.,pz=0.;
00356 double pInit = 0.;
00357
00358 if (thePrim != 0) {
00359 px = thePrim->GetPx();
00360 py = thePrim->GetPy();
00361 pz = thePrim->GetPz();
00362 pInit = sqrt(pow(px,2.)+pow(py,2.)+pow(pz,2.));
00363 if (pInit==0) {
00364 std::cout << "ZdcTest End Of Event ERR: primary has p=0 " << std::endl;
00365 }
00366 } else {
00367 std::cout << "ZdcTest End Of Event ERR: could not find primary "
00368 << std::endl;
00369 }
00370
00371 }
00372
00373
00374 }
00375
00376 int iEvt = (*evt)()->GetEventID();
00377 if (iEvt < 10)
00378 std::cout << " ZdcTest Event " << iEvt << std::endl;
00379 else if ((iEvt < 100) && (iEvt%10 == 0))
00380 std::cout << " ZdcTest Event " << iEvt << std::endl;
00381 else if ((iEvt < 1000) && (iEvt%100 == 0))
00382 std::cout << " ZdcTest Event " << iEvt << std::endl;
00383 else if ((iEvt < 10000) && (iEvt%1000 == 0))
00384 std::cout << " ZdcTest Event " << iEvt << std::endl;
00385 }
00386
00387 void ZdcTestAnalysis::update(const EndOfRun * run) {;}
00388
00389 void ZdcTestAnalysis::finish(){
00390 if (doNTzdcstep) {
00391 zdcOutputStepFile->cd();
00392 zdcstepntuple->Write();
00393 std::cout << "ZdcTestAnalysis: Ntuple step written for event: "<<eventIndex<<std::endl;
00394 zdcOutputStepFile->Close();
00395 std::cout << "ZdcTestAnalysis: Step file closed" << std::endl;
00396 }
00397
00398
00399 if (doNTzdcevent) {
00400 zdcOutputEventFile->cd();
00401 zdceventntuple->Write("",TObject::kOverwrite);
00402 std::cout << "ZdcTestAnalysis: Ntuple event written for event: "<<eventIndex<<std::endl;
00403 zdcOutputEventFile->Close();
00404 std::cout << "ZdcTestAnalysis: Event file closed" << std::endl;
00405 }
00406
00407 }