1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 ///////////////////////////////////////////////////////////////////////////////
20 // This class for running the Central Trigger Processor //
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 //
28 // Check the condition classes //
29 // Create the class mask //
33 ///////////////////////////////////////////////////////////////////////////////
36 #include <TObjString.h>
37 #include <TObjArray.h>
38 #include <TStopwatch.h>
44 #include "AliRunLoader.h"
45 #include "AliModule.h"
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"
56 ClassImp( AliCentralTrigger )
58 //_____________________________________________________________________________
59 AliCentralTrigger::AliCentralTrigger() :
63 // Default constructor
64 // LoadDescriptor("Pb-Pb");
67 //_____________________________________________________________________________
68 AliCentralTrigger::AliCentralTrigger( TString & descriptor ) :
72 // Default constructor
73 LoadDescriptor( descriptor );
76 //_____________________________________________________________________________
77 AliCentralTrigger::AliCentralTrigger( const AliCentralTrigger& ctp ):
79 fClassMask( ctp.fClassMask )
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 )) ) );
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 )) ) );
94 //_____________________________________________________________________________
95 AliCentralTrigger::~AliCentralTrigger()
101 //_____________________________________________________________________________
102 void AliCentralTrigger::DeleteDescriptors()
106 fDescriptors.SetOwner();
107 fDescriptors.Delete();
110 //_____________________________________________________________________________
111 void AliCentralTrigger::Reset()
113 // Reset Class Mask and conditions
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 );
126 //_____________________________________________________________________________
127 void AliCentralTrigger::MakeBranch( TString name, TTree * tree )
129 // Make a branch to store only trigger class mask event by event
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 );
140 AliDebug( 1, "Got Branch from Tree" );
141 branch->SetAddress( &(this->fClassMask) );
146 //_____________________________________________________________________________
147 Bool_t AliCentralTrigger::LoadDescriptor( TString & descriptor )
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"
153 // Delete any descriptor
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() );
163 fDescriptors.AddLast( des );
165 AliWarning( Form( "Descriptor (%s) not found", val->String().Data() ) );
167 Bool_t desfound = kTRUE;
168 if( fDescriptors.GetEntriesFast() == 0 ) desfound = kFALSE;
175 //_____________________________________________________________________________
176 TString AliCentralTrigger::GetDetectors()
178 // return TString with the detectors to be trigger
179 // merging detectors from all descriptors
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() );
198 //_____________________________________________________________________________
199 UChar_t AliCentralTrigger::GetClusterMask()
201 // Return the detector cluster mask following
202 // table 4.3 pag 60, TDR Trigger and DAQ
204 TString detStr = GetDetectors();
205 TObjArray* det = detStr.Tokenize(" ");
206 Int_t ndet = det->GetEntriesFast();
209 if( ndet >= 8 ) { // All detectors, should be 9 but CRT is not implemented yet
214 if( ndet >= 7 && !detStr.Contains("MUON") ) { // Central Barrel, All but MUON
219 if( detStr.Contains("MUON") && detStr.Contains("START") ) { // MUON arm
224 return idmask; // 0 something else!!!
226 //_____________________________________________________________________________
227 Bool_t AliCentralTrigger::RunTrigger( AliRunLoader* runLoader )
231 if( fDescriptors.GetEntriesFast() == 0 ) {
232 AliError( "not trigger descriptor loaded, skipping trigger" );
236 TTree *tree = runLoader->TreeCT();
238 AliError( "not folder with trigger loaded, skipping trigger" );
242 TStopwatch stopwatch;
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();
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) ) {
261 AliInfo( Form("Triggering from digits for %s", det->GetName() ) );
262 AliTriggerDetector* trgdet = det->CreateTriggerDetector();
263 trgdet->CreateInputs();
264 TStopwatch stopwatchDet;
265 stopwatchDet.Start();
267 AliInfo( Form("Execution time for %s: R:%.2fs C:%.2fs",
268 det->GetName(), stopwatchDet.RealTime(), stopwatchDet.CpuTime() ) );
271 TObjArray* detInp = trgdet->GetInputs();
272 for( Int_t i=0; i<detInp->GetEntriesFast(); i++ ) {
273 fInputs.AddLast( detInp->At(i) );
275 trgdetArray.AddLast( trgdet );
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 );
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() ) );
289 if( gFile && !gFile->IsWritable() ) {
290 gFile->ReOpen( "UPDATE" );
293 trgdet->Write( "Trigger", TObject::kOverwrite );
294 dataLoader->CloseFile();
296 else AliWarning( Form( "Not loader found for %s", det->GetName() ) );
300 // Check trigger conditions and create the trigger class mask
304 // Clear trigger detectors
305 trgdetArray.SetOwner();
306 trgdetArray.Delete();
308 if( (detStr.CompareTo( "ALL" ) != 0) && !detStr.IsNull() ) {
309 AliError( Form("the following detectors were not found: %s",
316 AliInfo( Form("**************** Central Trigger Class Mask:0x%X", fClassMask ) );
320 // cout << endl << " Print " << endl;
324 runLoader->WriteTrigger( "OVERWRITE" );
329 //_____________________________________________________________________________
330 Long_t AliCentralTrigger::CheckConditions()
332 // Check trigger conditions and create the trigger class mask
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 );
343 fClassMask |= cond->GetValue();
348 //_____________________________________________________________________________
349 TObjArray* AliCentralTrigger::GetResultConditions()
351 // return only the true conditions
353 TObjArray* result = new TObjArray();
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 );
368 //_____________________________________________________________________________
369 void AliCentralTrigger::Print( const Option_t* ) const
372 cout << "Central Trigger: " << endl;
373 cout << " Trigger Class Mask: 0x" << hex << fClassMask << dec << endl;
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();
384 //////////////////////////////////////////////////////////////////////////////
387 //_____________________________________________________________________________
388 Bool_t AliCentralTrigger::IsSelected( TString detName, TString& detectors ) const
390 // check whether detName is contained in detectors
391 // if yes, it is removed from detectors
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 ") ) {
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, "" );
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 );