CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/L1Trigger/RPCTechnicalTrigger/src/RPCTechnicalTrigger.cc

Go to the documentation of this file.
00001 // $Id: 
00002 
00003 //-----------------------------------------------------------------------------
00004 // Implementation file for class : RPCTechnicalTrigger
00005 //
00006 // 2008-10-15 : Andres Osorio
00007 //-----------------------------------------------------------------------------
00008 
00009 //=============================================================================
00010 // Standard constructor, initializes variables
00011 //=============================================================================
00012 
00013 // Include files
00014 
00015 // local
00016 #include "L1Trigger/RPCTechnicalTrigger/interface/RPCTechnicalTrigger.h"
00017 #include "L1Trigger/RPCTechnicalTrigger/interface/ProcessTestSignal.h"
00018 #include "L1Trigger/RPCTechnicalTrigger/interface/RBCProcessRPCDigis.h"
00019 #include "L1Trigger/RPCTechnicalTrigger/interface/RBCProcessRPCSimDigis.h"
00020 
00021 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00022 #include "FWCore/MessageLogger/interface/MessageDrop.h"
00023 
00024 //=============================================================================
00025 // Standard constructor, initializes variables
00026 //=============================================================================
00027 
00028 RPCTechnicalTrigger::RPCTechnicalTrigger(const edm::ParameterSet& iConfig) {
00029   
00030   //...........................................................................
00031   
00032   std::string configFile  = iConfig.getParameter<std::string>("ConfigFile");
00033   m_verbosity             = iConfig.getUntrackedParameter<int>("Verbosity", 0);
00034   m_rpcDigiLabel          = iConfig.getParameter<edm::InputTag>("RPCDigiLabel");
00035   m_ttBits                = iConfig.getParameter< std::vector<unsigned> >("BitNumbers");
00036   m_ttNames               = iConfig.getParameter< std::vector<std::string> >("BitNames");
00037   m_useEventSetup         = iConfig.getUntrackedParameter<int>("UseEventSetup", 0);
00038   m_useRPCSimLink         = iConfig.getUntrackedParameter<int>("UseRPCSimLink", 0);
00039   m_rpcSimLinkInstance    = iConfig.getParameter<edm::InputTag>("RPCSimLinkInstance");
00040 
00041   edm::FileInPath f1("L1Trigger/RPCTechnicalTrigger/data/" + configFile);
00042   m_configFile = f1.fullPath();
00043 
00044   if ( m_verbosity ) {
00045     LogTrace("RPCTechnicalTrigger")
00046       << m_rpcDigiLabel << '\n'
00047       << std::endl;
00048 
00049     LogTrace("RPCTechnicalTrigger")
00050       << "\nConfiguration file used for UseEventSetup = 0 \n" << m_configFile << '\n'
00051       << std::endl;
00052   }
00053   
00054   //...........................................................................
00055   //... There are three Technical Trigger Units Boards: 1 can handle 2 Wheels
00056   //... n_Wheels sets the number of wheels attached to board with index boardIndex
00057   
00058   m_boardIndex[0] = 1;
00059   m_boardIndex[1] = 2;
00060   m_boardIndex[2] = 3;
00061   
00062   m_nWheels[0]    = 2;
00063   m_nWheels[1]    = 1;
00064   m_nWheels[2]    = 2;
00065   
00066   m_ttu[0] = new TTUEmulator( m_boardIndex[0] , m_nWheels[0] );
00067   m_ttu[1] = new TTUEmulator( m_boardIndex[1] , m_nWheels[1] );
00068   m_ttu[2] = new TTUEmulator( m_boardIndex[2] , m_nWheels[2] );
00069 
00070   //... This is second line that delivers in parallel a second trigger
00071   m_ttuRbcLine[0] = new TTUEmulator( m_boardIndex[0] , m_nWheels[0] );
00072   m_ttuRbcLine[1] = new TTUEmulator( m_boardIndex[1] , m_nWheels[1] );
00073   m_ttuRbcLine[2] = new TTUEmulator( m_boardIndex[2] , m_nWheels[2] );
00074   
00075   m_WheelTtu[-2] = 3;
00076   m_WheelTtu[-1] = 3;
00077   m_WheelTtu[0 ] = 2;
00078   m_WheelTtu[1 ] = 1;
00079   m_WheelTtu[2 ] = 1;
00080   
00081   //...........................................................................
00082   //For the pointing Logic: declare here the first sector of each quadrant
00083   //
00084   m_quadrants.push_back(2);
00085   m_quadrants.push_back(3);
00086   m_quadrants.push_back(4);
00087   m_quadrants.push_back(5);
00088   m_quadrants.push_back(6);
00089   m_quadrants.push_back(7);
00090   m_quadrants.push_back(8);
00091   m_quadrants.push_back(9);
00092   m_quadrants.push_back(10);
00093   m_quadrants.push_back(11);
00094 
00095   //...........................................................................
00096   
00097   m_ievt = 0;
00098   m_cand = 0;
00099   m_maxTtuBoards = 3;
00100   m_maxBits = 5;
00101   m_hasConfig = false;
00102   m_readConfig = NULL;
00103   produces<L1GtTechnicalTriggerRecord>();
00104   
00105 }
00106 
00107 
00108 RPCTechnicalTrigger::~RPCTechnicalTrigger()
00109 {
00110   
00111   LogDebug("RPCTechnicalTrigger") << "RPCTechnicalTrigger: object starts deletion" << std::endl;
00112 
00113   if ( m_hasConfig ) {
00114     
00115     delete m_ttu[0];
00116     delete m_ttu[1];
00117     delete m_ttu[2];
00118     
00119     delete m_ttuRbcLine[0];
00120     delete m_ttuRbcLine[1];
00121     delete m_ttuRbcLine[2];
00122     
00123     if ( m_readConfig )
00124       delete m_readConfig;
00125     
00126   }
00127   
00128   m_WheelTtu.clear();
00129     
00130   LogDebug("RPCTechnicalTrigger") << "RPCTechnicalTrigger: object deleted" << '\n';
00131   
00132 }
00133 
00134 //=============================================================================
00135 void RPCTechnicalTrigger::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
00136 
00137 
00138   bool status(false);
00139   
00140   edm::Handle<RPCDigiCollection> pIn;
00141   
00142   edm::Handle<edm::DetSetVector<RPCDigiSimLink> > simIn;
00143   
00144   std::auto_ptr<L1GtTechnicalTriggerRecord> output(new L1GtTechnicalTriggerRecord());
00145   
00146   if ( m_useRPCSimLink == 0 ) {
00147 
00148     iEvent.getByLabel(m_rpcDigiLabel, pIn);
00149     if ( ! pIn.isValid() ) {
00150       edm::LogError("RPCTechnicalTrigger") << "can't find RPCDigiCollection with label: " 
00151                                            << m_rpcDigiLabel << '\n';
00152       iEvent.put(output);
00153       return;
00154     }
00155     m_signal  = dynamic_cast<ProcessInputSignal*>(new RBCProcessRPCDigis( m_rpcGeometry, pIn ));
00156     
00157   } else {
00158     
00159     iEvent.getByLabel("simMuonRPCDigis", "RPCDigiSimLink", simIn);
00160     
00161     if ( ! simIn.isValid() ) {
00162       edm::LogError("RPCTechnicalTrigger") << "can't find RPCDigiCollection with label: " 
00163                                            << m_rpcDigiLabel << '\n';
00164       iEvent.put(output);
00165       return;
00166     }
00167     m_signal  = dynamic_cast<ProcessInputSignal*>(new RBCProcessRPCSimDigis( m_rpcGeometry, simIn ));
00168   }
00169   
00170   LogDebug("RPCTechnicalTrigger") << "signal object created" << '\n';
00171   
00172   if ( ! m_hasConfig ) {
00173     edm::LogError("RPCTechnicalTrigger") << "cannot read hardware configuration \n";
00174     iEvent.put(output);
00175     return;
00176   }
00177   
00178   status = m_signal->next();
00179   
00180   if ( !status)  { 
00181     delete m_signal;
00182     iEvent.put(output);
00183     return;
00184   }
00185   
00186   m_input = m_signal->retrievedata();
00187   
00188   std::vector<L1GtTechnicalTrigger> ttVec( m_ttBits.size() );
00189   
00190   //. distribute data to different TTU emulator instances and process it
00191   
00192   m_triggerbits.reset();
00193   
00194   int indx(0);
00195   
00196   std::vector<TTUEmulator::TriggerResponse*>::const_iterator outItr;
00197   
00198   for(int k=0; k < m_maxTtuBoards; ++k) {
00199     
00200     indx=k*2;
00201     
00202     m_ttu[k]->processTtu( m_input );
00203     
00204     //work out Pointing Logic to Tracker
00205     for( m_firstSector = m_quadrants.begin(); m_firstSector != m_quadrants.end(); ++m_firstSector)
00206       m_ttuRbcLine[k]->processTtu( m_input , (*m_firstSector) );
00207     
00208     //...for trigger 1
00209     for( outItr  = m_ttu[k]->m_triggerBxVec.begin(); outItr != m_ttu[k]->m_triggerBxVec.end(); ++outItr )
00210       m_serializedInfoLine1.push_back( new TTUResults( k, (*outItr)->m_bx, (*outItr)->m_trigger[0], (*outItr)->m_trigger[1] ) );
00211     m_ttu[k]->clearTriggerResponse();
00212     
00213     //...for trigger 2
00214     for( outItr  = m_ttuRbcLine[k]->m_triggerBxVec.begin(); outItr != m_ttuRbcLine[k]->m_triggerBxVec.end(); ++outItr )
00215       m_serializedInfoLine2.push_back( new TTUResults( k, 
00216                                                        (*outItr)->m_bx, 
00217                                                        (*outItr)->m_trigger[0], 
00218                                                        (*outItr)->m_trigger[1], 
00219                                                        (*outItr)->m_wedge ) );
00220     
00221     m_ttuRbcLine[k]->clearTriggerResponse();
00222     
00223   }
00224   
00225   //.. write results to technical trigger bits
00226   int bx(0);
00227   int infoSize(0);
00228   
00229   infoSize = m_serializedInfoLine1.size();
00230 
00231   std::vector<RPCTechnicalTrigger::TTUResults*>::const_iterator ttuItr;
00232   
00233   std::sort( m_serializedInfoLine1.begin(), m_serializedInfoLine1.end(), sortByBx() );
00234   
00235   for( ttuItr = m_serializedInfoLine1.begin(); ttuItr != m_serializedInfoLine1.end(); ++ttuItr ) {
00236     if ( m_verbosity && abs( (*ttuItr)->m_bx ) <= 1 ) 
00237       std::cout << "RPCTechnicalTrigger> " 
00238                 << (*ttuItr)->m_ttuidx << '\t'
00239                 << (*ttuItr)->m_bx << '\t'
00240                 << (*ttuItr)->m_trigWheel1 << '\t'
00241                 << (*ttuItr)->m_trigWheel2 << '\n';
00242   }
00243   
00244   bool has_bx0 = false;
00245   
00246   for(int k = 0; k < infoSize; k+=m_maxTtuBoards) {
00247     
00248     bx = m_serializedInfoLine1[k]->m_bx;
00249     
00250     if ( bx == 0 ) {
00251       
00252       m_triggerbits.set(0, m_serializedInfoLine1[k]->m_trigWheel2);
00253       m_triggerbits.set(1, m_serializedInfoLine1[k]->m_trigWheel1);
00254       m_triggerbits.set(2, m_serializedInfoLine1[k+1]->m_trigWheel1);
00255       m_triggerbits.set(3, m_serializedInfoLine1[k+2]->m_trigWheel1);
00256       m_triggerbits.set(4, m_serializedInfoLine1[k+2]->m_trigWheel2);
00257       
00258       bool five_wheels_OR = m_triggerbits.any();
00259       
00260       ttVec.at(0)=L1GtTechnicalTrigger(m_ttNames.at(0), m_ttBits.at(0), bx, five_wheels_OR ) ;   // bit 24 = Or 5 wheels in TTU mode
00261       ttVec.at(2)=L1GtTechnicalTrigger(m_ttNames.at(2), m_ttBits.at(2), bx, m_triggerbits[0] ) ; // bit 26 
00262       ttVec.at(3)=L1GtTechnicalTrigger(m_ttNames.at(3), m_ttBits.at(3), bx, m_triggerbits[1] ) ; // bit 27 
00263       ttVec.at(4)=L1GtTechnicalTrigger(m_ttNames.at(4), m_ttBits.at(4), bx, m_triggerbits[2] ) ; // bit 28 
00264       ttVec.at(5)=L1GtTechnicalTrigger(m_ttNames.at(5), m_ttBits.at(5), bx, m_triggerbits[3] ) ; // bit 29
00265       ttVec.at(6)=L1GtTechnicalTrigger(m_ttNames.at(6), m_ttBits.at(6), bx, m_triggerbits[4] ) ; // bit 30
00266       
00267       m_triggerbits.reset();
00268       
00269       has_bx0 = true;
00270       
00271       break;
00272       
00273     } else continue;
00274     
00275   }
00276   
00277   infoSize = m_serializedInfoLine2.size();
00278   
00279   std::sort( m_serializedInfoLine2.begin(), m_serializedInfoLine2.end(), sortByBx() );
00280   
00281   for( ttuItr = m_serializedInfoLine2.begin(); ttuItr != m_serializedInfoLine2.end(); ++ttuItr ) {
00282     if ( m_verbosity && abs ( (*ttuItr)->m_bx ) <= 1 )
00283       std::cout << "RPCTechnicalTrigger> " 
00284                 << (*ttuItr)->m_ttuidx << '\t'
00285                 << (*ttuItr)->m_bx << '\t'
00286                 << (*ttuItr)->m_trigWheel1 << '\t'
00287                 << (*ttuItr)->m_trigWheel2 << '\t'
00288                 << (*ttuItr)->m_wedge << '\n';
00289   }
00290   
00291   infoSize = convertToMap( m_serializedInfoLine2 );
00292   
00293   std::bitset<8> triggerCoincidence;
00294   triggerCoincidence.reset();
00295   
00296   // searchCoincidence( W-2 , W0 )
00297   bool result = searchCoincidence( -2, 0 );
00298   triggerCoincidence.set(0, result );
00299   
00300   // searchCoincidence( W-2 , W+1 )
00301   result = searchCoincidence( -2, 1 );
00302   triggerCoincidence.set(1, result );
00303   
00304   // searchCoincidence( W-1 , W0  )
00305   result = searchCoincidence( -1, 0 );
00306   triggerCoincidence.set(2, result );
00307   
00308   // searchCoincidence( W-1 , W+1 )
00309   result = searchCoincidence( -1, 1 );
00310   triggerCoincidence.set(3, result );
00311   
00312   // searchCoincidence( W-1 , W+2 )
00313   result = searchCoincidence( -1, 2 );
00314   triggerCoincidence.set(4, result );
00315   
00316   // searchCoincidence( W0  , W0  )
00317   result = searchCoincidence( 0 , 0 );
00318   triggerCoincidence.set(5, result );
00319   
00320   // searchCoincidence( W+1 , W0  )
00321   result = searchCoincidence( 1, 0 );
00322   triggerCoincidence.set(6, result );
00323   
00324   // searchCoincidence( W+2 , W0  ) 
00325   result = searchCoincidence( 2, 0 );
00326   triggerCoincidence.set(7, result );
00327   
00328   bool five_wheels_OR = triggerCoincidence.any();
00329 
00330   if ( m_verbosity ) std::cout << "RPCTechnicalTrigger> pointing trigger: " << five_wheels_OR << '\n';
00331   
00332   ttVec.at(1)=L1GtTechnicalTrigger(m_ttNames.at(1), m_ttBits.at(1), bx, five_wheels_OR ) ; // bit 25 = Or 5 wheels in RBC mode
00333   
00334   triggerCoincidence.reset();
00335   
00336   //...check that data appeared at bx=0
00337   
00338   if ( ! has_bx0 ) {
00339     iEvent.put(output);
00340     status = Reset();
00341     ++m_ievt;
00342     LogDebug("RPCTechnicalTrigger") << "RPCTechnicalTrigger> end of event loop" << std::endl;
00343     return;
00344     
00345   }
00346   
00347   output->setGtTechnicalTrigger(ttVec);    
00348   iEvent.put(output);
00349   
00350   //.... all done
00351   
00352   status = Reset();
00353   ++m_ievt;
00354   LogDebug("RPCTechnicalTrigger") << "RPCTechnicalTrigger> end of event loop" << std::endl;
00355   
00356 }
00357 
00358 bool RPCTechnicalTrigger::Reset()
00359 {
00360   
00361   m_input->clear();
00362   m_triggerbits.reset();
00363   std::vector<TTUResults*>::iterator itrRes;
00364   
00365   for( itrRes=m_serializedInfoLine1.begin(); itrRes!=m_serializedInfoLine1.end(); ++itrRes)
00366     delete (*itrRes);
00367   
00368   for( itrRes=m_serializedInfoLine2.begin(); itrRes!=m_serializedInfoLine2.end(); ++itrRes)
00369     delete (*itrRes);
00370   
00371   m_serializedInfoLine1.clear();
00372   m_serializedInfoLine2.clear();
00373   m_ttuResultsByQuadrant.clear();
00374   
00375   delete m_signal; 
00376   
00377   return true;
00378   
00379 }
00380 
00381 // ------------ method called once each job just before starting event loop  ------------
00382 void RPCTechnicalTrigger::beginRun(edm::Run& iRun, const edm::EventSetup& evtSetup)
00383 {
00384   
00385   bool status(false);
00386   
00387   LogDebug("RPCTechnicalTrigger") << "RPCTechnicalTrigger::beginRun> starts" << std::endl;
00388   
00389   //.   Set up RPC geometry
00390   
00391   evtSetup.get<MuonGeometryRecord>().get( m_rpcGeometry );
00392   
00393   //..  Get Board Specifications (hardware configuration)
00394   
00395   if ( m_useEventSetup >= 1 ) {
00396     
00397     edm::ESHandle<RBCBoardSpecs> pRBCSpecs;
00398     evtSetup.get<RBCBoardSpecsRcd>().get(pRBCSpecs);
00399 
00400     edm::ESHandle<TTUBoardSpecs> pTTUSpecs;
00401     evtSetup.get<TTUBoardSpecsRcd>().get(pTTUSpecs);
00402     
00403     if ( !pRBCSpecs.isValid() ||  !pTTUSpecs.isValid() ) {
00404       edm::LogError("RPCTechnicalTrigger") << "can't find RBC/TTU BoardSpecsRcd" << '\n';
00405       m_hasConfig = false;
00406     }
00407     else  {
00408       m_rbcspecs = pRBCSpecs.product();
00409       m_ttuspecs = pTTUSpecs.product();
00410       m_hasConfig = true;
00411     }
00412     
00413   } else {
00414     
00415     // read hardware configuration from file
00416     m_readConfig = new TTUConfigurator( m_configFile );
00417     
00418     if ( m_readConfig->m_hasConfig ) {
00419       m_readConfig->process();
00420       m_rbcspecs = m_readConfig->getRbcSpecs();
00421       m_ttuspecs = m_readConfig->getTtuSpecs();
00422       m_hasConfig = true;
00423     }
00424     
00425     else m_hasConfig = false;
00426     
00427   }
00428   
00429   if ( m_hasConfig ) {
00430     
00431     //... Initialize all
00432     
00433     for (int k=0; k < m_maxTtuBoards; ++k ) {
00434 
00435       m_ttu[k]->SetLineId ( 1 );
00436       m_ttuRbcLine[k]->SetLineId( 2 );
00437       
00438       m_ttu[k]->setSpecifications( m_ttuspecs, m_rbcspecs );
00439       m_ttuRbcLine[k]->setSpecifications( m_ttuspecs, m_rbcspecs );
00440       
00441       status = m_ttu[k]->initialise();
00442       status = m_ttuRbcLine[k]->initialise();
00443       
00444       
00445     }
00446   
00447   }
00448     
00449 }
00450 
00451 //
00452 int RPCTechnicalTrigger::convertToMap( const std::vector<TTUResults*> & ttuResults )
00453 {
00454   
00455   std::vector<TTUResults*>::const_iterator itr = ttuResults.begin();
00456   
00457   while ( itr != ttuResults.end() ) {
00458     
00459     if ( (*itr)->m_bx != 0 ) {
00460       ++itr;
00461       continue;
00462     }
00463     
00464     int key(0);
00465     key = 1000 * ( (*itr)->m_ttuidx + 1 ) + 1*(*itr)->m_wedge;
00466     m_ttuResultsByQuadrant[ key ] = (*itr);
00467     ++itr;
00468     
00469   }
00470   
00471   return m_ttuResultsByQuadrant.size();
00472     
00473 }
00474 
00475 //...RBC pointing logic to tracker bit 25: hardwired
00476 bool RPCTechnicalTrigger::searchCoincidence( int wheel1, int wheel2 )
00477 {
00478   
00479   std::map<int, TTUResults*>::iterator itr;
00480   bool topRight(false);
00481   bool botLeft(false);
00482   
00483   int indxW1 = m_WheelTtu[wheel1];
00484   int indxW2 = m_WheelTtu[wheel2];
00485   
00486   int k(0);
00487   int key(0);
00488   bool finalTrigger(false);
00489   int maxTopQuadrants = 4;
00490   
00491   //work out Pointing Logic to Tracker
00492   
00493   for( m_firstSector = m_quadrants.begin(); m_firstSector != m_quadrants.end(); ++m_firstSector) {
00494     
00495     key = 1000 * ( indxW1 ) + (*m_firstSector);
00496     
00497     itr = m_ttuResultsByQuadrant.find( key );
00498     if ( itr != m_ttuResultsByQuadrant.end() )
00499       topRight  =  (*itr).second->getTriggerForWheel(wheel1);
00500 
00501     //std::cout << "W1: " << wheel1 << " " << "sec: " << (*m_firstSector) << " dec: " << topRight << '\n';
00502     
00503     key = 1000 * ( indxW2 ) + (*m_firstSector) + 5;
00504     
00505     itr = m_ttuResultsByQuadrant.find( key );
00506     
00507     if ( itr != m_ttuResultsByQuadrant.end() )
00508       botLeft   = (*itr).second->getTriggerForWheel(wheel2);
00509     
00510     //std::cout << "W2: " << wheel2 << " " << "sec: " << (*m_firstSector) + 5 << " dec: " << botLeft << '\n';
00511     
00512     finalTrigger |= ( topRight && botLeft );
00513     
00514     ++k;
00515     
00516     if ( k > maxTopQuadrants)
00517       break;
00518         
00519   }
00520   
00521   //Try the opposite now
00522 
00523   k=0;
00524   
00525   for( m_firstSector = m_quadrants.begin(); m_firstSector != m_quadrants.end(); ++m_firstSector) {
00526     
00527     key = 1000 * ( indxW2 ) + (*m_firstSector);
00528     
00529     itr = m_ttuResultsByQuadrant.find( key );
00530     if ( itr != m_ttuResultsByQuadrant.end() )
00531       topRight  =  (*itr).second->getTriggerForWheel(wheel1);
00532 
00533     //std::cout << "W1: " << wheel1 << " " << "sec: " << (*m_firstSector) << " dec: " << topRight << '\n';
00534     
00535     key = 1000 * ( indxW1 ) + (*m_firstSector) + 5;
00536     
00537     itr = m_ttuResultsByQuadrant.find( key );
00538     
00539     if ( itr != m_ttuResultsByQuadrant.end() )
00540       botLeft   = (*itr).second->getTriggerForWheel(wheel2);
00541     
00542     //std::cout << "W2: " << wheel2 << " " << "sec: " << (*m_firstSector) + 5 << " dec: " << botLeft << '\n';
00543     
00544     finalTrigger |= ( topRight && botLeft );
00545     
00546     ++k;
00547     
00548     if ( k > maxTopQuadrants)
00549       break;
00550         
00551   }
00552   
00553   return finalTrigger;
00554   
00555 }
00556 
00557 // ------------ method called once each job just after ending the event loop  ------------
00558 
00559 void RPCTechnicalTrigger::endJob() 
00560 {
00561   
00562   LogDebug("RPCTechnicalTrigger") << "RPCTechnicalTrigger::endJob>" << std::endl;
00563   
00564 }
00565 
00566 void RPCTechnicalTrigger::printinfo()
00567 {
00568   
00569   LogDebug("RPCTechnicalTrigger") << "RPCTechnicalTrigger::Printing TTU emulators info>" << std::endl;
00570   
00571   for (int k=0; k < m_maxTtuBoards; ++k ) {
00572     m_ttu[k]->printinfo();
00573     m_ttuRbcLine[k]->printinfo();
00574   }
00575   
00576     
00577 }
00578 
00579 
00580 //define this as a plug-in
00581 DEFINE_FWK_MODULE(RPCTechnicalTrigger);