1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
16 // ****************************************************************
18 // Trigger class for ZDC
20 // ****************************************************************
24 #include "AliLoader.h"
25 #include "AliRunLoader.h"
26 #include "AliTriggerInput.h"
29 #include "AliZDCDigit.h"
30 #include "AliZDCTrigger.h"
32 //________________________________________________________________
33 ClassImp(AliZDCTrigger)
35 //________________________________________________________________
36 AliZDCTrigger::AliZDCTrigger() : AliTriggerDetector()
45 SetZDCLeftEMDCuts(0,0);
46 SetZDCRightEMDCuts(0,0);
49 SetZDCLeftCentrCut(0);
50 SetZDCRightCentrCut(0);
51 SetZDCLeftSemiCentrCut(0);
52 SetZDCRightSemiCentrCut(0);
57 //________________________________________________________________
58 void AliZDCTrigger::CreateInputs()
62 // Do not create inputs again!!
63 if( fInputs.GetEntriesFast() > 0 ) return;
65 fInputs.AddLast(new AliTriggerInput("ZDC_1_L1", "ZDC Minimum Bias", 0x01));
66 fInputs.AddLast(new AliTriggerInput("ZDC_2_L1", "ZDC Semi-central", 0x02));
67 fInputs.AddLast(new AliTriggerInput("ZDC_3_L1", "ZDC Central", 0x04));
68 fInputs.AddLast(new AliTriggerInput("ZDC_EMD_L1", "ZDC EMD events", 0x08));
71 //________________________________________________________________
72 void AliZDCTrigger::Trigger()
77 AliRunLoader *runLoader = gAlice->GetRunLoader();
79 AliLoader *aZDCLoader = runLoader->GetLoader("ZDCLoader");
80 aZDCLoader->LoadDigits("READ");
82 AliZDCDigit* pdigit = &digit;
83 TTree* tD = aZDCLoader->TreeD();
85 cerr<<"AliZDCTrigger: digits tree not found\n";
88 tD->SetBranchAddress("ZDC", &pdigit);
90 Float_t signalZNLeft[]={0,0}, signalZPLeft[]={0,0}, signalZDCLeftSum[]={0,0};
91 Float_t signalZNRight[]={0,0}, signalZPRight[]={0,0}, signalZDCRightSum[]={0,0};
92 Float_t signalZEMSum[]={0,0};
93 for(Int_t iDigit=0; iDigit<tD->GetEntries(); iDigit++){
97 if(digit.GetSector(0)==1)
98 for(Int_t i=0; i<2; i++){ //0=high range; 1=low range
99 signalZNLeft[i] += digit.GetADCValue(i);
100 signalZDCLeftSum[i] += digit.GetADCValue(i);
102 else if(digit.GetSector(0)==2)
103 for(Int_t i=0; i<2; i++){
104 signalZPLeft[i] += digit.GetADCValue(i);
105 signalZDCLeftSum[i] += digit.GetADCValue(i);
107 else if(digit.GetSector(0)==3)
108 for(Int_t i=0; i<2; i++) signalZEMSum[i] += digit.GetADCValue(i);
110 else if(digit.GetSector(0)==4)
111 for(Int_t i=0; i<2; i++){ //0=high range; 1=low range
112 signalZNRight[i] += digit.GetADCValue(i);
113 signalZDCRightSum[i] += digit.GetADCValue(i);
115 else if(digit.GetSector(0)==5)
116 for(Int_t i=0; i<2; i++){
117 signalZPRight[i] += digit.GetADCValue(i);
118 signalZDCRightSum[i] += digit.GetADCValue(i);
121 // *******************************************************************
122 if(signalZDCLeftSum[1]>fZDCLeftMBCut && signalZDCRightSum[1]>fZDCRightMBCut)
123 // *** ZDC minimum bias trigger
124 SetInput("ZDC_1_L1");
125 // *******************************************************************
126 if(signalZDCLeftSum[1]>fZDCLeftCentrCut && signalZDCLeftSum[1]<fZDCLeftSemiCentrCut &&
127 signalZDCRightSum[1]>fZDCRightCentrCut && signalZDCRightSum[1]<fZDCRightSemiCentrCut
128 && signalZEMSum[1]>fZEMCentrCut)
129 // *** ZDC semi-central (10-40%)
130 SetInput("ZDC_2_L1");
131 // *******************************************************************
132 if(signalZDCLeftSum[1]>fZDCLeftMinCut && signalZDCLeftSum[1]<fZDCLeftCentrCut &&
133 signalZDCRightSum[1]>fZDCRightMinCut && signalZDCRightSum[1]<fZDCRightCentrCut &&
134 signalZEMSum[1]>fZEMCentrCut)
135 // *** ZDC central (0-10%)
136 SetInput("ZDC_3_L1");
137 // *******************************************************************
138 if(signalZNLeft[0]>fZDCLeftEMDCuts[0] && signalZNLeft[0]<fZDCLeftEMDCuts[1] &&
139 signalZNRight[0]>fZDCRightEMDCuts[0] && signalZNRight[0]<fZDCRightEMDCuts[1] &&
140 signalZEMSum[1]<fZEMMinCut){ // *** 1n EMD trigger
141 SetInput("ZDC_EMD_L1");
146 //________________________________________________________________
147 void AliZDCTrigger::SetZDCLeftMinCut(Float_t ZDCLeftMinCut)
149 // Set default cut values for ZDC trigger
151 if(ZDCLeftMinCut) fZDCLeftMinCut = ZDCLeftMinCut;
152 else fZDCLeftMinCut = 800.;
154 //________________________________________________________________
155 void AliZDCTrigger::SetZDCRightMinCut(Float_t ZDCRightMinCut)
157 // Set default cut values for ZDC trigger
159 if(ZDCRightMinCut) fZDCRightMinCut = ZDCRightMinCut;
160 else fZDCRightMinCut = 800.;
163 //________________________________________________________________
164 void AliZDCTrigger::SetZEMMinCut(Float_t ZEMMinCut)
166 // Set default cut values for ZDC trigger
168 if(ZEMMinCut) fZEMMinCut = ZEMMinCut;
169 else fZEMMinCut = 80.;
171 //________________________________________________________________
172 void AliZDCTrigger::SetZDCLeftEMDCuts(Float_t* ZDCLeftEMDCuts)
174 // Set default cut values for ZDC trigger
176 if(ZDCLeftEMDCuts) for(int j=0; j<2; j++) fZDCLeftEMDCuts[j] = ZDCLeftEMDCuts[j];
178 fZDCLeftEMDCuts[0] = 600.;
179 fZDCLeftEMDCuts[1] = 1000.;
182 //________________________________________________________________
183 void AliZDCTrigger::SetZDCLeftEMDCuts(Float_t ZDCLeftEMDCutInf,
184 Float_t ZDCLeftEMDCutSup)
186 // Set default cut values for ZDC trigger
188 if(ZDCLeftEMDCutInf && ZDCLeftEMDCutSup){
189 fZDCLeftEMDCuts[0]=ZDCLeftEMDCutInf;
190 fZDCLeftEMDCuts[1]=ZDCLeftEMDCutSup;
193 fZDCLeftEMDCuts[0] = 600.;
194 fZDCLeftEMDCuts[1] = 1000.;
197 //________________________________________________________________
198 void AliZDCTrigger::SetZDCRightEMDCuts(Float_t* ZDCRightEMDCuts)
200 // Set default cut values for ZDC trigger
202 if(ZDCRightEMDCuts) for(int j=0; j<2; j++) fZDCRightEMDCuts[j] = ZDCRightEMDCuts[j];
204 fZDCRightEMDCuts[0] = 600.;
205 fZDCRightEMDCuts[1] = 1000.;
208 //________________________________________________________________
209 void AliZDCTrigger::SetZDCRightEMDCuts(Float_t ZDCRightEMDCutInf,
210 Float_t ZDCRightEMDCutSup)
212 // Set default cut values for ZDC trigger
214 if(ZDCRightEMDCutInf && ZDCRightEMDCutSup){
215 fZDCRightEMDCuts[0]=ZDCRightEMDCutInf;
216 fZDCRightEMDCuts[1]=ZDCRightEMDCutSup;
219 fZDCRightEMDCuts[0] = 600.;
220 fZDCRightEMDCuts[1] = 1000.;
223 //________________________________________________________________
224 void AliZDCTrigger::SetZDCLeftMBCut(Float_t ZDCLeftMBCut)
226 // Set default cut values for ZDC trigger
228 if(ZDCLeftMBCut) fZDCLeftMBCut = ZDCLeftMBCut;
229 else fZDCLeftMBCut = 800.;
231 //________________________________________________________________
232 void AliZDCTrigger::SetZDCRightMBCut(Float_t ZDCRightMBCut)
234 // Set default cut values for ZDC trigger
236 if(ZDCRightMBCut) fZDCRightMBCut = ZDCRightMBCut;
237 else fZDCRightMBCut = 800.;
239 //________________________________________________________________
240 void AliZDCTrigger::SetZDCLeftCentrCut(Float_t ZDCLeftCentrCut)
242 // Set default cut values for ZDC trigger
244 if(ZDCLeftCentrCut) fZDCLeftCentrCut = ZDCLeftCentrCut;
245 else fZDCLeftCentrCut = 10000.;
247 //________________________________________________________________
248 void AliZDCTrigger::SetZDCRightCentrCut(Float_t ZDCRightCentrCut)
250 // Set default cut values for ZDC trigger
252 if(ZDCRightCentrCut) fZDCRightCentrCut = ZDCRightCentrCut;
253 else fZDCRightCentrCut = 10000.;
255 //________________________________________________________________
256 void AliZDCTrigger::SetZDCLeftSemiCentrCut(Float_t ZDCLeftSemiCentrCut)
258 // Set default cut values for ZDC trigger
260 if(ZDCLeftSemiCentrCut) fZDCLeftSemiCentrCut = ZDCLeftSemiCentrCut;
261 else fZDCLeftSemiCentrCut = 18500.;
263 //________________________________________________________________
264 void AliZDCTrigger::SetZDCRightSemiCentrCut(Float_t ZDCRightSemiCentrCut)
266 // Set default cut values for ZDC trigger
268 if(ZDCRightSemiCentrCut) fZDCRightSemiCentrCut = ZDCRightSemiCentrCut;
269 else fZDCRightSemiCentrCut = 18500.;
271 //________________________________________________________________
272 void AliZDCTrigger::SetZEMCentrCut(Float_t ZEMCentrCut)
274 // Set default cut values for ZDC trigger
276 if(ZEMCentrCut) fZEMCentrCut = ZEMCentrCut;
277 else fZEMCentrCut = 210.;