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