Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:27:58

0001 #include "RecoTBCalo/ZDCTBAnalysis/interface/ZdcTBAnalysis.h"
0002 #include <sstream>
0003 #include <iostream>
0004 #include <vector>
0005 
0006 ZdcTBAnalysis::ZdcTBAnalysis() { ; }
0007 
0008 void ZdcTBAnalysis::setup(const std::string& outFileName) {
0009   TString outName = outFileName;
0010   outFile = new TFile(outName, "RECREATE");
0011   ZdcAnalize = new TTree("ZdcAnaTree", "ZdcAnaTree");
0012   ZdcAnalize->Branch("Trigger",
0013                      0,
0014                      "run/I:event/I:beamTrigger/I:fakeTrigger/I:"
0015                      "pedestalTrigger/I:outSpillPedestalTrigger/I:inSpillPedestalTrigger/I:"
0016                      "laserTrigger/I:laserTrigger/I:ledTrigger/I:spillTrigger/I");
0017   ZdcAnalize->Branch("TDC",
0018                      0,
0019                      "trigger/D:ttcL1/D:beamCoincidence/D:laserFlash/D:qiePhase/D:"
0020                      "TTOF1/D:TTOF2/D:m1[5]/D:m2[5]/D:m3[5]/D:"
0021                      "s1[5]/D:s2[5]/D:s3[5]/D:s4[5]/D:"
0022                      "bh1[5]/D:bh2[5]/D:bh3[5]/D:bh4[5]/D");
0023   ZdcAnalize->Branch("ADC",
0024                      0,
0025                      "VM/D:V3/D:V6/D:VH1/D:VH2/D:VH3/D:VH4/D:Ecal7x7/D:"
0026                      "Sci521/D:Sci528/D:CK1/D:CK2/D:CK3/D:SciVLE/D:S1/D:S2/D:S3/D:S4/D:"
0027                      "VMF/D:VMB/D:VM1/D:VM2/D:VM3/D:VM4/D:VM5/D:VM6/D:VM7/D:VM8/D:"
0028                      "TOF1/D:TOF2/D:BH1/D:BH2/D:BH3/BH4/D");
0029   ZdcAnalize->Branch("Chamb",
0030                      0,
0031                      "WCAx[5]/D:WCAy[5]/D:WCBx[5]/D:WCBy[5]/D:"
0032                      "WCCx[5]/D:WCCy[5]/D:WCDx[5]/D:WCDy[5]/D:WCEx[5]/D:WCEy[5]/D:"
0033                      "WCFx[5]/D:WCFy[5]/D:WCGx[5]/D:WCGy[5]/D:WCHx[5]/D:WCHy[5]/D");
0034   ZdcAnalize->Branch("ZDCP",
0035                      0,
0036                      "zdcHAD1/D:zdcHAD2/D:zdcHAD3/D:zdcHAD4/D:"
0037                      "zdcEM1/D:zdcEM2/D:zdcEM3/D:zdcEM4/D:zdcEM5/D:"
0038                      "zdcScint1/D:zdcScint2/D:"
0039                      "zdcExtras[7]/D");
0040   ZdcAnalize->Branch("ZDCN",
0041                      0,
0042                      "zdcHAD1/D:zdcHAD2/D:zdcHAD3/D:zdcHAD4/D:"
0043                      "zdcEM1/D:zdcEM2/D:zdcEM3/D:zdcEM4/D:zdcEM5/D:"
0044                      "zdcScint1/D:zdcScint2/D:"
0045                      "zdcExtras[7]/D");
0046   ZdcAnalize->GetBranch("Trigger")->SetAddress(&trigger);
0047   ZdcAnalize->GetBranch("TDC")->SetAddress(&tdc);
0048   ZdcAnalize->GetBranch("ADC")->SetAddress(&adc);
0049   ZdcAnalize->GetBranch("Chamb")->SetAddress(&chamb);
0050   ZdcAnalize->GetBranch("ZDCP")->SetAddress(&zdcp);
0051   ZdcAnalize->GetBranch("ZDCN")->SetAddress(&zdcn);
0052   ZdcAnalize->SetAutoSave();
0053 }
0054 
0055 void ZdcTBAnalysis::analyze(const HcalTBTriggerData& trg) {
0056   // trigger
0057   trigger.runNum = runNumber = trg.runNumber();
0058   trigger.eventNum = eventNumber = trg.eventNumber();
0059   isBeamTrigger = trg.wasBeamTrigger();
0060   isFakeTrigger = trg.wasFakeTrigger();
0061   isCalibTrigger = trg.wasSpillIgnorantPedestalTrigger();
0062   isOutSpillPedestalTrigger = trg.wasOutSpillPedestalTrigger();
0063   isInSpillPedestalTrigger = trg.wasInSpillPedestalTrigger();
0064   isLaserTrigger = trg.wasLaserTrigger();
0065   isLedTrigger = trg.wasLEDTrigger();
0066   isSpillTrigger = trg.wasInSpill();
0067 
0068   trigger.beamTrigger = trigger.fakeTrigger = trigger.calibTrigger = trigger.outSpillPedestalTrigger =
0069       trigger.inSpillPedestalTrigger = trigger.laserTrigger = trigger.ledTrigger = trigger.spillTrigger = 0;
0070 
0071   if (isBeamTrigger)
0072     trigger.beamTrigger = 1;
0073   if (isFakeTrigger)
0074     trigger.fakeTrigger = 1;
0075   if (isCalibTrigger)
0076     trigger.calibTrigger = 1;
0077   if (isOutSpillPedestalTrigger)
0078     trigger.outSpillPedestalTrigger = 1;
0079   if (isInSpillPedestalTrigger)
0080     trigger.inSpillPedestalTrigger = 1;
0081   if (isLaserTrigger)
0082     trigger.laserTrigger = 1;
0083   if (isLedTrigger)
0084     trigger.ledTrigger = 1;
0085   if (isSpillTrigger)
0086     trigger.spillTrigger = 1;
0087 }
0088 
0089 void ZdcTBAnalysis::analyze(const HcalTBTiming& times) {
0090   //times
0091   tdc.trigger = trigger_time = times.triggerTime();
0092   tdc.ttcL1 = ttc_L1a_time = times.ttcL1Atime();
0093   tdc.laserFlash = laser_flash = times.laserFlash();
0094   tdc.qiePhase = qie_phase = times.qiePhase();
0095   tdc.TOF1 = TOF1_time = times.TOF1Stime();
0096   tdc.TOF2 = TOF2_time = times.TOF2Stime();
0097 
0098   // just take 5 first hits of multihit tdc (5 tick cycles)
0099   int indx = 0;
0100   int indTop = 5;
0101   for (indx = 0; indx < times.BeamCoincidenceCount(); indx++)
0102     if (indx < indTop)
0103       tdc.beamCoincidence[indx] = beam_coincidence[indx] = times.BeamCoincidenceHits(indx);
0104   for (indx = 0; indx < times.M1Count(); indx++)
0105     if (indx < indTop)
0106       tdc.m1[indx] = m1hits[indx] = times.M1Hits(indx);
0107   for (indx = 0; indx < times.M2Count(); indx++)
0108     if (indx < indTop)
0109       tdc.m2[indx] = m2hits[indx] = times.M2Hits(indx);
0110   for (indx = 0; indx < times.M3Count(); indx++)
0111     if (indx < indTop)
0112       tdc.m3[indx] = m3hits[indx] = times.M3Hits(indx);
0113   for (indx = 0; indx < times.S1Count(); indx++)
0114     if (indx < indTop)
0115       tdc.s1[indx] = s1hits[indx] = times.S1Hits(indx);
0116   for (indx = 0; indx < times.S2Count(); indx++)
0117     if (indx < indTop)
0118       tdc.s2[indx] = s2hits[indx] = times.S2Hits(indx);
0119   for (indx = 0; indx < times.S3Count(); indx++)
0120     if (indx < indTop)
0121       tdc.s3[indx] = s3hits[indx] = times.S3Hits(indx);
0122   for (indx = 0; indx < times.S4Count(); indx++)
0123     if (indx < indTop)
0124       tdc.s4[indx] = s4hits[indx] = times.S4Hits(indx);
0125   for (indx = 0; indx < times.BH1Count(); indx++)
0126     if (indx < indTop)
0127       tdc.bh1[indx] = bh1hits[indx] = times.BH1Hits(indx);
0128   for (indx = 0; indx < times.BH2Count(); indx++)
0129     if (indx < indTop)
0130       tdc.bh2[indx] = bh2hits[indx] = times.BH2Hits(indx);
0131   for (indx = 0; indx < times.BH3Count(); indx++)
0132     if (indx < indTop)
0133       tdc.bh3[indx] = bh3hits[indx] = times.BH3Hits(indx);
0134   for (indx = 0; indx < times.BH4Count(); indx++)
0135     if (indx < indTop)
0136       tdc.bh4[indx] = bh4hits[indx] = times.BH4Hits(indx);
0137 }
0138 
0139 void ZdcTBAnalysis::analyze(const HcalTBBeamCounters& bc) {
0140   //beam counters
0141   adc.VM = VMadc = bc.VMadc();
0142   adc.V3 = V3adc = bc.V3adc();
0143   adc.V6 = V6adc = bc.V6adc();
0144   adc.VH1 = VH1adc = bc.VH1adc();
0145   adc.VH2 = VH2adc = bc.VH2adc();
0146   adc.VH3 = VH3adc = bc.VH3adc();
0147   adc.VH4 = VH4adc = bc.VH4adc();
0148   adc.Ecal7x7 = Ecal7x7adc = bc.Ecal7x7();
0149   adc.Sci521 = Sci521adc = bc.Sci521adc();
0150   adc.Sci528 = Sci528adc = bc.Sci528adc();
0151   adc.CK1 = CK1adc = bc.CK1adc();
0152   adc.CK2 = CK2adc = bc.CK2adc();
0153   adc.CK3 = CK3adc = bc.CK3adc();
0154   adc.SciVLE = SciVLEadc = bc.SciVLEadc();
0155   adc.S1 = S1adc = bc.S1adc();
0156   adc.S2 = S2adc = bc.S2adc();
0157   adc.S3 = S3adc = bc.S3adc();
0158   adc.S4 = S4adc = bc.S4adc();
0159   adc.VMF = VMFadc = bc.VMFadc();
0160   adc.VMB = VMBadc = bc.VMBadc();
0161   adc.VM1 = VM1adc = bc.VM1adc();
0162   adc.VM2 = VM2adc = bc.VM2adc();
0163   adc.VM3 = VM3adc = bc.VM3adc();
0164   adc.VM4 = VM4adc = bc.VM4adc();
0165   adc.VM5 = VM5adc = bc.VM5adc();
0166   adc.VM6 = VM6adc = bc.VM6adc();
0167   adc.VM7 = VM7adc = bc.VM7adc();
0168   adc.VM8 = VM8adc = bc.VM8adc();
0169   adc.TOF1 = TOF1adc = bc.TOF1Sadc();
0170   adc.TOF2 = TOF2adc = bc.TOF2Sadc();
0171   adc.BH1 = BH1adc = bc.BH1adc();
0172   adc.BH2 = BH2adc = bc.BH2adc();
0173   adc.BH3 = BH3adc = bc.BH3adc();
0174   adc.BH4 = BH4adc = bc.BH4adc();
0175 }
0176 
0177 void ZdcTBAnalysis::analyze(const HcalTBEventPosition& chpos) {
0178   //chambers position
0179   chpos.getChamberHits('A', wcax, wcay);
0180   chpos.getChamberHits('B', wcbx, wcby);
0181   chpos.getChamberHits('C', wccx, wccy);
0182   chpos.getChamberHits('D', wcdx, wcdy);
0183   chpos.getChamberHits('E', wcex, wcey);
0184   chpos.getChamberHits('F', wcfx, wcfy);
0185   chpos.getChamberHits('G', wcgx, wcgy);
0186   chpos.getChamberHits('H', wchx, wchy);
0187 
0188   // just take 5 first hits of chambers (5 tick cycles)
0189   unsigned int indTop = 5;
0190   unsigned int indx = 0;
0191   for (indx = 0; indx < wcax.size(); indx++)
0192     if (indx < indTop)
0193       chamb.WCAx[indx] = wcax[indx];
0194   for (indx = 0; indx < wcay.size(); indx++)
0195     if (indx < indTop)
0196       chamb.WCAy[indx] = wcay[indx];
0197   for (indx = 0; indx < wcbx.size(); indx++)
0198     if (indx < indTop)
0199       chamb.WCBx[indx] = wcbx[indx];
0200   for (indx = 0; indx < wcby.size(); indx++)
0201     if (indx < indTop)
0202       chamb.WCBy[indx] = wcby[indx];
0203   for (indx = 0; indx < wccx.size(); indx++)
0204     if (indx < indTop)
0205       chamb.WCCx[indx] = wccx[indx];
0206   for (indx = 0; indx < wccy.size(); indx++)
0207     if (indx < indTop)
0208       chamb.WCCy[indx] = wccy[indx];
0209   for (indx = 0; indx < wcdx.size(); indx++)
0210     if (indx < indTop)
0211       chamb.WCDx[indx] = wcdx[indx];
0212   for (indx = 0; indx < wcdy.size(); indx++)
0213     if (indx < indTop)
0214       chamb.WCDy[indx] = wcdy[indx];
0215   for (indx = 0; indx < wcdx.size(); indx++)
0216     if (indx < indTop)
0217       chamb.WCEx[indx] = wcex[indx];
0218   for (indx = 0; indx < wcey.size(); indx++)
0219     if (indx < indTop)
0220       chamb.WCEy[indx] = wcey[indx];
0221   for (indx = 0; indx < wcfx.size(); indx++)
0222     if (indx < indTop)
0223       chamb.WCFx[indx] = wcfx[indx];
0224   for (indx = 0; indx < wcfy.size(); indx++)
0225     if (indx < indTop)
0226       chamb.WCFy[indx] = wcfy[indx];
0227   for (indx = 0; indx < wcgx.size(); indx++)
0228     if (indx < indTop)
0229       chamb.WCGx[indx] = wcgx[indx];
0230   for (indx = 0; indx < wcgy.size(); indx++)
0231     if (indx < indTop)
0232       chamb.WCGy[indx] = wcgy[indx];
0233   for (indx = 0; indx < wchx.size(); indx++)
0234     if (indx < indTop)
0235       chamb.WCHx[indx] = wchx[indx];
0236   for (indx = 0; indx < wchy.size(); indx++)
0237     if (indx < indTop)
0238       chamb.WCHy[indx] = wchy[indx];
0239 }
0240 
0241 void ZdcTBAnalysis::analyze(const ZDCRecHitCollection& zdcHits) {
0242   // zdc hits
0243   std::cout << "****************************************************" << std::endl;
0244   ZDCRecHitCollection::const_iterator i;
0245   for (i = zdcHits.begin(); i != zdcHits.end(); i++) {
0246     energy = i->energy();
0247     detID = i->id();
0248     iside = detID.zside();
0249     isection = detID.section();
0250     ichannel = detID.channel();
0251     idepth = detID.depth();
0252     std::cout << "energy: " << energy << " detID: " << detID << " side: " << iside << " section: " << isection
0253               << " channel: " << ichannel << " depth: " << idepth << std::endl;
0254 
0255     if (iside > 0) {
0256       if (ichannel == 1 && isection == 1)
0257         zdcp.zdcEMMod1 = energy;
0258       if (ichannel == 2 && isection == 1)
0259         zdcp.zdcEMMod2 = energy;
0260       if (ichannel == 3 && isection == 1)
0261         zdcp.zdcEMMod3 = energy;
0262       if (ichannel == 4 && isection == 1)
0263         zdcp.zdcEMMod4 = energy;
0264       if (ichannel == 5 && isection == 1)
0265         zdcp.zdcEMMod5 = energy;
0266       if (ichannel == 1 && isection == 2)
0267         zdcp.zdcHADMod1 = energy;
0268       if (ichannel == 2 && isection == 2)
0269         zdcp.zdcHADMod2 = energy;
0270       if (ichannel == 3 && isection == 2)
0271         zdcp.zdcHADMod3 = energy;
0272       if (ichannel == 4 && isection == 2)
0273         zdcp.zdcHADMod4 = energy;
0274       if (ichannel == 1 && isection == 3)
0275         zdcp.zdcScint1 = energy;
0276     }
0277     if (iside < 0) {
0278       if (ichannel == 1 && isection == 1)
0279         zdcn.zdcEMMod1 = energy;
0280       if (ichannel == 2 && isection == 1)
0281         zdcn.zdcEMMod2 = energy;
0282       if (ichannel == 3 && isection == 1)
0283         zdcn.zdcEMMod3 = energy;
0284       if (ichannel == 4 && isection == 1)
0285         zdcn.zdcEMMod4 = energy;
0286       if (ichannel == 5 && isection == 1)
0287         zdcn.zdcEMMod5 = energy;
0288       if (ichannel == 1 && isection == 2)
0289         zdcn.zdcHADMod1 = energy;
0290       if (ichannel == 2 && isection == 2)
0291         zdcn.zdcHADMod2 = energy;
0292       if (ichannel == 3 && isection == 2)
0293         zdcn.zdcHADMod3 = energy;
0294       if (ichannel == 4 && isection == 2)
0295         zdcn.zdcHADMod4 = energy;
0296       if (ichannel == 1 && isection == 3)
0297         zdcn.zdcScint1 = energy;
0298     }
0299   }
0300 }
0301 
0302 void ZdcTBAnalysis::fillTree() { ZdcAnalize->Fill(); }
0303 
0304 void ZdcTBAnalysis::done() {
0305   ZdcAnalize->Print();
0306   outFile->cd();
0307   ZdcAnalize->Write();
0308   outFile->Close();
0309 }