Go to the documentation of this file.00001 #include "RecoTBCalo/ZDCTBAnalysis/interface/ZdcTBAnalysis.h"
00002 #include <sstream>
00003 #include <iostream>
00004 #include <vector>
00005
00006 ZdcTBAnalysis::ZdcTBAnalysis(){;}
00007
00008 void ZdcTBAnalysis::setup(const std::string& outFileName) {
00009 TString outName = outFileName;
00010 outFile = new TFile(outName,"RECREATE");
00011 ZdcAnalize = new TTree("ZdcAnaTree","ZdcAnaTree");
00012 ZdcAnalize->Branch("Trigger",0,"run/I:event/I:beamTrigger/I:fakeTrigger/I:"
00013 "pedestalTrigger/I:outSpillPedestalTrigger/I:inSpillPedestalTrigger/I:"
00014 "laserTrigger/I:laserTrigger/I:ledTrigger/I:spillTrigger/I");
00015 ZdcAnalize->Branch("TDC",0,"trigger/D:ttcL1/D:beamCoincidence/D:laserFlash/D:qiePhase/D:"
00016 "TTOF1/D:TTOF2/D:m1[5]/D:m2[5]/D:m3[5]/D:"
00017 "s1[5]/D:s2[5]/D:s3[5]/D:s4[5]/D:"
00018 "bh1[5]/D:bh2[5]/D:bh3[5]/D:bh4[5]/D");
00019 ZdcAnalize->Branch("ADC",0,"VM/D:V3/D:V6/D:VH1/D:VH2/D:VH3/D:VH4/D:Ecal7x7/D:"
00020 "Sci521/D:Sci528/D:CK1/D:CK2/D:CK3/D:SciVLE/D:S1/D:S2/D:S3/D:S4/D:"
00021 "VMF/D:VMB/D:VM1/D:VM2/D:VM3/D:VM4/D:VM5/D:VM6/D:VM7/D:VM8/D:"
00022 "TOF1/D:TOF2/D:BH1/D:BH2/D:BH3/BH4/D");
00023 ZdcAnalize->Branch("Chamb",0,"WCAx[5]/D:WCAy[5]/D:WCBx[5]/D:WCBy[5]/D:"
00024 "WCCx[5]/D:WCCy[5]/D:WCDx[5]/D:WCDy[5]/D:WCEx[5]/D:WCEy[5]/D:"
00025 "WCFx[5]/D:WCFy[5]/D:WCGx[5]/D:WCGy[5]/D:WCHx[5]/D:WCHy[5]/D");
00026 ZdcAnalize->Branch("ZDCP",0,"zdcHAD1/D:zdcHAD2/D:zdcHAD3/D:zdcHAD4/D:"
00027 "zdcEM1/D:zdcEM2/D:zdcEM3/D:zdcEM4/D:zdcEM5/D:"
00028 "zdcScint1/D:zdcScint2/D:"
00029 "zdcExtras[7]/D");
00030 ZdcAnalize->Branch("ZDCN",0,"zdcHAD1/D:zdcHAD2/D:zdcHAD3/D:zdcHAD4/D:"
00031 "zdcEM1/D:zdcEM2/D:zdcEM3/D:zdcEM4/D:zdcEM5/D:"
00032 "zdcScint1/D:zdcScint2/D:"
00033 "zdcExtras[7]/D");
00034 ZdcAnalize->GetBranch("Trigger")->SetAddress(&trigger);
00035 ZdcAnalize->GetBranch("TDC")->SetAddress(&tdc);
00036 ZdcAnalize->GetBranch("ADC")->SetAddress(&adc);
00037 ZdcAnalize->GetBranch("Chamb")->SetAddress(&chamb);
00038 ZdcAnalize->GetBranch("ZDCP")->SetAddress(&zdcp);
00039 ZdcAnalize->GetBranch("ZDCN")->SetAddress(&zdcn);
00040 ZdcAnalize->SetAutoSave();
00041 }
00042
00043
00044 void ZdcTBAnalysis::analyze(const HcalTBTriggerData& trg){
00045
00046 trigger.runNum = runNumber = trg.runNumber();
00047 trigger.eventNum = eventNumber = trg.eventNumber();
00048 isBeamTrigger = trg.wasBeamTrigger();
00049 isFakeTrigger = trg.wasFakeTrigger();
00050 isCalibTrigger = trg.wasSpillIgnorantPedestalTrigger();
00051 isOutSpillPedestalTrigger = trg.wasOutSpillPedestalTrigger();
00052 isInSpillPedestalTrigger = trg.wasInSpillPedestalTrigger();
00053 isLaserTrigger = trg.wasLaserTrigger();
00054 isLedTrigger = trg.wasLEDTrigger();
00055 isSpillTrigger = trg.wasInSpill();
00056
00057 trigger.beamTrigger = trigger.fakeTrigger = trigger.calibTrigger =
00058 trigger.outSpillPedestalTrigger = trigger.inSpillPedestalTrigger =
00059 trigger.laserTrigger = trigger.ledTrigger = trigger.spillTrigger = 0;
00060
00061 if(isBeamTrigger)trigger.beamTrigger = 1;
00062 if(isFakeTrigger)trigger.fakeTrigger = 1;
00063 if(isCalibTrigger)trigger.calibTrigger = 1;
00064 if(isOutSpillPedestalTrigger)trigger.outSpillPedestalTrigger = 1;
00065 if(isInSpillPedestalTrigger)trigger.inSpillPedestalTrigger = 1;
00066 if(isLaserTrigger)trigger.laserTrigger = 1;
00067 if(isLedTrigger)trigger.ledTrigger = 1;
00068 if(isSpillTrigger)trigger.spillTrigger = 1;
00069 }
00070
00071
00072 void ZdcTBAnalysis::analyze(const HcalTBTiming& times){
00073
00074 tdc.trigger =trigger_time = times.triggerTime();
00075 tdc.ttcL1 = ttc_L1a_time = times.ttcL1Atime();
00076 tdc.laserFlash = laser_flash = times.laserFlash();
00077 tdc.qiePhase = qie_phase = times.qiePhase();
00078 tdc.TOF1 = TOF1_time = times.TOF1Stime();
00079 tdc.TOF2 = TOF2_time = times.TOF2Stime();
00080
00081
00082 int indx = 0;
00083 int indTop = 5;
00084 for (indx=0; indx<times.BeamCoincidenceCount(); indx++)
00085 if (indx < indTop)tdc.beamCoincidence[indx] = beam_coincidence[indx] = times.BeamCoincidenceHits(indx);
00086 for (indx=0; indx<times.M1Count(); indx++)
00087 if (indx < indTop)tdc.m1[indx] = m1hits[indx] = times.M1Hits(indx);
00088 for (indx=0; indx<times.M2Count(); indx++)
00089 if (indx < indTop) tdc.m2[indx] = m2hits[indx] = times.M2Hits(indx);
00090 for (indx=0; indx<times.M3Count(); indx++)
00091 if (indx < indTop) tdc.m3[indx] = m3hits[indx] = times.M3Hits(indx);
00092 for (indx=0; indx<times.S1Count(); indx++)
00093 if (indx < indTop)tdc.s1[indx] = s1hits[indx]= times.S1Hits(indx);
00094 for (indx=0; indx<times.S2Count(); indx++)
00095 if (indx < indTop) tdc.s2[indx] = s2hits[indx] = times.S2Hits(indx);
00096 for (indx=0; indx<times.S3Count() ; indx++)
00097 if (indx < indTop) tdc.s3[indx] = s3hits[indx] = times.S3Hits(indx);
00098 for (indx=0; indx<times.S4Count(); indx++)
00099 if (indx < indTop) tdc.s4[indx] = s4hits[indx] = times.S4Hits(indx);
00100 for (indx=0; indx<times.BH1Count(); indx++)
00101 if (indx < indTop) tdc.bh1[indx] = bh1hits[indx] = times.BH1Hits(indx);
00102 for (indx=0; indx<times.BH2Count(); indx++)
00103 if (indx < indTop) tdc.bh2[indx] = bh2hits[indx] = times.BH2Hits(indx);
00104 for (indx=0; indx<times.BH3Count() ; indx++)
00105 if (indx < indTop) tdc.bh3[indx] = bh3hits[indx] = times.BH3Hits(indx);
00106 for (indx=0; indx<times.BH4Count(); indx++)
00107 if (indx < indTop) tdc.bh4[indx] = bh4hits[indx] = times.BH4Hits(indx);
00108 }
00109
00110 void ZdcTBAnalysis::analyze(const HcalTBBeamCounters& bc){
00111
00112 adc.VM = VMadc = bc.VMadc();
00113 adc.V3 = V3adc = bc.V3adc();
00114 adc.V6 = V6adc = bc.V6adc();
00115 adc.VH1 = VH1adc = bc.VH1adc();
00116 adc.VH2 = VH2adc = bc.VH2adc();
00117 adc.VH3 = VH3adc = bc.VH3adc();
00118 adc.VH4 = VH4adc = bc.VH4adc();
00119 adc.Ecal7x7 = Ecal7x7adc = bc.Ecal7x7();
00120 adc.Sci521 = Sci521adc = bc.Sci521adc();
00121 adc.Sci528 = Sci528adc = bc.Sci528adc();
00122 adc.CK1 = CK1adc = bc.CK1adc();
00123 adc.CK2 = CK2adc = bc.CK2adc();
00124 adc.CK3 = CK3adc = bc.CK3adc();
00125 adc.SciVLE = SciVLEadc = bc.SciVLEadc();
00126 adc.S1 = S1adc = bc.S1adc();
00127 adc.S2 = S2adc = bc.S2adc();
00128 adc.S3 = S3adc = bc.S3adc();
00129 adc.S4 = S4adc = bc.S4adc();
00130 adc.VMF = VMFadc = bc.VMFadc();
00131 adc.VMB = VMBadc = bc.VMBadc();
00132 adc.VM1 = VM1adc = bc.VM1adc();
00133 adc.VM2 = VM2adc = bc.VM2adc();
00134 adc.VM3 = VM3adc = bc.VM3adc();
00135 adc.VM4 = VM4adc = bc.VM4adc();
00136 adc.VM5 = VM5adc = bc.VM5adc();
00137 adc.VM6 = VM6adc = bc.VM6adc();
00138 adc.VM7 = VM7adc = bc.VM7adc();
00139 adc.VM8 = VM8adc = bc.VM8adc();
00140 adc.TOF1 = TOF1adc = bc.TOF1Sadc();
00141 adc.TOF2 = TOF2adc = bc.TOF2Sadc();
00142 adc.BH1 = BH1adc = bc.BH1adc();
00143 adc.BH2 = BH2adc = bc.BH2adc();
00144 adc.BH3 = BH3adc = bc.BH3adc();
00145 adc.BH4 = BH4adc = bc.BH4adc();
00146 }
00147
00148 void ZdcTBAnalysis::analyze(const HcalTBEventPosition& chpos){
00149
00150 chpos.getChamberHits('A',wcax,wcay);
00151 chpos.getChamberHits('B',wcbx,wcby);
00152 chpos.getChamberHits('C',wccx,wccy);
00153 chpos.getChamberHits('D',wcdx,wcdy);
00154 chpos.getChamberHits('E',wcex,wcey);
00155 chpos.getChamberHits('F',wcfx,wcfy);
00156 chpos.getChamberHits('G',wcgx,wcgy);
00157 chpos.getChamberHits('H',wchx,wchy);
00158
00159
00160 unsigned int indTop = 5;
00161 unsigned int indx = 0;
00162 for (indx = 0; indx < wcax.size(); indx++)
00163 if (indx < indTop)chamb.WCAx[indx] = wcax[indx];
00164 for (indx = 0; indx < wcay.size(); indx++)
00165 if (indx < indTop)chamb.WCAy[indx] = wcay[indx];
00166 for (indx = 0; indx < wcbx.size(); indx++)
00167 if (indx < indTop)chamb.WCBx[indx] = wcbx[indx];
00168 for (indx = 0; indx < wcby.size(); indx++)
00169 if (indx < indTop)chamb.WCBy[indx] = wcby[indx];
00170 for (indx = 0; indx < wccx.size(); indx++)
00171 if (indx < indTop)chamb.WCCx[indx] = wccx[indx];
00172 for (indx = 0; indx < wccy.size(); indx++)
00173 if (indx < indTop)chamb.WCCy[indx] = wccy[indx];
00174 for (indx = 0; indx < wcdx.size(); indx++)
00175 if (indx < indTop)chamb.WCDx[indx] = wcdx[indx];
00176 for (indx = 0; indx < wcdy.size(); indx++)
00177 if (indx < indTop)chamb.WCDy[indx] = wcdy[indx];
00178 for (indx = 0; indx < wcdx.size(); indx++)
00179 if (indx < indTop)chamb.WCEx[indx] = wcex[indx];
00180 for (indx = 0; indx < wcey.size(); indx++)
00181 if (indx < indTop)chamb.WCEy[indx] = wcey[indx];
00182 for (indx = 0; indx < wcfx.size(); indx++)
00183 if (indx < indTop)chamb.WCFx[indx] = wcfx[indx];
00184 for (indx = 0; indx < wcfy.size(); indx++)
00185 if (indx < indTop)chamb.WCFy[indx] = wcfy[indx];
00186 for (indx = 0; indx < wcgx.size(); indx++)
00187 if (indx < indTop)chamb.WCGx[indx] = wcgx[indx];
00188 for (indx = 0; indx < wcgy.size(); indx++)
00189 if (indx < indTop)chamb.WCGy[indx] = wcgy[indx];
00190 for (indx = 0; indx < wchx.size(); indx++)
00191 if (indx < indTop)chamb.WCHx[indx] = wchx[indx];
00192 for (indx = 0; indx < wchy.size(); indx++)
00193 if (indx < indTop)chamb.WCHy[indx] = wchy[indx];
00194 }
00195
00196 void ZdcTBAnalysis::analyze(const ZDCRecHitCollection& zdcHits){
00197
00198 std::cout<<"****************************************************"<<std::endl;
00199 ZDCRecHitCollection::const_iterator i;
00200 for(i = zdcHits.begin(); i!=zdcHits.end(); i++){
00201 energy = i->energy();
00202 detID = i->id();
00203 iside= detID.zside();
00204 isection = detID.section();
00205 ichannel = detID.channel();
00206 idepth = detID.depth();
00207 std::cout<<"energy: "<<energy<<" detID: "<<detID
00208 <<" side: "<<iside<<" section: "<<isection
00209 <<" channel: "<<ichannel<< " depth: "<<idepth<<std::endl;
00210
00211 if(iside>0){
00212 if(ichannel ==1 && isection ==1)zdcp.zdcEMMod1 = energy;
00213 if(ichannel ==2 && isection ==1)zdcp.zdcEMMod2 = energy;
00214 if(ichannel ==3 && isection ==1)zdcp.zdcEMMod3 = energy;
00215 if(ichannel ==4 && isection ==1)zdcp.zdcEMMod4 = energy;
00216 if(ichannel ==5 && isection ==1)zdcp.zdcEMMod5 = energy;
00217 if(ichannel ==1 && isection ==2)zdcp.zdcHADMod1 = energy;
00218 if(ichannel ==2 && isection ==2)zdcp.zdcHADMod2 = energy;
00219 if(ichannel ==3 && isection ==2)zdcp.zdcHADMod3 = energy;
00220 if(ichannel ==4 && isection ==2)zdcp.zdcHADMod4 = energy;
00221 if(ichannel ==1 && isection ==3)zdcp.zdcScint1 = energy;
00222 }
00223 if(iside<0){
00224 if(ichannel ==1 && isection ==1)zdcn.zdcEMMod1 = energy;
00225 if(ichannel ==2 && isection ==1)zdcn.zdcEMMod2 = energy;
00226 if(ichannel ==3 && isection ==1)zdcn.zdcEMMod3 = energy;
00227 if(ichannel ==4 && isection ==1)zdcn.zdcEMMod4 = energy;
00228 if(ichannel ==5 && isection ==1)zdcn.zdcEMMod5 = energy;
00229 if(ichannel ==1 && isection ==2)zdcn.zdcHADMod1 = energy;
00230 if(ichannel ==2 && isection ==2)zdcn.zdcHADMod2 = energy;
00231 if(ichannel ==3 && isection ==2)zdcn.zdcHADMod3 = energy;
00232 if(ichannel ==4 && isection ==2)zdcn.zdcHADMod4 = energy;
00233 if(ichannel ==1 && isection ==3)zdcn.zdcScint1 = energy;
00234 }
00235 }
00236 }
00237
00238 void ZdcTBAnalysis::fillTree(){
00239 ZdcAnalize->Fill();
00240 }
00241
00242 void ZdcTBAnalysis::done(){
00243 ZdcAnalize->Print();
00244 outFile->cd();
00245 ZdcAnalize->Write();
00246 outFile->Close();
00247 }
00248
00249