]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliCentralTrigger.cxx
Corrected LUT for the alignable volume paths in TRD (R.Grosso)
[u/mrichter/AliRoot.git] / STEER / AliCentralTrigger.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 /* $Id$ */
17
18 ///////////////////////////////////////////////////////////////////////////////
19 //                                                                           //
20 // This class for running the Central Trigger Processor                      //
21 //                                                                           //
22 //                                                                           //
23 //    Load Descriptors                                                       //
24 //    Make a list the trigger detectors involve from the descriptors         //
25 //    For the each event                                                     //
26 //           Run the Trigger for the each detector                           //
27 //           Get the inputs                                                  //
28 //           Check the condition classes                                     //
29 //           Create the class mask                                           //
30 //           Save result                                                     //
31 //                                                                           //
32 //                                                                           //
33 ///////////////////////////////////////////////////////////////////////////////
34
35 #include <TString.h>
36 #include <TObjString.h>
37 #include <TObjArray.h>
38 #include <TStopwatch.h>
39
40 #include "AliLog.h"
41 #include "AliRun.h"
42 #include "AliRunLoader.h"
43 #include "AliModule.h"
44
45 #include "AliTriggerInput.h"
46 #include "AliTriggerDetector.h"
47 #include "AliTriggerCondition.h"
48 #include "AliTriggerDescriptor.h"
49 #include "AliCentralTrigger.h"
50
51 ClassImp( AliCentralTrigger )
52
53 //_____________________________________________________________________________
54 AliCentralTrigger::AliCentralTrigger() :
55    TObject(),
56    fClassMask(0)
57 {
58    // Default constructor
59 //   LoadDescriptor("Pb-Pb");
60 }
61
62 //_____________________________________________________________________________
63 AliCentralTrigger::AliCentralTrigger( TString & descriptor ) :
64    TObject(),
65    fClassMask(0)
66 {
67    // Default constructor
68    LoadDescriptor( descriptor );
69 }
70
71 //_____________________________________________________________________________
72 AliCentralTrigger::~AliCentralTrigger()
73 {
74    // Destructor
75    fDescriptors.SetOwner();
76    fDescriptors.Delete();
77 }
78
79 //_____________________________________________________________________________
80 Bool_t AliCentralTrigger::LoadDescriptor( TString & descriptor )
81 {
82    // Load one or more pre-created Descriptors from database/file that match
83    // with the input string 'descriptor'
84    // Ej: "Pb-Pb" or "p-p-DIMUON CALIBRATION-CENTRAL-BARREL"
85
86    // Delete any descriptor
87    fDescriptors.Delete();
88
89    // Load the selected descriptors
90    TObjArray* desArray = descriptor.Tokenize( " " );
91    Int_t ndes = desArray->GetEntriesFast();
92    for( Int_t i=0; i<ndes; i++ ) {
93       TObjString* val = (TObjString*)desArray->At( i );
94       AliTriggerDescriptor* des = AliTriggerDescriptor::LoadDescriptor( val->String() );
95       if( des ) 
96          fDescriptors.AddLast( des );
97       else
98          AliWarning( Form( "Descriptor (%s) not found", val->String().Data() ) );
99    }
100    Bool_t desfound = kTRUE;
101    if( fDescriptors.GetEntriesFast() == 0 ) desfound = kFALSE;
102
103    delete desArray;
104
105    return desfound;
106 }
107
108 //_____________________________________________________________________________
109 TString AliCentralTrigger::GetDetectors()
110 {
111    // return TString with the detectors to be trigger
112    // merging detectors from all descriptors
113
114    TString result;
115
116    Int_t ndes = fDescriptors.GetEntriesFast();
117    for( Int_t i=0; i<ndes; i++ ) {
118       TString detStr = ((AliTriggerDescriptor*)fDescriptors.At( i ))->GetDetectorCluster();
119       TObjArray* det = detStr.Tokenize(" ");
120       Int_t ndet = det->GetEntriesFast();
121       for( Int_t j=0; j<ndet; j++ ) {
122          if( result.Contains( ((TObjString*)det->At(j))->String() ) )continue;
123          result.Append( " " );
124          result.Append( ((TObjString*)det->At(j))->String() );
125       }
126    }
127
128    return result;
129 }
130
131 //_____________________________________________________________________________
132 Bool_t AliCentralTrigger::RunTrigger( AliRunLoader* runLoader )
133 {
134    // run the trigger
135    
136    if( fDescriptors.GetEntriesFast() == 0 ) {
137       AliError("not trigger descriptor loaded, skipping trigger");
138       return kFALSE;
139    }
140
141    TStopwatch stopwatch;
142    stopwatch.Start();
143    
144    // Process each event
145    for( Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++ ) {
146       AliInfo( Form("\n ***** Processing event %d *****\n", iEvent) );
147       runLoader->GetEvent( iEvent );       
148       // Get detectors involve
149       TString detStr = GetDetectors();
150       TObjArray* detArray = runLoader->GetAliRun()->Detectors();
151       // Reset Mask
152       fClassMask = 0;
153       TObjArray trgdetArray;
154       for( Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++ ) {
155          AliModule* det = (AliModule*) detArray->At( iDet );
156          if( !det || !det->IsActive() ) continue;
157          if( IsSelected(det->GetName(), detStr) ) {
158             AliInfo( Form("triggering from digits for %s", det->GetName() ) );
159             AliTriggerDetector* trgdet = det->CreateTriggerDetector();
160             trgdet->CreateInputs();
161             trgdetArray.AddLast( trgdet );
162             TStopwatch stopwatchDet;
163             stopwatchDet.Start();
164             trgdet->Trigger();
165             AliInfo( Form("Execution time for %s: R:%.2fs C:%.2fs",
166                      det->GetName(), stopwatchDet.RealTime(), stopwatchDet.CpuTime() ) );
167             // Get the inputs
168             TObjArray* detInp = trgdet->GetInputs();
169             for( Int_t i=0; i<detInp->GetEntriesFast(); i++ ) {
170                fInputs.AddLast( detInp->At(i) );
171             }
172          }
173       }
174
175       // Check trigger conditions and create the trigger class mask
176       CheckConditions();
177
178       fInputs.Clear();
179       
180       // Clear trigger detectors
181       trgdetArray.SetOwner();
182       trgdetArray.Delete();
183
184       if( (detStr.CompareTo( "ALL" ) != 0) && !detStr.IsNull() ) {
185          AliError( Form("the following detectors were not found: %s",
186                    detStr.Data()));
187          return kFALSE;
188       }
189
190       // Write trigger ???? -> Event Header
191       // Write();
192       Print();
193
194    } // end event loop
195    return kTRUE;
196 }
197
198
199
200 //_____________________________________________________________________________
201 Long_t AliCentralTrigger::CheckConditions()
202 {
203    // Check trigger conditions and create the trigger class mask
204
205    Int_t ndes = fDescriptors.GetEntriesFast();
206    for( Int_t i=0; i<ndes; i++ ) {
207       TObjArray* condArray = ((AliTriggerDescriptor*)fDescriptors.At( i ))->GetTriggerConditions();
208       Int_t ncond = condArray->GetEntriesFast();
209       for( Int_t j=0; j<ncond; j++ ) {
210          AliTriggerCondition* cond = (AliTriggerCondition*)condArray->At( j );
211          if( !cond->CheckInputs( fInputs ) ) continue;
212          cond->Trigger( fInputs );
213     //     cond->Print();
214          fClassMask |= cond->GetValue();
215       }
216    }
217    return fClassMask;
218 }
219 //_____________________________________________________________________________
220 TObjArray* AliCentralTrigger::GetResultConditions()
221 {
222    // return only the true conditions
223
224    TObjArray* result = new TObjArray();
225
226    Int_t ndes = fDescriptors.GetEntriesFast();
227    for( Int_t i=0; i<ndes; i++ ) {
228       TObjArray* condArray = ((AliTriggerDescriptor*)fDescriptors.At( i ))->GetTriggerConditions();
229       Int_t ncond = condArray->GetEntriesFast();
230       for( Int_t j=0; j<ncond; j++ ) {
231          AliTriggerCondition* cond = (AliTriggerCondition*)condArray->At( j );
232          if( cond->GetStatus() ) result->AddLast( cond );
233       }
234    }
235
236    return result;
237 }
238
239
240 //_____________________________________________________________________________
241 void AliCentralTrigger::Print( const Option_t*  ) const
242 {
243    // Print
244    cout << "Central Trigger: " << endl;
245    cout << "  Trigger Class Mask: 0x" << hex << fClassMask << dec << endl;
246
247    Int_t ndes = fDescriptors.GetEntriesFast();
248    for( Int_t i=0; i<ndes; i++ ) {
249       AliTriggerDescriptor* des = (AliTriggerDescriptor*)fDescriptors.At( i );
250       if( des ) des->Print();
251    }
252    cout << endl;
253 }
254
255
256 //////////////////////////////////////////////////////////////////////////////
257 // Helper method
258
259 //_____________________________________________________________________________
260 Bool_t AliCentralTrigger::IsSelected( TString detName, TString& detectors ) const
261 {
262    // check whether detName is contained in detectors
263    // if yes, it is removed from detectors
264
265    // check if all detectors are selected
266    if( (detectors.CompareTo("ALL") == 0 ) ||
267         detectors.BeginsWith("ALL ") ||
268         detectors.EndsWith(" ALL") ||
269         detectors.Contains(" ALL ") ) {
270       detectors = "ALL";
271       return kTRUE;
272    }
273
274    // search for the given detector
275    Bool_t result = kFALSE;
276    if( (detectors.CompareTo( detName ) == 0) ||
277         detectors.BeginsWith( detName+" " ) ||
278         detectors.EndsWith( " "+detName ) ||
279         detectors.Contains( " "+detName+" " ) ) {
280       detectors.ReplaceAll( detName, "" );
281       result = kTRUE;
282    }
283
284    // clean up the detectors string
285    while( detectors.Contains("  ") )  detectors.ReplaceAll( "  ", " " );
286    while( detectors.BeginsWith(" ") ) detectors.Remove( 0, 1 );
287    while( detectors.EndsWith(" ") )   detectors.Remove( detectors.Length()-1, 1 );
288
289    return result;
290 }