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() :
45 fZDCLeftSemiCentrCut(0),
46 fZDCRightSemiCentrCut(0),
53 SetZDCLeftEMDCuts(0,0);
54 SetZDCRightEMDCuts(0,0);
58 //________________________________________________________________
59 void AliZDCTrigger::CreateInputs()
63 // Do not create inputs again!!
64 if( fInputs.GetEntriesFast() > 0 ) return;
66 fInputs.AddLast(new AliTriggerInput("ZDC_1_L1", "ZDC", 1));
67 fInputs.AddLast(new AliTriggerInput("ZDC_2_L1", "ZDC", 1));
68 fInputs.AddLast(new AliTriggerInput("ZDC_3_L1", "ZDC", 1));
69 fInputs.AddLast(new AliTriggerInput("ZDC_EMD_L1", "ZDC", 1));
72 //________________________________________________________________
73 void AliZDCTrigger::Trigger()
78 AliRunLoader *runLoader = AliRunLoader::Instance();
80 AliLoader *aZDCLoader = runLoader->GetLoader("ZDCLoader");
82 aZDCLoader->LoadDigits("READ");
84 AliZDCDigit* pdigit = &digit;
85 TTree* tD = aZDCLoader->TreeD();
87 cerr<<"AliZDCTrigger: digits tree not found\n";
90 tD->SetBranchAddress("ZDC", &pdigit);
92 Float_t signalZNLeft[]={0,0}, signalZPLeft[]={0,0}, signalZDCLeftSum[]={0,0};
93 Float_t signalZNRight[]={0,0}, signalZPRight[]={0,0}, signalZDCRightSum[]={0,0};
94 Float_t signalZEMSum[]={0,0};
95 for(Int_t iDigit=0; iDigit<tD->GetEntries(); iDigit++){
99 if(digit.GetSector(0)==1)
100 for(Int_t i=0; i<2; i++){ //0=high range; 1=low range
101 signalZNLeft[i] += digit.GetADCValue(i);
102 signalZDCLeftSum[i] += digit.GetADCValue(i);
104 else if(digit.GetSector(0)==2)
105 for(Int_t i=0; i<2; i++){
106 signalZPLeft[i] += digit.GetADCValue(i);
107 signalZDCLeftSum[i] += digit.GetADCValue(i);
109 else if(digit.GetSector(0)==3)
110 for(Int_t i=0; i<2; i++) signalZEMSum[i] += digit.GetADCValue(i);
112 else if(digit.GetSector(0)==4)
113 for(Int_t i=0; i<2; i++){ //0=high range; 1=low range
114 signalZNRight[i] += digit.GetADCValue(i);
115 signalZDCRightSum[i] += digit.GetADCValue(i);
117 else if(digit.GetSector(0)==5)
118 for(Int_t i=0; i<2; i++){
119 signalZPRight[i] += digit.GetADCValue(i);
120 signalZDCRightSum[i] += digit.GetADCValue(i);
123 // *******************************************************************
124 if(signalZDCLeftSum[1]>fZDCLeftMBCut && signalZDCRightSum[1]>fZDCRightMBCut)
125 // *** ZDC minimum bias trigger
126 SetInput("ZDC_1_L1");
127 // *******************************************************************
128 if(signalZDCLeftSum[1]>fZDCLeftCentrCut && signalZDCLeftSum[1]<fZDCLeftSemiCentrCut &&
129 signalZDCRightSum[1]>fZDCRightCentrCut && signalZDCRightSum[1]<fZDCRightSemiCentrCut
130 && signalZEMSum[1]>fZEMCentrCut)
131 // *** ZDC semi-central (10-40%)
132 SetInput("ZDC_2_L1");
133 // *******************************************************************
134 if(signalZDCLeftSum[1]>fZDCLeftMinCut && signalZDCLeftSum[1]<fZDCLeftCentrCut &&
135 signalZDCRightSum[1]>fZDCRightMinCut && signalZDCRightSum[1]<fZDCRightCentrCut &&
136 signalZEMSum[1]>fZEMCentrCut)
137 // *** ZDC central (0-10%)
138 SetInput("ZDC_3_L1");
139 // *******************************************************************
140 if(signalZNLeft[0]>fZDCLeftEMDCuts[0] && signalZNLeft[0]<fZDCLeftEMDCuts[1] &&
141 signalZNRight[0]>fZDCRightEMDCuts[0] && signalZNRight[0]<fZDCRightEMDCuts[1] &&
142 signalZEMSum[1]<fZEMMinCut){ // *** 1n EMD trigger
143 SetInput("ZDC_EMD_L1");
148 //________________________________________________________________
149 void AliZDCTrigger::SetZDCLeftMinCut(Float_t ZDCLeftMinCut)
151 // Set default cut values for ZDC trigger
153 if(ZDCLeftMinCut) fZDCLeftMinCut = ZDCLeftMinCut;
154 else fZDCLeftMinCut = 800.;
156 //________________________________________________________________
157 void AliZDCTrigger::SetZDCRightMinCut(Float_t ZDCRightMinCut)
159 // Set default cut values for ZDC trigger
161 if(ZDCRightMinCut) fZDCRightMinCut = ZDCRightMinCut;
162 else fZDCRightMinCut = 800.;
165 //________________________________________________________________
166 void AliZDCTrigger::SetZEMMinCut(Float_t ZEMMinCut)
168 // Set default cut values for ZDC trigger
170 if(ZEMMinCut) fZEMMinCut = ZEMMinCut;
171 else fZEMMinCut = 80.;
173 //________________________________________________________________
174 void AliZDCTrigger::SetZDCLeftEMDCuts(Float_t* ZDCLeftEMDCuts)
176 // Set default cut values for ZDC trigger
178 if(ZDCLeftEMDCuts) for(int j=0; j<2; j++) fZDCLeftEMDCuts[j] = ZDCLeftEMDCuts[j];
180 fZDCLeftEMDCuts[0] = 600.;
181 fZDCLeftEMDCuts[1] = 1000.;
184 //________________________________________________________________
185 void AliZDCTrigger::SetZDCLeftEMDCuts(Float_t ZDCLeftEMDCutInf,
186 Float_t ZDCLeftEMDCutSup)
188 // Set default cut values for ZDC trigger
190 if(ZDCLeftEMDCutInf && ZDCLeftEMDCutSup){
191 fZDCLeftEMDCuts[0]=ZDCLeftEMDCutInf;
192 fZDCLeftEMDCuts[1]=ZDCLeftEMDCutSup;
195 fZDCLeftEMDCuts[0] = 600.;
196 fZDCLeftEMDCuts[1] = 1000.;
199 //________________________________________________________________
200 void AliZDCTrigger::SetZDCRightEMDCuts(Float_t* ZDCRightEMDCuts)
202 // Set default cut values for ZDC trigger
204 if(ZDCRightEMDCuts) for(int j=0; j<2; j++) fZDCRightEMDCuts[j] = ZDCRightEMDCuts[j];
206 fZDCRightEMDCuts[0] = 600.;
207 fZDCRightEMDCuts[1] = 1000.;
210 //________________________________________________________________
211 void AliZDCTrigger::SetZDCRightEMDCuts(Float_t ZDCRightEMDCutInf,
212 Float_t ZDCRightEMDCutSup)
214 // Set default cut values for ZDC trigger
216 if(ZDCRightEMDCutInf && ZDCRightEMDCutSup){
217 fZDCRightEMDCuts[0]=ZDCRightEMDCutInf;
218 fZDCRightEMDCuts[1]=ZDCRightEMDCutSup;
221 fZDCRightEMDCuts[0] = 600.;
222 fZDCRightEMDCuts[1] = 1000.;
225 //________________________________________________________________
226 void AliZDCTrigger::SetZDCLeftMBCut(Float_t ZDCLeftMBCut)
228 // Set default cut values for ZDC trigger
230 if(ZDCLeftMBCut) fZDCLeftMBCut = ZDCLeftMBCut;
231 else fZDCLeftMBCut = 800.;
233 //________________________________________________________________
234 void AliZDCTrigger::SetZDCRightMBCut(Float_t ZDCRightMBCut)
236 // Set default cut values for ZDC trigger
238 if(ZDCRightMBCut) fZDCRightMBCut = ZDCRightMBCut;
239 else fZDCRightMBCut = 800.;
241 //________________________________________________________________
242 void AliZDCTrigger::SetZDCLeftCentrCut(Float_t ZDCLeftCentrCut)
244 // Set default cut values for ZDC trigger
246 if(ZDCLeftCentrCut) fZDCLeftCentrCut = ZDCLeftCentrCut;
247 else fZDCLeftCentrCut = 10000.;
249 //________________________________________________________________
250 void AliZDCTrigger::SetZDCRightCentrCut(Float_t ZDCRightCentrCut)
252 // Set default cut values for ZDC trigger
254 if(ZDCRightCentrCut) fZDCRightCentrCut = ZDCRightCentrCut;
255 else fZDCRightCentrCut = 10000.;
257 //________________________________________________________________
258 void AliZDCTrigger::SetZDCLeftSemiCentrCut(Float_t ZDCLeftSemiCentrCut)
260 // Set default cut values for ZDC trigger
262 if(ZDCLeftSemiCentrCut) fZDCLeftSemiCentrCut = ZDCLeftSemiCentrCut;
263 else fZDCLeftSemiCentrCut = 18500.;
265 //________________________________________________________________
266 void AliZDCTrigger::SetZDCRightSemiCentrCut(Float_t ZDCRightSemiCentrCut)
268 // Set default cut values for ZDC trigger
270 if(ZDCRightSemiCentrCut) fZDCRightSemiCentrCut = ZDCRightSemiCentrCut;
271 else fZDCRightSemiCentrCut = 18500.;
273 //________________________________________________________________
274 void AliZDCTrigger::SetZEMCentrCut(Float_t ZEMCentrCut)
276 // Set default cut values for ZDC trigger
278 if(ZEMCentrCut) fZEMCentrCut = ZEMCentrCut;
279 else fZEMCentrCut = 210.;