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