]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliCentralTrigger.cxx
The default thickness of the chips is set to 150 mkm (D.Elia). Removing some obsolete...
[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 UChar_t AliCentralTrigger::GetClusterMask()
200 {
201    // Return the detector cluster mask following
202    // table 4.3 pag 60, TDR Trigger and DAQ
203
204    TString detStr = GetDetectors();
205    TObjArray* det = detStr.Tokenize(" ");
206    Int_t ndet = det->GetEntriesFast();
207
208    UInt_t idmask = 0;
209    if( ndet >= 8 ) {  // All detectors, should be 9 but CRT is not implemented yet
210       idmask = 1;
211       return idmask;
212    }
213
214    if( ndet >= 7 && !detStr.Contains("MUON") ) {  // Central Barrel, All but MUON
215       idmask = 2;
216       return idmask;
217    }
218
219    if( detStr.Contains("MUON") && detStr.Contains("START") ) {  // MUON arm
220       idmask = 4;
221       return idmask;
222    }
223
224    return idmask; // 0 something else!!!
225 }
226 //_____________________________________________________________________________
227 Bool_t AliCentralTrigger::RunTrigger( AliRunLoader* runLoader )
228 {
229    // run the trigger
230
231    if( fDescriptors.GetEntriesFast() == 0 ) {
232       AliError( "not trigger descriptor loaded, skipping trigger" );
233       return kFALSE;
234    }
235
236    TTree *tree = runLoader->TreeCT();
237    if( !tree ) {
238       AliError( "not folder with trigger loaded, skipping trigger" );
239       return kFALSE;
240    }
241
242    TStopwatch stopwatch;
243    stopwatch.Start();
244
245    // Process each event
246    for( Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++ ) {
247       AliInfo( Form("  ***** Processing event %d *****\n", iEvent) );
248       runLoader->GetEvent( iEvent );
249       // Get detectors involve
250       TString detStr = GetDetectors();
251       AliInfo( Form(" Cluster Detectors %s \n", detStr.Data() ) );
252       TObjArray* detArray = runLoader->GetAliRun()->Detectors();
253       // Reset Mask
254       fClassMask = 0;
255       TObjArray trgdetArray; // use as garbage collector
256       for( Int_t iDet = 0; iDet < detArray->GetEntriesFast(); iDet++ ) {
257          AliModule* det = (AliModule*) detArray->At( iDet );
258          if( !det || !det->IsActive() ) continue;
259          if( IsSelected(det->GetName(), detStr) ) {
260
261             AliInfo( Form("Triggering from digits for %s", det->GetName() ) );
262             AliTriggerDetector* trgdet = det->CreateTriggerDetector();
263             trgdet->CreateInputs();
264             TStopwatch stopwatchDet;
265             stopwatchDet.Start();
266             trgdet->Trigger();
267             AliInfo( Form("Execution time for %s: R:%.2fs C:%.2fs",
268                      det->GetName(), stopwatchDet.RealTime(), stopwatchDet.CpuTime() ) );
269
270             // Get the inputs
271             TObjArray* detInp = trgdet->GetInputs();
272             for( Int_t i=0; i<detInp->GetEntriesFast(); i++ ) {
273                fInputs.AddLast( detInp->At(i) );
274             }
275             trgdetArray.AddLast( trgdet );
276
277             // Write trigger detector in Event folder in Digits file
278             TString loadername = det->GetName();
279             loadername.Append( "Loader" );
280             AliLoader * loader = runLoader->GetLoader( loadername );
281             if( loader ) {
282                AliDataLoader * dataLoader = loader->GetDigitsDataLoader();
283                if( !dataLoader->IsFileOpen() ) {
284                   if( dataLoader->OpenFile( "UPDATE" ) ) {
285                      AliWarning( Form( "\n\nCan't write trigger for %s\n", det->GetName() ) );
286                   }
287                }
288                dataLoader->Cd();
289                if( gFile && !gFile->IsWritable() ) {
290                   gFile->ReOpen( "UPDATE" );
291                   dataLoader->Cd();
292                }
293                trgdet->Write( "Trigger", TObject::kOverwrite );
294                dataLoader->CloseFile();
295             }
296             else  AliWarning( Form( "Not loader found for %s", det->GetName() ) );
297          }
298       }
299
300       // Check trigger conditions and create the trigger class mask
301       CheckConditions();
302       fInputs.Clear();
303
304       // Clear trigger detectors
305       trgdetArray.SetOwner();
306       trgdetArray.Delete();
307
308       if( (detStr.CompareTo( "ALL" ) != 0) && !detStr.IsNull() ) {
309          AliError( Form("the following detectors were not found: %s",
310                    detStr.Data()));
311          return kFALSE;
312       }
313
314       // Save trigger mask
315       tree->Fill();
316       AliInfo( Form("**************** Central Trigger Class Mask:0x%X", fClassMask ) );
317    } // end event loop
318
319    Reset();
320 //   cout << endl <<  " Print " << endl;
321 //   Print();
322
323    // Write
324    runLoader->WriteTrigger( "OVERWRITE" );
325
326    return kTRUE;
327 }
328
329 //_____________________________________________________________________________
330 Long_t AliCentralTrigger::CheckConditions()
331 {
332    // Check trigger conditions and create the trigger class mask
333
334    Int_t ndes = fDescriptors.GetEntriesFast();
335    for( Int_t i=0; i<ndes; i++ ) {
336       TObjArray* condArray = ((AliTriggerDescriptor*)fDescriptors.At( i ))->GetTriggerConditions();
337       Int_t ncond = condArray->GetEntriesFast();
338       for( Int_t j=0; j<ncond; j++ ) {
339          AliTriggerCondition* cond = (AliTriggerCondition*)condArray->At( j );
340          if( !cond->CheckInputs( fInputs ) ) continue;
341          cond->Trigger( fInputs );
342     //     cond->Print();
343          fClassMask |= cond->GetValue();
344       }
345    }
346    return fClassMask;
347 }
348 //_____________________________________________________________________________
349 TObjArray* AliCentralTrigger::GetResultConditions()
350 {
351    // return only the true conditions
352
353    TObjArray* result = new TObjArray();
354
355    Int_t ndes = fDescriptors.GetEntriesFast();
356    for( Int_t i=0; i<ndes; i++ ) {
357       TObjArray* condArray = ((AliTriggerDescriptor*)fDescriptors.At( i ))->GetTriggerConditions();
358       Int_t ncond = condArray->GetEntriesFast();
359       for( Int_t j=0; j<ncond; j++ ) {
360          AliTriggerCondition* cond = (AliTriggerCondition*)condArray->At( j );
361          if( cond->GetStatus() ) result->AddLast( cond );
362       }
363    }
364
365    return result;
366 }
367
368 //_____________________________________________________________________________
369 void AliCentralTrigger::Print( const Option_t*  ) const
370 {
371    // Print
372    cout << "Central Trigger: " << endl;
373    cout << "  Trigger Class Mask: 0x" << hex << fClassMask << dec << endl;
374
375    Int_t ndes = fDescriptors.GetEntriesFast();
376    for( Int_t i=0; i<ndes; i++ ) {
377       AliTriggerDescriptor* des = (AliTriggerDescriptor*)fDescriptors.At( i );
378       if( des ) des->Print();
379    }
380    cout << endl;
381 }
382
383
384 //////////////////////////////////////////////////////////////////////////////
385 // Helper method
386
387 //_____________________________________________________________________________
388 Bool_t AliCentralTrigger::IsSelected( TString detName, TString& detectors ) const
389 {
390    // check whether detName is contained in detectors
391    // if yes, it is removed from detectors
392
393    // check if all detectors are selected
394    if( (detectors.CompareTo("ALL") == 0 ) ||
395         detectors.BeginsWith("ALL ") ||
396         detectors.EndsWith(" ALL") ||
397         detectors.Contains(" ALL ") ) {
398       detectors = "ALL";
399       return kTRUE;
400    }
401
402    // search for the given detector
403    Bool_t result = kFALSE;
404    if( (detectors.CompareTo( detName ) == 0) ||
405         detectors.BeginsWith( detName+" " ) ||
406         detectors.EndsWith( " "+detName ) ||
407         detectors.Contains( " "+detName+" " ) ) {
408       detectors.ReplaceAll( detName, "" );
409       result = kTRUE;
410    }
411
412    // clean up the detectors string
413    while( detectors.Contains("  ") )  detectors.ReplaceAll( "  ", " " );
414    while( detectors.BeginsWith(" ") ) detectors.Remove( 0, 1 );
415    while( detectors.EndsWith(" ") )   detectors.Remove( detectors.Length()-1, 1 );
416
417    return result;
418 }