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