]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliCentralTrigger.cxx
Return to the original methods from AliLoader (Marco)
[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 #include <TFile.h>
40 #include <TTree.h>
41
42 #include "AliLog.h"
43 #include "AliRun.h"
44 #include "AliRunLoader.h"
45 #include "AliModule.h"
46
47 #include "AliTriggerInput.h"
48 #include "AliTriggerDetector.h"
49 #include "AliTriggerCondition.h"
50 #include "AliTriggerDescriptor.h"
51 #include "AliCentralTrigger.h"
52 #include "AliDetectorEventHeader.h"
53 #include "AliHeader.h"
54
55
56 ClassImp( AliCentralTrigger )
57
58 //_____________________________________________________________________________
59 AliCentralTrigger::AliCentralTrigger() :
60    TObject(),
61    fClassMask(0)
62 {
63    // Default constructor
64 //   LoadDescriptor("Pb-Pb");
65 }
66
67 //_____________________________________________________________________________
68 AliCentralTrigger::AliCentralTrigger( TString & descriptor ) :
69    TObject(),
70    fClassMask(0)
71 {
72    // Default constructor
73    LoadDescriptor( descriptor );
74 }
75
76 //_____________________________________________________________________________
77 AliCentralTrigger::AliCentralTrigger( const AliCentralTrigger& ctp ):
78    TObject( ctp ),
79    fClassMask( ctp.fClassMask )
80 {
81    // Copy constructor
82
83    Int_t ndes = ctp.fDescriptors.GetEntriesFast();
84    for( Int_t j=0; j<ndes; j++ ) {
85       fDescriptors.AddLast( new AliTriggerDescriptor( *(AliTriggerDescriptor*)(ctp.fDescriptors.At( j )) ) );
86    }
87
88    Int_t nInp = ctp.fInputs.GetEntriesFast();
89    for( Int_t j=0; j<nInp; j++ ) {
90       fInputs.AddLast( new AliTriggerInput( *(AliTriggerInput*)(ctp.fInputs.At( j )) ) );
91    }
92 }
93
94 //_____________________________________________________________________________
95 AliCentralTrigger::~AliCentralTrigger()
96 {
97    // Destructor
98    DeleteDescriptors();
99 }
100
101 //_____________________________________________________________________________
102 void AliCentralTrigger::DeleteDescriptors()
103 {
104    // Reset Descriptors
105    fClassMask = 0;
106    fDescriptors.SetOwner();
107    fDescriptors.Delete();
108 }
109
110 //_____________________________________________________________________________
111 void AliCentralTrigger::Reset()
112 {
113    // Reset Class Mask and conditions
114    fClassMask = 0;
115    Int_t ndes = fDescriptors.GetEntriesFast();
116    for( Int_t i=0; i<ndes; i++ ) {
117       TObjArray* condArray = ((AliTriggerDescriptor*)fDescriptors.At( i ))->GetTriggerConditions();
118       Int_t ncond = condArray->GetEntriesFast();
119       for( Int_t j=0; j<ncond; j++ ) {
120          AliTriggerCondition* cond = (AliTriggerCondition*)condArray->At( j );
121          cond->Reset();
122       }
123    }
124 }
125
126 //_____________________________________________________________________________
127 void AliCentralTrigger::MakeBranch( TString name, TTree * tree )
128 {
129    // Make a branch to store only trigger class mask event by event
130
131    if( tree )  {
132       AliDebug( 1, "Got Tree from folder." );
133       TBranch* branch = tree->GetBranch( name );
134       if( branch == 0x0 ) {
135          AliDebug( 1, "Creating new branch" );
136          branch = tree->Branch( name, &(this->fClassMask), "fClassMask/l" );
137          branch->SetAutoDelete( kFALSE );
138       }
139       else {
140          AliDebug( 1, "Got Branch from Tree" );
141          branch->SetAddress( &(this->fClassMask) );
142       }
143    }
144 }
145
146 //_____________________________________________________________________________
147 Bool_t AliCentralTrigger::LoadDescriptor( TString & descriptor )
148 {
149    // Load one or more pre-created Descriptors from database/file that match
150    // with the input string 'descriptor'
151    // Ej: "Pb-Pb" or "p-p-DIMUON CALIBRATION-CENTRAL-BARREL"
152
153    // Delete any descriptor
154    Reset();
155
156    // Load the selected descriptors
157    TObjArray* desArray = descriptor.Tokenize( " " );
158    Int_t ndes = desArray->GetEntriesFast();
159    for( Int_t i=0; i<ndes; i++ ) {
160       TObjString* val = (TObjString*)desArray->At( i );
161       AliTriggerDescriptor* des = AliTriggerDescriptor::LoadDescriptor( val->String() );
162       if( des ) 
163          fDescriptors.AddLast( des );
164       else
165          AliWarning( Form( "Descriptor (%s) not found", val->String().Data() ) );
166    }
167    Bool_t desfound = kTRUE;
168    if( fDescriptors.GetEntriesFast() == 0 ) desfound = kFALSE;
169
170    delete desArray;
171
172    return desfound;
173 }
174
175 //_____________________________________________________________________________
176 TString AliCentralTrigger::GetDetectors()
177 {
178    // return TString with the detectors to be trigger
179    // merging detectors from all descriptors
180
181    TString result;
182
183    Int_t ndes = fDescriptors.GetEntriesFast();
184    for( Int_t i=0; i<ndes; i++ ) {
185       TString detStr = ((AliTriggerDescriptor*)fDescriptors.At( i ))->GetDetectorCluster();
186       TObjArray* det = detStr.Tokenize(" ");
187       Int_t ndet = det->GetEntriesFast();
188       for( Int_t j=0; j<ndet; j++ ) {
189          if( result.Contains( ((TObjString*)det->At(j))->String() ) )continue;
190          result.Append( " " );
191          result.Append( ((TObjString*)det->At(j))->String() );
192       }
193    }
194
195    return result;
196 }
197
198 //_____________________________________________________________________________
199 Bool_t AliCentralTrigger::RunTrigger( AliRunLoader* runLoader )
200 {
201    // run the trigger
202
203    if( fDescriptors.GetEntriesFast() == 0 ) {
204       AliError( "not trigger descriptor loaded, skipping trigger" );
205       return kFALSE;
206    }
207
208    TTree *tree = runLoader->TreeCT();
209    if( !tree ) {
210       AliError( "not folder with trigger loaded, skipping trigger" );
211       return kFALSE;
212    }
213
214    TStopwatch stopwatch;
215    stopwatch.Start();
216
217    // Process each event
218    for( Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++ ) {
219       AliInfo( Form("  ***** Processing event %d *****\n", iEvent) );
220       runLoader->GetEvent( iEvent );
221       // Get detectors involve
222       TString detStr = GetDetectors();
223       AliInfo( Form(" Cluster Detectors %s \n", detStr.Data() ) );
224       TObjArray* detArray = runLoader->GetAliRun()->Detectors();
225       // Reset Mask
226       fClassMask = 0;
227       TObjArray trgdetArray; // use as garbage collector
228       for( Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++ ) {
229          AliModule* det = (AliModule*) detArray->At( iDet );
230          if( !det || !det->IsActive() ) continue;
231          if( IsSelected(det->GetName(), detStr) ) {
232
233             AliInfo( Form("Triggering from digits for %s", det->GetName() ) );
234             AliTriggerDetector* trgdet = det->CreateTriggerDetector();
235             trgdet->CreateInputs();
236             TStopwatch stopwatchDet;
237             stopwatchDet.Start();
238             trgdet->Trigger();
239             AliInfo( Form("Execution time for %s: R:%.2fs C:%.2fs",
240                      det->GetName(), stopwatchDet.RealTime(), stopwatchDet.CpuTime() ) );
241
242             // Get the inputs
243             TObjArray* detInp = trgdet->GetInputs();
244             for( Int_t i=0; i<detInp->GetEntriesFast(); i++ ) {
245                fInputs.AddLast( detInp->At(i) );
246             }
247             trgdetArray.AddLast( trgdet );
248
249             // Write trigger detector in Event folder in Digits file
250             TString loadername = det->GetName();
251             loadername.Append( "Loader" );
252             AliLoader * loader = runLoader->GetLoader( loadername );
253             if( loader ) {
254                AliDataLoader * dataLoader = loader->GetDigitsDataLoader();
255                if( !dataLoader->IsFileOpen() ) {
256                   if( dataLoader->OpenFile( "UPDATE" ) ) {
257                      AliWarning( Form( "\n\nCan't write trigger for %s\n", det->GetName() ) );
258                   }
259                }
260                dataLoader->Cd();
261                if( gFile && !gFile->IsWritable() ) {
262                   gFile->ReOpen( "UPDATE" );
263                   dataLoader->Cd();
264                }
265                trgdet->Write( "Trigger", TObject::kOverwrite );
266                dataLoader->CloseFile();
267             }
268             else  AliWarning( Form( "Not loader found for %s", det->GetName() ) );
269          }
270       }
271
272       // Check trigger conditions and create the trigger class mask
273       CheckConditions();
274       fInputs.Clear();
275
276       // Clear trigger detectors
277       trgdetArray.SetOwner();
278       trgdetArray.Delete();
279
280       if( (detStr.CompareTo( "ALL" ) != 0) && !detStr.IsNull() ) {
281          AliError( Form("the following detectors were not found: %s",
282                    detStr.Data()));
283          return kFALSE;
284       }
285
286       // Save trigger mask
287       tree->Fill();
288       AliInfo( Form("**************** Central Trigger Class Mask:0x%X", fClassMask ) );
289    } // end event loop
290
291    Reset();
292 //   cout << endl <<  " Print " << endl;
293 //   Print();
294
295    // Write
296    runLoader->WriteTrigger( "OVERWRITE" );
297
298    return kTRUE;
299 }
300
301 //_____________________________________________________________________________
302 Long_t AliCentralTrigger::CheckConditions()
303 {
304    // Check trigger conditions and create the trigger class mask
305
306    Int_t ndes = fDescriptors.GetEntriesFast();
307    for( Int_t i=0; i<ndes; i++ ) {
308       TObjArray* condArray = ((AliTriggerDescriptor*)fDescriptors.At( i ))->GetTriggerConditions();
309       Int_t ncond = condArray->GetEntriesFast();
310       for( Int_t j=0; j<ncond; j++ ) {
311          AliTriggerCondition* cond = (AliTriggerCondition*)condArray->At( j );
312          if( !cond->CheckInputs( fInputs ) ) continue;
313          cond->Trigger( fInputs );
314     //     cond->Print();
315          fClassMask |= cond->GetValue();
316       }
317    }
318    return fClassMask;
319 }
320 //_____________________________________________________________________________
321 TObjArray* AliCentralTrigger::GetResultConditions()
322 {
323    // return only the true conditions
324
325    TObjArray* result = new TObjArray();
326
327    Int_t ndes = fDescriptors.GetEntriesFast();
328    for( Int_t i=0; i<ndes; i++ ) {
329       TObjArray* condArray = ((AliTriggerDescriptor*)fDescriptors.At( i ))->GetTriggerConditions();
330       Int_t ncond = condArray->GetEntriesFast();
331       for( Int_t j=0; j<ncond; j++ ) {
332          AliTriggerCondition* cond = (AliTriggerCondition*)condArray->At( j );
333          if( cond->GetStatus() ) result->AddLast( cond );
334       }
335    }
336
337    return result;
338 }
339
340 //_____________________________________________________________________________
341 void AliCentralTrigger::Print( const Option_t*  ) const
342 {
343    // Print
344    cout << "Central Trigger: " << endl;
345    cout << "  Trigger Class Mask: 0x" << hex << fClassMask << dec << endl;
346
347    Int_t ndes = fDescriptors.GetEntriesFast();
348    for( Int_t i=0; i<ndes; i++ ) {
349       AliTriggerDescriptor* des = (AliTriggerDescriptor*)fDescriptors.At( i );
350       if( des ) des->Print();
351    }
352    cout << endl;
353 }
354
355
356 //////////////////////////////////////////////////////////////////////////////
357 // Helper method
358
359 //_____________________________________________________________________________
360 Bool_t AliCentralTrigger::IsSelected( TString detName, TString& detectors ) const
361 {
362    // check whether detName is contained in detectors
363    // if yes, it is removed from detectors
364
365    // check if all detectors are selected
366    if( (detectors.CompareTo("ALL") == 0 ) ||
367         detectors.BeginsWith("ALL ") ||
368         detectors.EndsWith(" ALL") ||
369         detectors.Contains(" ALL ") ) {
370       detectors = "ALL";
371       return kTRUE;
372    }
373
374    // search for the given detector
375    Bool_t result = kFALSE;
376    if( (detectors.CompareTo( detName ) == 0) ||
377         detectors.BeginsWith( detName+" " ) ||
378         detectors.EndsWith( " "+detName ) ||
379         detectors.Contains( " "+detName+" " ) ) {
380       detectors.ReplaceAll( detName, "" );
381       result = kTRUE;
382    }
383
384    // clean up the detectors string
385    while( detectors.Contains("  ") )  detectors.ReplaceAll( "  ", " " );
386    while( detectors.BeginsWith(" ") ) detectors.Remove( 0, 1 );
387    while( detectors.EndsWith(" ") )   detectors.Remove( detectors.Length()-1, 1 );
388
389    return result;
390 }