]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ZDC/AliZDCTrigger.cxx
Corrected input definition of the ZDC trigger (Chiara)
[u/mrichter/AliRoot.git] / ZDC / AliZDCTrigger.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                         *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
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  **************************************************************************/
15
16 // ****************************************************************
17 //
18 //      Trigger class for ZDC
19 //
20 // ****************************************************************
21
22 #include "AliLog.h"
23 #include "AliRun.h"
24 #include "AliLoader.h"
25 #include "AliRunLoader.h"
26 #include "AliTriggerInput.h"
27
28 #include "AliZDC.h"
29 #include "AliZDCDigit.h"
30 #include "AliZDCTrigger.h"
31
32 //________________________________________________________________
33 ClassImp(AliZDCTrigger)
34
35 //________________________________________________________________
36 AliZDCTrigger::AliZDCTrigger() : AliTriggerDetector() 
37 {  
38    // Constructor
39    SetName("ZDC");
40    CreateInputs();
41    //
42    SetZDCLeftMinCut(0);
43    SetZDCRightMinCut(0);
44    SetZEMMinCut(0);
45    SetZDCLeftEMDCuts(0,0);
46    SetZDCRightEMDCuts(0,0);
47    SetZDCLeftMBCut(0);
48    SetZDCRightMBCut(0);
49    SetZDCLeftCentrCut(0);
50    SetZDCRightCentrCut(0);
51    SetZDCLeftSemiCentrCut(0);
52    SetZDCRightSemiCentrCut(0);
53    SetZEMCentrCut(0);
54
55 }
56
57 //________________________________________________________________
58 void AliZDCTrigger::CreateInputs()
59 {
60    // Trigger inputs
61    
62    // Do not create inputs again!!
63    if( fInputs.GetEntriesFast() > 0 ) return;
64    
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));
69 }
70
71 //________________________________________________________________
72 void AliZDCTrigger::Trigger()
73 {
74
75    // Trigger selection
76    //
77    AliRunLoader *runLoader = gAlice->GetRunLoader();
78
79    AliLoader *aZDCLoader = runLoader->GetLoader("ZDCLoader");
80    aZDCLoader->LoadDigits("READ");
81    AliZDCDigit digit;
82    AliZDCDigit* pdigit = &digit;
83    TTree* tD = aZDCLoader->TreeD();
84    if (!tD) cerr<<"AliZDCTrigger: digits tree not found\n";
85    tD->SetBranchAddress("ZDC", &pdigit);
86    //
87    Float_t signalZNLeft[2], signalZPLeft[2], signalZDCLeftSum[2];
88    Float_t signalZNRight[2], signalZPRight[2], signalZDCRightSum[2];
89    Float_t signalZEMSum[2];
90    for(Int_t iDigit=0; iDigit<tD->GetEntries(); iDigit++){
91       tD->GetEntry(iDigit);
92       //
93       // *** ZDC LEFT
94       if(digit.GetSector(0)==1)
95          for(Int_t i=0; i<2; i++){ //0=high range; 1=low range
96             signalZNLeft[i] += digit.GetADCValue(i);
97             signalZDCLeftSum[i] += digit.GetADCValue(i);
98           }
99       else if(digit.GetSector(0)==2)
100          for(Int_t i=0; i<2; i++){
101             signalZPLeft[i] += digit.GetADCValue(i);
102             signalZDCLeftSum[i] += digit.GetADCValue(i);
103           }
104       else if(digit.GetSector(0)==3)
105          for(Int_t i=0; i<2; i++) signalZEMSum[i] += digit.GetADCValue(i);
106       // *** ZDC RIGHT
107       else if(digit.GetSector(0)==4)
108          for(Int_t i=0; i<2; i++){ //0=high range; 1=low range
109             signalZNRight[i] += digit.GetADCValue(i);
110             signalZDCRightSum[i] += digit.GetADCValue(i);
111           }
112       else if(digit.GetSector(0)==5)
113          for(Int_t i=0; i<2; i++){
114             signalZPRight[i] += digit.GetADCValue(i);
115             signalZDCRightSum[i] += digit.GetADCValue(i);
116           }
117    }
118    // *******************************************************************
119    if(signalZDCLeftSum[1]>fZDCLeftMBCut && signalZDCRightSum[1]>fZDCRightMBCut) 
120        // *** ZDC minimum bias trigger
121        SetInput("ZDC_1_L1");
122    // *******************************************************************
123    if(signalZDCLeftSum[1]>fZDCLeftCentrCut && signalZDCLeftSum[1]<fZDCLeftSemiCentrCut &&
124       signalZDCRightSum[1]>fZDCRightCentrCut && signalZDCRightSum[1]<fZDCRightSemiCentrCut
125       && signalZEMSum[1]>fZEMCentrCut) 
126        // *** ZDC semi-central (10-40%)
127        SetInput("ZDC_2_L1");
128    // *******************************************************************
129    if(signalZDCLeftSum[1]>fZDCLeftMinCut && signalZDCLeftSum[1]<fZDCLeftCentrCut &&
130       signalZDCRightSum[1]>fZDCRightMinCut && signalZDCRightSum[1]<fZDCRightCentrCut &&
131       signalZEMSum[1]>fZEMCentrCut) 
132        // *** ZDC central (0-10%)
133        SetInput("ZDC_3_L1");
134    // *******************************************************************
135    if(signalZNLeft[0]>fZDCLeftEMDCuts[0] && signalZNLeft[0]<fZDCLeftEMDCuts[1] && 
136       signalZNRight[0]>fZDCRightEMDCuts[0] && signalZNRight[0]<fZDCRightEMDCuts[1] &&
137       signalZEMSum[1]<fZEMMinCut){ // *** 1n EMD trigger
138         SetInput("ZDC_EMD_L1");
139    }
140    
141 }
142
143 //________________________________________________________________
144 void AliZDCTrigger::SetZDCLeftMinCut(Float_t ZDCLeftMinCut) 
145 {
146   // Set default cut values for ZDC trigger
147   //
148   if(ZDCLeftMinCut)  fZDCLeftMinCut = ZDCLeftMinCut;
149   else  fZDCLeftMinCut = 800.;
150 }
151 //________________________________________________________________
152 void AliZDCTrigger::SetZDCRightMinCut(Float_t ZDCRightMinCut) 
153 {
154   // Set default cut values for ZDC trigger
155   //
156   if(ZDCRightMinCut)  fZDCRightMinCut = ZDCRightMinCut;
157   else  fZDCRightMinCut = 800.;
158 }
159
160 //________________________________________________________________
161 void AliZDCTrigger::SetZEMMinCut(Float_t ZEMMinCut) 
162 {
163   // Set default cut values for ZDC trigger
164   //
165   if(ZEMMinCut)  fZEMMinCut = ZEMMinCut;
166   else  fZEMMinCut = 80.;
167 }
168 //________________________________________________________________
169 void AliZDCTrigger::SetZDCLeftEMDCuts(Float_t* ZDCLeftEMDCuts) 
170 {
171   // Set default cut values for ZDC trigger
172   //
173   if(ZDCLeftEMDCuts) for(int j=0; j<2; j++) fZDCLeftEMDCuts[j] = ZDCLeftEMDCuts[j];
174   else{
175     fZDCLeftEMDCuts[0] = 600.;
176     fZDCLeftEMDCuts[1] = 1000.;
177   }
178 }
179 //________________________________________________________________
180 void AliZDCTrigger::SetZDCLeftEMDCuts(Float_t ZDCLeftEMDCutInf, 
181         Float_t ZDCLeftEMDCutSup) 
182 {
183   // Set default cut values for ZDC trigger
184   //
185   if(ZDCLeftEMDCutInf && ZDCLeftEMDCutSup){
186     fZDCLeftEMDCuts[0]=ZDCLeftEMDCutInf; 
187     fZDCLeftEMDCuts[1]=ZDCLeftEMDCutSup;
188   }     
189   else{
190     fZDCLeftEMDCuts[0] = 600.;
191     fZDCLeftEMDCuts[1] = 1000.;
192   }
193 }
194 //________________________________________________________________
195 void AliZDCTrigger::SetZDCRightEMDCuts(Float_t* ZDCRightEMDCuts) 
196 {
197   // Set default cut values for ZDC trigger
198   //
199   if(ZDCRightEMDCuts) for(int j=0; j<2; j++) fZDCRightEMDCuts[j] = ZDCRightEMDCuts[j];
200   else{
201     fZDCRightEMDCuts[0] = 600.;
202     fZDCRightEMDCuts[1] = 1000.;
203   }
204 }
205 //________________________________________________________________
206 void AliZDCTrigger::SetZDCRightEMDCuts(Float_t ZDCRightEMDCutInf, 
207         Float_t ZDCRightEMDCutSup) 
208 {
209   // Set default cut values for ZDC trigger
210   //
211   if(ZDCRightEMDCutInf && ZDCRightEMDCutSup){
212     fZDCRightEMDCuts[0]=ZDCRightEMDCutInf; 
213     fZDCRightEMDCuts[1]=ZDCRightEMDCutSup;
214   }     
215   else{
216     fZDCRightEMDCuts[0] = 600.;
217     fZDCRightEMDCuts[1] = 1000.;
218   }
219 }
220 //________________________________________________________________
221 void AliZDCTrigger::SetZDCLeftMBCut(Float_t ZDCLeftMBCut) 
222 {
223   // Set default cut values for ZDC trigger
224   //
225   if(ZDCLeftMBCut) fZDCLeftMBCut = ZDCLeftMBCut;
226   else fZDCLeftMBCut = 800.;
227 }
228 //________________________________________________________________
229 void AliZDCTrigger::SetZDCRightMBCut(Float_t ZDCRightMBCut) 
230 {
231   // Set default cut values for ZDC trigger
232   //
233   if(ZDCRightMBCut) fZDCRightMBCut = ZDCRightMBCut;
234   else fZDCRightMBCut = 800.;
235 }
236 //________________________________________________________________
237 void AliZDCTrigger::SetZDCLeftCentrCut(Float_t ZDCLeftCentrCut) 
238 {
239   // Set default cut values for ZDC trigger
240   //
241   if(ZDCLeftCentrCut) fZDCLeftCentrCut = ZDCLeftCentrCut;
242   else fZDCLeftCentrCut = 10000.;
243 }
244 //________________________________________________________________
245 void AliZDCTrigger::SetZDCRightCentrCut(Float_t ZDCRightCentrCut) 
246 {
247   // Set default cut values for ZDC trigger
248   //
249   if(ZDCRightCentrCut) fZDCRightCentrCut = ZDCRightCentrCut;
250   else fZDCRightCentrCut = 10000.;
251 }
252 //________________________________________________________________
253 void AliZDCTrigger::SetZDCLeftSemiCentrCut(Float_t ZDCLeftSemiCentrCut) 
254 {
255   // Set default cut values for ZDC trigger
256   //
257   if(ZDCLeftSemiCentrCut) fZDCLeftSemiCentrCut = ZDCLeftSemiCentrCut;
258   else fZDCLeftSemiCentrCut = 18500.;
259 }
260 //________________________________________________________________
261 void AliZDCTrigger::SetZDCRightSemiCentrCut(Float_t ZDCRightSemiCentrCut) 
262 {
263   // Set default cut values for ZDC trigger
264   //
265   if(ZDCRightSemiCentrCut) fZDCRightSemiCentrCut = ZDCRightSemiCentrCut;
266   else fZDCRightSemiCentrCut = 18500.;
267 }
268 //________________________________________________________________
269 void AliZDCTrigger::SetZEMCentrCut(Float_t ZEMCentrCut) 
270 {
271   // Set default cut values for ZDC trigger
272   //
273   if(ZEMCentrCut) fZEMCentrCut = ZEMCentrCut;
274   else fZEMCentrCut = 210.;
275 }