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 //________________________________________________________________
34 ClassImp(AliZDCTrigger)
36 //________________________________________________________________
37 AliZDCTrigger::AliZDCTrigger() :
46 fZDCLeftSemiCentrCut(0),
47 fZDCRightSemiCentrCut(0),
54 SetZDCLeftEMDCuts(0,0);
55 SetZDCRightEMDCuts(0,0);
59 //________________________________________________________________
60 void AliZDCTrigger::CreateInputs()
64 // Do not create inputs again!!
65 if( fInputs.GetEntriesFast() > 0 ) return;
67 fInputs.AddLast(new AliTriggerInput("ZDC_1_L1", "ZDC", 1));
68 fInputs.AddLast(new AliTriggerInput("ZDC_2_L1", "ZDC", 1));
69 fInputs.AddLast(new AliTriggerInput("ZDC_3_L1", "ZDC", 1));
70 fInputs.AddLast(new AliTriggerInput("ZDC_EMD_L1", "ZDC", 1));
73 //________________________________________________________________
74 void AliZDCTrigger::Trigger()
79 AliRunLoader *runLoader = AliRunLoader::Instance();
81 AliLoader *aZDCLoader = runLoader->GetLoader("ZDCLoader");
83 aZDCLoader->LoadDigits("READ");
85 AliZDCDigit* pdigit = &digit;
86 TTree* tD = aZDCLoader->TreeD();
88 cerr<<"AliZDCTrigger: digits tree not found\n";
91 tD->SetBranchAddress("ZDC", &pdigit);
93 Float_t signalZNLeft[]={0,0}, signalZPLeft[]={0,0}, signalZDCLeftSum[]={0,0};
94 Float_t signalZNRight[]={0,0}, signalZPRight[]={0,0}, signalZDCRightSum[]={0,0};
95 Float_t signalZEMSum[]={0,0};
96 for(Int_t iDigit=0; iDigit<tD->GetEntries(); iDigit++){
100 if(digit.GetSector(0)==1)
101 for(Int_t i=0; i<2; i++){ //0=high range; 1=low range
102 signalZNLeft[i] += digit.GetADCValue(i);
103 signalZDCLeftSum[i] += digit.GetADCValue(i);
105 else if(digit.GetSector(0)==2)
106 for(Int_t i=0; i<2; i++){
107 signalZPLeft[i] += digit.GetADCValue(i);
108 signalZDCLeftSum[i] += digit.GetADCValue(i);
110 else if(digit.GetSector(0)==3)
111 for(Int_t i=0; i<2; i++) signalZEMSum[i] += digit.GetADCValue(i);
113 else if(digit.GetSector(0)==4)
114 for(Int_t i=0; i<2; i++){ //0=high range; 1=low range
115 signalZNRight[i] += digit.GetADCValue(i);
116 signalZDCRightSum[i] += digit.GetADCValue(i);
118 else if(digit.GetSector(0)==5)
119 for(Int_t i=0; i<2; i++){
120 signalZPRight[i] += digit.GetADCValue(i);
121 signalZDCRightSum[i] += digit.GetADCValue(i);
124 // *******************************************************************
125 if(signalZDCLeftSum[1]>fZDCLeftMBCut && signalZDCRightSum[1]>fZDCRightMBCut)
126 // *** ZDC minimum bias trigger
127 SetInput("ZDC_1_L1");
128 // *******************************************************************
129 if(signalZDCLeftSum[1]>fZDCLeftCentrCut && signalZDCLeftSum[1]<fZDCLeftSemiCentrCut &&
130 signalZDCRightSum[1]>fZDCRightCentrCut && signalZDCRightSum[1]<fZDCRightSemiCentrCut
131 && signalZEMSum[1]>fZEMCentrCut)
132 // *** ZDC semi-central (10-40%)
133 SetInput("ZDC_2_L1");
134 // *******************************************************************
135 if(signalZDCLeftSum[1]>fZDCLeftMinCut && signalZDCLeftSum[1]<fZDCLeftCentrCut &&
136 signalZDCRightSum[1]>fZDCRightMinCut && signalZDCRightSum[1]<fZDCRightCentrCut &&
137 signalZEMSum[1]>fZEMCentrCut)
138 // *** ZDC central (0-10%)
139 SetInput("ZDC_3_L1");
140 // *******************************************************************
141 if(signalZNLeft[0]>fZDCLeftEMDCuts[0] && signalZNLeft[0]<fZDCLeftEMDCuts[1] &&
142 signalZNRight[0]>fZDCRightEMDCuts[0] && signalZNRight[0]<fZDCRightEMDCuts[1] &&
143 signalZEMSum[1]<fZEMMinCut){ // *** 1n EMD trigger
144 SetInput("ZDC_EMD_L1");
149 //________________________________________________________________
150 void AliZDCTrigger::SetZDCLeftMinCut(Float_t ZDCLeftMinCut)
152 // Set default cut values for ZDC trigger
154 if(ZDCLeftMinCut) fZDCLeftMinCut = ZDCLeftMinCut;
155 else fZDCLeftMinCut = 800.;
157 //________________________________________________________________
158 void AliZDCTrigger::SetZDCRightMinCut(Float_t ZDCRightMinCut)
160 // Set default cut values for ZDC trigger
162 if(ZDCRightMinCut) fZDCRightMinCut = ZDCRightMinCut;
163 else fZDCRightMinCut = 800.;
166 //________________________________________________________________
167 void AliZDCTrigger::SetZEMMinCut(Float_t ZEMMinCut)
169 // Set default cut values for ZDC trigger
171 if(ZEMMinCut) fZEMMinCut = ZEMMinCut;
172 else fZEMMinCut = 80.;
174 //________________________________________________________________
175 void AliZDCTrigger::SetZDCLeftEMDCuts(Float_t* ZDCLeftEMDCuts)
177 // Set default cut values for ZDC trigger
179 if(ZDCLeftEMDCuts) for(int j=0; j<2; j++) fZDCLeftEMDCuts[j] = ZDCLeftEMDCuts[j];
181 fZDCLeftEMDCuts[0] = 600.;
182 fZDCLeftEMDCuts[1] = 1000.;
185 //________________________________________________________________
186 void AliZDCTrigger::SetZDCLeftEMDCuts(Float_t ZDCLeftEMDCutInf,
187 Float_t ZDCLeftEMDCutSup)
189 // Set default cut values for ZDC trigger
191 if(ZDCLeftEMDCutInf && ZDCLeftEMDCutSup){
192 fZDCLeftEMDCuts[0]=ZDCLeftEMDCutInf;
193 fZDCLeftEMDCuts[1]=ZDCLeftEMDCutSup;
196 fZDCLeftEMDCuts[0] = 600.;
197 fZDCLeftEMDCuts[1] = 1000.;
200 //________________________________________________________________
201 void AliZDCTrigger::SetZDCRightEMDCuts(Float_t* ZDCRightEMDCuts)
203 // Set default cut values for ZDC trigger
205 if(ZDCRightEMDCuts) for(int j=0; j<2; j++) fZDCRightEMDCuts[j] = ZDCRightEMDCuts[j];
207 fZDCRightEMDCuts[0] = 600.;
208 fZDCRightEMDCuts[1] = 1000.;
211 //________________________________________________________________
212 void AliZDCTrigger::SetZDCRightEMDCuts(Float_t ZDCRightEMDCutInf,
213 Float_t ZDCRightEMDCutSup)
215 // Set default cut values for ZDC trigger
217 if(ZDCRightEMDCutInf && ZDCRightEMDCutSup){
218 fZDCRightEMDCuts[0]=ZDCRightEMDCutInf;
219 fZDCRightEMDCuts[1]=ZDCRightEMDCutSup;
222 fZDCRightEMDCuts[0] = 600.;
223 fZDCRightEMDCuts[1] = 1000.;
226 //________________________________________________________________
227 void AliZDCTrigger::SetZDCLeftMBCut(Float_t ZDCLeftMBCut)
229 // Set default cut values for ZDC trigger
231 if(ZDCLeftMBCut) fZDCLeftMBCut = ZDCLeftMBCut;
232 else fZDCLeftMBCut = 800.;
234 //________________________________________________________________
235 void AliZDCTrigger::SetZDCRightMBCut(Float_t ZDCRightMBCut)
237 // Set default cut values for ZDC trigger
239 if(ZDCRightMBCut) fZDCRightMBCut = ZDCRightMBCut;
240 else fZDCRightMBCut = 800.;
242 //________________________________________________________________
243 void AliZDCTrigger::SetZDCLeftCentrCut(Float_t ZDCLeftCentrCut)
245 // Set default cut values for ZDC trigger
247 if(ZDCLeftCentrCut) fZDCLeftCentrCut = ZDCLeftCentrCut;
248 else fZDCLeftCentrCut = 10000.;
250 //________________________________________________________________
251 void AliZDCTrigger::SetZDCRightCentrCut(Float_t ZDCRightCentrCut)
253 // Set default cut values for ZDC trigger
255 if(ZDCRightCentrCut) fZDCRightCentrCut = ZDCRightCentrCut;
256 else fZDCRightCentrCut = 10000.;
258 //________________________________________________________________
259 void AliZDCTrigger::SetZDCLeftSemiCentrCut(Float_t ZDCLeftSemiCentrCut)
261 // Set default cut values for ZDC trigger
263 if(ZDCLeftSemiCentrCut) fZDCLeftSemiCentrCut = ZDCLeftSemiCentrCut;
264 else fZDCLeftSemiCentrCut = 18500.;
266 //________________________________________________________________
267 void AliZDCTrigger::SetZDCRightSemiCentrCut(Float_t ZDCRightSemiCentrCut)
269 // Set default cut values for ZDC trigger
271 if(ZDCRightSemiCentrCut) fZDCRightSemiCentrCut = ZDCRightSemiCentrCut;
272 else fZDCRightSemiCentrCut = 18500.;
274 //________________________________________________________________
275 void AliZDCTrigger::SetZEMCentrCut(Float_t ZEMCentrCut)
277 // Set default cut values for ZDC trigger
279 if(ZEMCentrCut) fZEMCentrCut = ZEMCentrCut;
280 else fZEMCentrCut = 210.;