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 // ****************************************************************
25 #include "AliLoader.h"
26 #include "AliRunLoader.h"
27 #include "AliTriggerInput.h"
30 #include "AliZDCDigit.h"
31 #include "AliZDCTrigger.h"
33 //________________________________________________________________
35 ClassImp(AliZDCTrigger)
37 //________________________________________________________________
38 AliZDCTrigger::AliZDCTrigger() :
47 fZDCLeftSemiCentrCut(0),
48 fZDCRightSemiCentrCut(0),
55 SetZDCLeftEMDCuts(0,0);
56 SetZDCRightEMDCuts(0,0);
60 //________________________________________________________________
61 void AliZDCTrigger::CreateInputs()
65 // Do not create inputs again!!
66 if( fInputs.GetEntriesFast() > 0 ) return;
68 fInputs.AddLast(new AliTriggerInput("ZDC_1_L1", "ZDC", 1));
69 fInputs.AddLast(new AliTriggerInput("ZDC_2_L1", "ZDC", 1));
70 fInputs.AddLast(new AliTriggerInput("ZDC_3_L1", "ZDC", 1));
71 fInputs.AddLast(new AliTriggerInput("ZDC_EMD_L1", "ZDC", 1));
74 //________________________________________________________________
75 void AliZDCTrigger::Trigger()
80 AliRunLoader *runLoader = AliRunLoader::Instance();
82 AliLoader *aZDCLoader = runLoader->GetLoader("ZDCLoader");
84 aZDCLoader->LoadDigits("READ");
86 AliZDCDigit* pdigit = &digit;
87 TTree* tD = aZDCLoader->TreeD();
89 cerr<<"AliZDCTrigger: digits tree not found\n";
92 tD->SetBranchAddress("ZDC", &pdigit);
94 Float_t signalZNLeft[]={0,0}, signalZPLeft[]={0,0}, signalZDCLeftSum[]={0,0};
95 Float_t signalZNRight[]={0,0}, signalZPRight[]={0,0}, signalZDCRightSum[]={0,0};
96 Float_t signalZEMSum[]={0,0};
97 for(Int_t iDigit=0; iDigit<tD->GetEntries(); iDigit++){
101 if(digit.GetSector(0)==1)
102 for(Int_t i=0; i<2; i++){ //0=high range; 1=low range
103 signalZNLeft[i] += digit.GetADCValue(i);
104 signalZDCLeftSum[i] += digit.GetADCValue(i);
106 else if(digit.GetSector(0)==2)
107 for(Int_t i=0; i<2; i++){
108 signalZPLeft[i] += digit.GetADCValue(i);
109 signalZDCLeftSum[i] += digit.GetADCValue(i);
111 else if(digit.GetSector(0)==3)
112 for(Int_t i=0; i<2; i++) signalZEMSum[i] += digit.GetADCValue(i);
114 else if(digit.GetSector(0)==4)
115 for(Int_t i=0; i<2; i++){ //0=high range; 1=low range
116 signalZNRight[i] += digit.GetADCValue(i);
117 signalZDCRightSum[i] += digit.GetADCValue(i);
119 else if(digit.GetSector(0)==5)
120 for(Int_t i=0; i<2; i++){
121 signalZPRight[i] += digit.GetADCValue(i);
122 signalZDCRightSum[i] += digit.GetADCValue(i);
125 // *******************************************************************
126 if(signalZDCLeftSum[1]>fZDCLeftMBCut && signalZDCRightSum[1]>fZDCRightMBCut)
127 // *** ZDC minimum bias trigger
128 SetInput("ZDC_1_L1");
129 // *******************************************************************
130 if(signalZDCLeftSum[1]>fZDCLeftCentrCut && signalZDCLeftSum[1]<fZDCLeftSemiCentrCut &&
131 signalZDCRightSum[1]>fZDCRightCentrCut && signalZDCRightSum[1]<fZDCRightSemiCentrCut
132 && signalZEMSum[1]>fZEMCentrCut)
133 // *** ZDC semi-central (10-40%)
134 SetInput("ZDC_2_L1");
135 // *******************************************************************
136 if(signalZDCLeftSum[1]>fZDCLeftMinCut && signalZDCLeftSum[1]<fZDCLeftCentrCut &&
137 signalZDCRightSum[1]>fZDCRightMinCut && signalZDCRightSum[1]<fZDCRightCentrCut &&
138 signalZEMSum[1]>fZEMCentrCut)
139 // *** ZDC central (0-10%)
140 SetInput("ZDC_3_L1");
141 // *******************************************************************
142 if(signalZNLeft[0]>fZDCLeftEMDCuts[0] && signalZNLeft[0]<fZDCLeftEMDCuts[1] &&
143 signalZNRight[0]>fZDCRightEMDCuts[0] && signalZNRight[0]<fZDCRightEMDCuts[1] &&
144 signalZEMSum[1]<fZEMMinCut){ // *** 1n EMD trigger
145 SetInput("ZDC_EMD_L1");
150 //________________________________________________________________
151 void AliZDCTrigger::SetZDCLeftMinCut(Float_t ZDCLeftMinCut)
153 // Set default cut values for ZDC trigger
155 if(ZDCLeftMinCut) fZDCLeftMinCut = ZDCLeftMinCut;
156 else fZDCLeftMinCut = 800.;
158 //________________________________________________________________
159 void AliZDCTrigger::SetZDCRightMinCut(Float_t ZDCRightMinCut)
161 // Set default cut values for ZDC trigger
163 if(ZDCRightMinCut) fZDCRightMinCut = ZDCRightMinCut;
164 else fZDCRightMinCut = 800.;
167 //________________________________________________________________
168 void AliZDCTrigger::SetZEMMinCut(Float_t ZEMMinCut)
170 // Set default cut values for ZDC trigger
172 if(ZEMMinCut) fZEMMinCut = ZEMMinCut;
173 else fZEMMinCut = 80.;
175 //________________________________________________________________
176 void AliZDCTrigger::SetZDCLeftEMDCuts(Float_t* ZDCLeftEMDCuts)
178 // Set default cut values for ZDC trigger
180 if(ZDCLeftEMDCuts) for(int j=0; j<2; j++) fZDCLeftEMDCuts[j] = ZDCLeftEMDCuts[j];
182 fZDCLeftEMDCuts[0] = 600.;
183 fZDCLeftEMDCuts[1] = 1000.;
186 //________________________________________________________________
187 void AliZDCTrigger::SetZDCLeftEMDCuts(Float_t ZDCLeftEMDCutInf,
188 Float_t ZDCLeftEMDCutSup)
190 // Set default cut values for ZDC trigger
192 if(ZDCLeftEMDCutInf && ZDCLeftEMDCutSup){
193 fZDCLeftEMDCuts[0]=ZDCLeftEMDCutInf;
194 fZDCLeftEMDCuts[1]=ZDCLeftEMDCutSup;
197 fZDCLeftEMDCuts[0] = 600.;
198 fZDCLeftEMDCuts[1] = 1000.;
201 //________________________________________________________________
202 void AliZDCTrigger::SetZDCRightEMDCuts(Float_t* ZDCRightEMDCuts)
204 // Set default cut values for ZDC trigger
206 if(ZDCRightEMDCuts) for(int j=0; j<2; j++) fZDCRightEMDCuts[j] = ZDCRightEMDCuts[j];
208 fZDCRightEMDCuts[0] = 600.;
209 fZDCRightEMDCuts[1] = 1000.;
212 //________________________________________________________________
213 void AliZDCTrigger::SetZDCRightEMDCuts(Float_t ZDCRightEMDCutInf,
214 Float_t ZDCRightEMDCutSup)
216 // Set default cut values for ZDC trigger
218 if(ZDCRightEMDCutInf && ZDCRightEMDCutSup){
219 fZDCRightEMDCuts[0]=ZDCRightEMDCutInf;
220 fZDCRightEMDCuts[1]=ZDCRightEMDCutSup;
223 fZDCRightEMDCuts[0] = 600.;
224 fZDCRightEMDCuts[1] = 1000.;
227 //________________________________________________________________
228 void AliZDCTrigger::SetZDCLeftMBCut(Float_t ZDCLeftMBCut)
230 // Set default cut values for ZDC trigger
232 if(ZDCLeftMBCut) fZDCLeftMBCut = ZDCLeftMBCut;
233 else fZDCLeftMBCut = 800.;
235 //________________________________________________________________
236 void AliZDCTrigger::SetZDCRightMBCut(Float_t ZDCRightMBCut)
238 // Set default cut values for ZDC trigger
240 if(ZDCRightMBCut) fZDCRightMBCut = ZDCRightMBCut;
241 else fZDCRightMBCut = 800.;
243 //________________________________________________________________
244 void AliZDCTrigger::SetZDCLeftCentrCut(Float_t ZDCLeftCentrCut)
246 // Set default cut values for ZDC trigger
248 if(ZDCLeftCentrCut) fZDCLeftCentrCut = ZDCLeftCentrCut;
249 else fZDCLeftCentrCut = 10000.;
251 //________________________________________________________________
252 void AliZDCTrigger::SetZDCRightCentrCut(Float_t ZDCRightCentrCut)
254 // Set default cut values for ZDC trigger
256 if(ZDCRightCentrCut) fZDCRightCentrCut = ZDCRightCentrCut;
257 else fZDCRightCentrCut = 10000.;
259 //________________________________________________________________
260 void AliZDCTrigger::SetZDCLeftSemiCentrCut(Float_t ZDCLeftSemiCentrCut)
262 // Set default cut values for ZDC trigger
264 if(ZDCLeftSemiCentrCut) fZDCLeftSemiCentrCut = ZDCLeftSemiCentrCut;
265 else fZDCLeftSemiCentrCut = 18500.;
267 //________________________________________________________________
268 void AliZDCTrigger::SetZDCRightSemiCentrCut(Float_t ZDCRightSemiCentrCut)
270 // Set default cut values for ZDC trigger
272 if(ZDCRightSemiCentrCut) fZDCRightSemiCentrCut = ZDCRightSemiCentrCut;
273 else fZDCRightSemiCentrCut = 18500.;
275 //________________________________________________________________
276 void AliZDCTrigger::SetZEMCentrCut(Float_t ZEMCentrCut)
278 // Set default cut values for ZDC trigger
280 if(ZEMCentrCut) fZEMCentrCut = ZEMCentrCut;
281 else fZEMCentrCut = 210.;