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