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