LoKi's Array Particle Functions



ACHI2DOCA

It is an alias for ADOCACHI2

ACHI2V

It is an alias for AVCHI2

ACHI2VX

It is an alias for AVCHI2

ACHILD, C++ type LoKi::AParticles::ChildFun

The function which evaluates another function using the child particle. Get the DLL(K-pi) for the first particle:

 1000// Get the DLL(K-pi) for the first particle 
 1010AFun fun = ACHILD ( PIDK , 1 ) ;
 1020
 1030const LHCb::Particle::ConstVector& combination = ... ;
 1040
 1050const double kaonDLL = fun ( combination ) ;

ACHILDCUT, C++ type LoKi::AParticles::ChildCut

The predicate which evaluates another predicate using the child particle. Check the DLL(K-pi) for the first particle:

 1000// check the DLL(K-pi) for the first particle 
 1010ACut cut = ACHILDCUT ( PIDK > 0 , 1 ) ;
 1020
 1030const LHCb::Particle::ConstVector& combination = ... ;
 1040
 1050const book goodKaon = cut ( combination ) ;

ACHILDFUN

It is an alias for ACHILD

ACOUNTER, C++ type LoKi::Monitoring::Counter<LoKi::ATypes::Combination>

The special monitoring predicate, which could be used for monitoring of other predicates. The generic Gaudi counter StatEntity is used for monitoring:

 1000// the predicate to be monitored
 1010ACut cut = ... ;
 1020
 1030// the generic counter:
 1040StatEntity* counter = ... ;
 1050
 1060// create the monitoring predicate:
 1070ACut mon = ACOUNTER ( cut , counter ) ;
 1080
 1090// use the counter:
 1100for ( .... ) 
 1110{
 1120   ...
 1130   const LHCb::Particle::ConstVector& comb = ... ;
 1140   const bool result = mon ( comb ) ; ///< use the monitored predicate
 1150   ...
 1160}
 1170...
 1180// inspect the counter:
 1190info () 
 1200   << " Monitoring results : "                       << endreq 
 1210   << " NEntries:  #" << counter->entries ()         << endreq 
 1220   << " Efficiency:(" << counter->eff     ()   
 1230   << "+="            << counter->effErr  () << ")%" << endreq ;
 1240 
The alternative and recommended way for creation of of monitored predicates:
 1000// the predicate to be monitored 
 1010ACut cut = ... ;
 1020
 1030// the monitored predicate 
 1040ACut mon = monitor ( cut , counter ( "SomeEffCounter" ) ) ;
 1050     

The substitution of the predicate by the monitored predicate could be done "on-flight" without the disturbing of the actual processing:

 1000// The predicate to be monitored
 1010ACut cut = ... ;
 1020     
 1030if ( monitoring ) 
 1040 {
 1050    cut = monitor ( cut , counter ( "Efficiency1" ) ) ;
 1060 }
 1070    
Attention:
  • The string representation of the object is delegated to the underlying predicate, therefore the object is NOT reconstructable from its string representations. It is done on-purpose to avoid the disturbing of unique function IDs.

ACUTCHILD

It is an alias for ACHILDCUT

ACUTDOCA, C++ type LoKi::AParticles::MaxDOCACut

The predicate which checks the maximal distance of the closest approach for all two-particle subcombinations, using IGeomDispCalculator tool with respect some threshold:

 1000IGeomDispCalculator* doca = ... ;
 1010ACut ok = ACUTDOCA ( doca , 0.2 * mm ) ;
 1020     
 1030const LHCb::Particle::ConstVector& combination = ... ;
 1040     
 1050const bool good = ok ( combination ) ;
 1060     
It should be faster than AMAXDOCA<0.2*mm

ACUTDOCACHI2, C++ type LoKi::AParticles::MaxDOCAChi2Cut

Simple predicate which check if all two-particle subcombinations have a -distance of the closest approach below some threshold, using IGeomDispCalculator tool:

 1000IGeomDispCalculator* doca = ... ;
 1010ACut cut = ACUTDOCACHI2 ( doca , 3*3 ) ;
 1020     
 1030const LHCb::Particle::ConstVector& combination = ... ;
 1040     
 1050const bool good = cut ( combination ) ;
 1060    
It should be more efficient than ADOCACHI2<9

ADOCACHI2, C++ type LoKi::AParticles::MaxDOCAChi2

The evaluator of the maximal -distance of the closest approach for all two-particle subcombinations, using IGeomDispCalculator tool:

 1000IGeomDispCalculator* doca = ... ;
 1010AFun docaMAX = ADOCACHI2 ( doca ) ;
 1020     
 1030const LHCb::Particle::ConstVector& combination = ... ;
 1040     
 1050const double maxDOCA = fun ( combination ) ;
 1060     

ADOCACHI2CUT

It is an alias for ACUTDOCACHI2

ADOCACUT

It is an alias for ACUTDOCA

ADOCAMAX

It is an alias for AMAXDOCA

AETA, the instance of LoKi::AParticles::Eta(0)

Simple evaluator of the pseudorapidity, ,for the combination

 1000// get the combination:
 1010const LHCb::Paricle::ConstVector& comb = ... ;
 1020
 1030// use the functor 
 1040const double eta = AETA ( comb ) ;
See AETAP

AETA0, the instance of LoKi::AParticles::Eta(0)

Simple evaluator of the pseudorapidity, ,for the combination

 1000// get the combination:
 1010const LHCb::Paricle::ConstVector& comb = ... ;
 1020
 1030// use the functor 
 1040const double eta = AETA0 ( comb ) ;
See AETAP

AETAP, C++ type LoKi::AParticles::Eta

Simple evaluator of the pseudorapidity,, for the various subcombinations of particles:

 1000// get the combination:
 1010const LHCb::Paricle::ConstVector& comb = ... ;
 1020
 1030// pseudorapidity of 1st and 4th particles:
 1040AFun eta14 = AETAP(1,4) ;
 1050
 1060// use the functor 
 1070const double e14 = eta14 ( comb ) ;

AFUNCHILD

It is an alias for ACHILD

AHASCHILD, C++ type LoKi::AParticles::HasChild

Check the existence of daughter particle, which satisfy certain criteria.

Check the existnce of "badKaon" in the combination:

 1000// check the existence of "badKaon" in the combination:
 1010ACut cut = AHASCHILD ( "K+" = ABSID && PIDK < -5 ) ;
 1020
 1030const LHCb::Particle::ConstVector& combination = ... ;
 1040
 1050const bool hasBadKaon = cut ( combination ) ;

Check the presence of photons, including the photons, used for recontruction of daughter particles:

 1000// Check the presence of photons, including teh photons,
 1010// used for recontruction of daughter particles:
 1020ACut cut = AHASCHILD ( INTREE ( "gamma" == ID )  ) ;
 1030     
 1040const LHCb::Particle::ConstVector& combination = ... ;
 1050
 1060const bool hasGamma = cut ( combination ) ;

ALV0, C++ type LoKi::AParticles::DecayAngle

The function evaluates the cosine of the angle between daughter's momentum and the flight direction of the whole combination in the rest system of the combination. For 2-body decays it is just a polarization angle of mother particle.

 1000// get the combination:
 1010const LHCb::Paricle::ConstVector& comb = ... ;
 1020
 1030// make functor:
 1040AFun lv01 = ALV0 ( 1 ) ;
 1050
 1060// use the fucntor
 1070const double  value = lv01 ( comb ) ; 
The name comes from KAL language by H.Albrecht used in ARGUS collaboration.

ALV01, the instance of LoKi::AParticles::DecayAngle(1)

The function evaluates the cosine of the angle between the momentum of the first daughter particle and the flight direction of the whole combination in the rest system of the combination. For 2-body decays it is just a polarization angle.

 1000// get the combination:
 1010const LHCb::Paricle::ConstVector& comb = ... ;
 1020
 1030// use the fucntor
 1040const double  value = ALV01 ( comb ) ; 
The name comes from KAL language by H.Albrecht used in ARGUS collaboration.

ALV02, the instance of LoKi::AParticles::DecayAngle(2)

The function evaluates the cosine of the angle between the momentum of the second daughter particle and the flight direction of the whole combination in the rest system of the combination. For 2-body decays it is just a polarization angle.

 1000// get the combination:
 1010const LHCb::Paricle::ConstVector& comb = ... ;
 1020
 1030// use the fucntor
 1040const double  value = ALV02 ( comb ) ; 
The name comes from KAL language by H.Albrecht used in ARGUS collaboration.

ALV03, the instance of LoKi::AParticles::DecayAngle(3)

The function evaluates the cosine of the angle between the momentum of the third daughter particle and the flight direction of the whole combination in the rest system of the combination.

 1000// get the combination:
 1010const LHCb::Paricle::ConstVector& comb = ... ;
 1020
 1030// use the fucntor
 1040const double  value = ALV03 ( comb ) ; 
The name comes from KAL language by H.Albrecht used in ARGUS collaboration.

ALV04, the instance of LoKi::AParticles::DecayAngle(4)

The function evaluates the cosine of the angle between the momentum of the fourth daughter particle and the flight direction of the whole combination in the rest system of the combination.

 1000// get the combination:
 1010const LHCb::Paricle::ConstVector& comb = ... ;
 1020
 1030// use the fucntor
 1040const double  value = ALV04 ( comb ) ; 
The name comes from KAL language by H.Albrecht used in ARGUS collaboration.

AM, the instance of LoKi::AParticles::InvariantMass(0)

Simple evaluator of the invariant mass for the whole combination

 1000// get the combination:
 1010const LHCb::Paricle::ConstVector& comb = ... ;
 1020
 1030// use the functor AM      
 1040const double mass = AM ( comb ) ;
See AMASS

AM0, the instance of LoKi::AParticles::InvariantMass(0)

Simple evaluator of the invariant mass for the whole combination

 1000// get the combination:
 1010const LHCb::Paricle::ConstVector& comb = ... ;
 1020
 1030// use the functor AM0      
 1040const double mass = AM0 ( comb ) ;
See AMASS

AM1, the instance of LoKi::AParticles::InvariantMass(1)

Simple evaluator of the invariant mass of the first daughter particle:

 1000// get the combination:
 1010const LHCb::Paricle::ConstVector& comb = ... ;
 1020
 1030// use the functor AM1     
 1040const double m1 = AM1 ( comb ) ;
See AMASS

AM12, the instance of LoKi::AParticles::InvariantMass(1,2)

Simple evaluator of the invariant mass for the 1st and 2nd particles

 1000// get the combination:
 1010const LHCb::Paricle::ConstVector& comb = ... ;
 1020
 1030// use the functor AM12      
 1040const double m12 = AM12 ( comb ) ;
See AMASS

AM13, the instance of LoKi::AParticles::InvariantMass(1,3)

Simple evaluator of the invariant mass for the 1st and 3rd particles

 1000// get the combination:
 1010const LHCb::Paricle::ConstVector& comb = ... ;
 1020
 1030// use the functor AM13      
 1040const double m13 = AM13 ( comb ) ;
See AMASS

AM14, the instance of LoKi::AParticles::InvariantMass(1,4)

Simple evaluator of the invariant mass for the 1st and 3rd particles

 1000// get the combination:
 1010const LHCb::Paricle::ConstVector& comb = ... ;
 1020
 1030// use the functor AM14      
 1040const double m14 = AM14 ( comb ) ;
See AMASS

AM2, the instance of LoKi::AParticles::InvariantMass(2)

Simple evaluator of the invariant mass of the second daughter particle:

 1000// get the combination:
 1010const LHCb::Paricle::ConstVector& comb = ... ;
 1020
 1030// use the functor AM2     
 1040const double m2 = AM2 ( comb ) ;
See AMASS

AM23, the instance of LoKi::AParticles::InvariantMass(2,3)

Simple evaluator of the invariant mass for the 2nd and 3rd particles

 1000// get the combination:
 1010const LHCb::Paricle::ConstVector& comb = ... ;
 1020
 1030// use the functor AM23      
 1040const double m23 = AM23 ( comb ) ;
See AMASS

AM24, the instance of LoKi::AParticles::InvariantMass(2,4)

Simple evaluator of the invariant mass for the 2nd and 4th particles

 1000// get the combination:
 1010const LHCb::Paricle::ConstVector& comb = ... ;
 1020
 1030// use the functor AM24      
 1040const double m24 = AM24 ( comb ) ;
See AMASS

AM3, the instance of LoKi::AParticles::InvariantMass(3)

Simple evaluator of the invariant mass of the third daughter particle:

 1000// get the combination:
 1010const LHCb::Paricle::ConstVector& comb = ... ;
 1020
 1030// use the functor AM3     
 1040const double m3 = AM3 ( comb ) ;
See AMASS

AM34, the instance of LoKi::AParticles::InvariantMass(3,4)

Simple evaluator of the invariant mass for the 3rd and 4th particles

 1000// get the combination:
 1010const LHCb::Paricle::ConstVector& comb = ... ;
 1020
 1030// use the functor AM34      
 1040const double m34 = AM34 ( comb ) ;
See AMASS

AM4, the instance of LoKi::AParticles::InvariantMass(4)

Simple evaluator of the invariant mass of the fourth daughter particle:

 1000// get the combination:
 1010const LHCb::Paricle::ConstVector& comb = ... ;
 1020
 1030// use the functor AM4     
 1040const double m1 = AM4 ( comb ) ;
See AMASS

AMASS, C++ type LoKi::AParticles::InvariantMass

Simple evaluator of the invariant mass for the subcombination from the combination itself:

 1000// get the combination:
 1010const LHCb::Paricle::ConstVector& comb = ... ;
 1020
 1030// the functor which evaluates the invariant mass of 1st and 2nd particles:
 1040const AFun M12 = AMASS(1,2) ;
 1050
 1060// use it!      
 1070const double m12 = M12( comb ) ;

AMAXCHILD, C++ type LoKi::AParticles::MaxChild

Simple function which evaluates the maximal value of another function over the combination components:

Search for the maximal transverse momentum:

 1000AFun maxPT = AMAXCHILD(PT) ;
 1010    
 1020const LHCb::Particle::ConstVector& combination = ... ;
 1030   
 1040const double ptMax = maxPT ( combination ) ;
Search for the maximal transverse momentum of pion:
 1000AFun maxPT = AMAXCHILD ( PT , "pi+" == ABSID ) ;
 1010
 1020const LHCb::Particle::ConstVector& combination = ... ;
 1030
 1040const double ptMax = maxPT ( combination ) ;
Search for the maximal transverse momentum of pion including the pions, used for recontruction of daughter particles:
 1000AFun maxPT = AMAXCHILD ( MAXTREE ( PT , "pi+" == ABSID )  ) ;
 1010
 1020const LHCb::Particle::ConstVector& combination = ... ;
 1030const double ptMax = maxPT ( combination ) ;

AMAXDOCA, C++ type LoKi::AParticles::MaxDOCA

The evaluator of the maximal distance of the closest approach for all two-particle subcombinations, using IGeomDispCalculator tool:

 1000IGeomDispCalculator* doca = ... ;
 1010AFun docaMAX = AMAXDOCA ( doca ) ;
 1020     
 1030const LHCb::Particle::ConstVector& combination = ... ;
 1040     
 1050const double maxDOCA = fun ( combination ) ;
 1060     

AMINCHILD, C++ type LoKi::AParticles::MinChild

Simple function which evaluates the minimal value of another function over the combination components:

Search for the minimum transverse momentum:

 1000AFun minPT = AMINCHILD(PT) ;
 1010    
 1020const LHCb::Particle::ConstVector& combination = ... ;
 1030   
 1040const double ptMin = minPT ( combination ) ;
Search for the minimum transverse momentum of pion:
 1000AFun minPT = AMINCHILD ( PT , "pi+" == ABSID ) ;
 1010
 1020const LHCb::Particle::ConstVector& combination = ... ;
 1030
 1040const double ptMin = minPT ( combination ) ;
Search for the minimum transverse momentum of pion including the pions, used for recontruction of daughter particles:
 1000AFun minPT = AMINCHILD ( MINTREE ( PT , "pi+" == ABSID )  ) ;
 1010
 1020const LHCb::Particle::ConstVector& combination = ... ;
 1030const double ptMin = minPT ( combination ) ;

AMOM, C++ type LoKi::AParticles::Momentum

Simple evaluator of the scalar momentum for the various subcombinations from thecombination itself:

 1000// get the combination:
 1010const LHCb::Paricle::ConstVector& comb = ... ;
 1020
 1030// the functor which evaluates the momentum of 1st 2nd and 4th particles:
 1040const AFun p124 = AMOM(1,2,4) ;
 1050
 1060// use it!      
 1070const double mom = p124( comb ) ;

AMOMT, C++ type LoKi::AParticles::TransverseMomentum

Simple evaluator of the transverse momentum for the various subcombinations from thecombination itself:

 1000// get the combination:
 1010const LHCb::Paricle::ConstVector& comb = ... ;
 1020
 1030// the functor which evaluates the transverse momentum of 1st 2nd and 4th particles:
 1040const AFun pt124 = AMOMT(1,2,4) ;
 1050
 1060// use it!      
 1070const double pt = pt124( comb ) ;

ANUM, C++ type LoKi::AParticles::Count

Simple functor to count number of particles in a combination, which satisfy certain criteria. It uses the algorithm LoKi::Algs::count_if

 1000const LHCb::Particle::ConstVector& combination = ... ;
 1010// count number of positive kaons: 
 1020//    1) create the functor:
 1030const AFun nK = ANUM ( "K+" == ID ) ;
 1040//    2) use it!
 1050const double nKplus = nK( combination ) ;

AP, the instance of LoKi::AParticles::Momentum(0)

Simple evaluator of the scalar momentum for the combination itself:

 1000// get the combination:
 1010const LHCb::Paricle::ConstVector& comb = ... ;
 1020
 1030const double mom = AP ( comb ) ;
See AMOM

AP0, the instance of LoKi::AParticles::Momentum(0)

Simple evaluator of the scalar momentum for the combination itself:

 1000// get the combination:
 1010const LHCb::Paricle::ConstVector& comb = ... ;
 1020
 1030const double mom = AP0 ( comb ) ;
See AMOM

AP1, the instance of LoKi::AParticles::Momentum(1)

Simple evaluator of the scalar momentum for the first daughter particle:

 1000// get the combination:
 1010const LHCb::Paricle::ConstVector& comb = ... ;
 1020
 1030const double p1 = AP1 ( comb ) ;
See AMOM

AP2, the instance of LoKi::AParticles::Momentum(2)

Simple evaluator of the scalar momentum for the second daughter:

 1000// get the combination:
 1010const LHCb::Paricle::ConstVector& comb = ... ;
 1020
 1030const double p2 = AP2 ( comb ) ;
See AMOM

AP3, the instance of LoKi::AParticles::Momentum(3)

Simple evaluator of the scalar momentum for the third daughter particle:

 1000// get the combination:
 1010const LHCb::Paricle::ConstVector& comb = ... ;
 1020
 1030const double p3 = AP3 ( comb ) ;
See AMOM

AP4, the instance of LoKi::AParticles::Momentum(4)

Simple evaluator of the scalar momentum for the fourth daughter particle:

 1000// get the combination:
 1010const LHCb::Paricle::ConstVector& comb = ... ;
 1020
 1030const double p4 = AP4 ( comb ) ;
See AMOM

APLOT, C++ type LoKi::Monitoring::Plot<LoKi::ATypes::Combination>

The special monitoring function, which could be used for monitoring of other functions. The histogram (through its AIDA::IHistogram1D interface) is filled with the results of the function evaluation.

 1000// the function to be monitored
 1010AFun fun = ... ;
 1020
 1030// the histogram
 1040AIDA::IHistogram1D* histo = ... ;
 1050
 1060// create the monitoring function:
 1070AFun mon = APLOT ( fun , histo ) ;
 1080
 1090// use the function:
 1100for ( .... ) 
 1110{
 1120   ...
 1130   const LHCb::Particle::ConstVector& comb = ... ;
 1140   const double value = mon ( comb ) ; ///< use the monitored function
 1150   ...
 1160}
The alternative and recommended way for creation of of monitored functions:
 1000// the function to be monitored 
 1010AFun fun = ... ;
 1020
 1030// the monitored function:
 1040AIDA::IHistogram1D* histo = ... ; 
 1050AFun mon = monitor ( fun , histo ) ;
 1060     

The substitution of the function by the monitored function could be done "on-flight" without the disturbing of the actual processing:

 1000// The function to be monitored
 1010AFun fun = ... ;
 1020     
 1030if ( monitoring ) 
 1040 {
 1050    AIDA::IHistogram1D* histo = ... ;
 1060    fun = monitor ( fun , counter ( "Efficiency1" ) ) ;
 1070 }
 1080    
Attention:
  • The string representation of the object is delegated to the underlying predicate, therefore the object is NOT reconstructable from its string representations. It is done on-purpose to avoid the disturbing of unique function IDs.

APT, the instance of LoKi::AParticles::TransverseMomentum(0)

Simple evaluator of the transverse momentum for the combination itself:

 1000// get the combination:
 1010const LHCb::Paricle::ConstVector& comb = ... ;
 1020
 1030const double pt = APT ( comb ) ;
See AMOMT

#LoKiCutsAPT0

APT0, the instance of LoKi::AParticles::TransverseMomentum(0)

Simple evaluator of the transverse momentum for the combination itself:

 1000// get the combination:
 1010const LHCb::Paricle::ConstVector& comb = ... ;
 1020
 1030const double pt = APT0 ( comb ) ;
See AMOMT

APT1, the instance of LoKi::AParticles::TransverseMomentum(1)

Simple evaluator of the transverse momentum for the first daughter particle:

 1000// get the combination:
 1010const LHCb::Paricle::ConstVector& comb = ... ;
 1020
 1030const double pt1 = APT1 ( comb ) ;
See AMOMT

APT2, the instance of LoKi::AParticles::TransverseMomentum(2)

Simple evaluator of the transverse momentum for the second daughter particle:

 1000// get the combination:
 1010const LHCb::Paricle::ConstVector& comb = ... ;
 1020
 1030const double pt2 = APT2 ( comb ) ;
See AMOMT

APT3, the instance of LoKi::AParticles::TransverseMomentum(3)

Simple evaluator of the transverse momentum for the third daughter particle:

 1000// get the combination:
 1010const LHCb::Paricle::ConstVector& comb = ... ;
 1020
 1030const double pt3 = APT3 ( comb ) ;
See AMOMT

APT4, the instance of LoKi::AParticles::TransverseMomentum(4)

Simple evaluator of the transverse momentum for the fourth daughter particle:

 1000// get the combination:
 1010const LHCb::Paricle::ConstVector& comb = ... ;
 1020
 1030const double pt4 = APT4 ( comb ) ;
See AMOMT

ASIZE, the instance of LoKi::AParticles::Size()

Simple functor to get the size of the combination

 1000const LHCb::Particle::ConstVector& combination = ... ;
 1010/// empty ???
 1020const bool empty = 1 >  ASIZE ( combination ) ;

TrSSWITCH, C++ type LoKi::SimpleSwitch<LoKi::ATypes::Combination>

The simple function whcuh acts according to the symbolic rule:

result = condition ( argument ) ? value1 : value2
Essentially it is the simplest way to convert predicate into function:
 1000// some complicated condition/predicate:
 1010ACut cut = ... ;
 1020
 1030// effectively "convert" it into function:
 1040AFun fun = ASSWITCH ( cut , 1 , -1 ) ; 

TrSTAT, C++ type LoKi::Monitoring::Stat<LoKi::ATypes::Combination>

The special monitoring function, which could be used for monitoring of other functions. The generic Gaudi counter StatEntity is used for monitoring:

 1000// the function to be monitored
 1010AFun fun = ... ;
 1020
 1030// the generic counter:
 1040StatEntity* counter = ... ;
 1050
 1060// create the monitoring function:
 1070AFun mon = ASTAT ( fun , counter ) ;
 1080
 1090// use the function:
 1100for ( .... ) 
 1110{
 1120   ...
 1130   const LHCb::Particle::ConstVector& comb = ... ;
 1140   const double value = mon ( comb ) ; ///< use the monitored function
 1150   ...
 1160}
 1170...
 1180// inspect the counter:
 1190info () 
 1200   << " Monitoring results : "                 << endreq 
 1210   << " NEntries:  #" << counter->entries  ()  << endreq 
 1220   << " TotalSum:   " << counter->flag     ()  << endreq 
 1230   << " Mean+-RMS:  " << counter->flagMean () 
 1240   << "+="            << counter->flagRMS  ()  << endreq      
 1250   << " Min/Max:    " << counter->flagMin  ()  
 1260   << "/"             << counter->flagMax  ()  << endreq ;
The alternative and recommended way for creation of of monitored functions:
 1000// the function to be monitored 
 1010AFun fun = ... ;
 1020
 1030// the monitored function: 
 1040AFun mon = monitor ( fun , counter ( "SomeEffCounter" ) ) ;
 1050     

The substitution of the function by the monitored function could be done "on-flight" without the disturbing of the actual processing:

 1000// The function to be monitored
 1010AFun fun = ... ;
 1020     
 1030if ( monitoring ) 
 1040 {
 1050    fun = monitor ( fun , counter ( "Efficiency1" ) ) ;
 1060 }
 1070    
Attention:
  • The string representation of the object is delegated to the underlying predicate, therefore the object is NOT reconstructable from its string representations. It is done on-purpose to avoid the disturbing of unique function IDs.

ASWITCH, C++ type LoKi::Switch<LoKi::ATypes::Combination>

The simple function which acts according to the symbolic rule:

result = condition ( argument ) ? function1 ( argument )  : function2 ( argument )
 1000// use negative Pt for invalid combinations 
 1010AFun pt = TrSWITCH ( AVALID , APT , -1 * TeV  ) ; 

AUNIQUE, C++ type LoKi::AParticles::Unique

Simple predicate to check the overlaps, using = ICheckOverlap tool

 1000const LHCb::Particle::ConstVector& combination = ... ;
 1010// get the tool:
 1020ICheckOverlap* tool = ... ;
 1030// create the predicate:
 1040const ACut unique = AUNIQUE( tool ) ;
 1050     
 1060// check the uniqueness:
 1070const bool OK = unique ( combination ) ;     

AVALID, the instance of LoKi::Valid<LoKi::ATypes::Combination>()

 1000const LHCb::Particle::ConstVector& combination = ... ;
 1010/// empty ???
 1020const bool empty = ! AVALID ( combination ) ;

AVCHI2, C++ type LoKi::AParticles::VertexChi2

The function which evaluates of the vertex fit for the combination, using IVertexFit tool:

 1000IVertexFit* fitter = ... ;
 1010
 1020AFun vChi2 = AVCHI2 ( fitter ) ;
 1030
 1040const LHCb::Particle::ConstVector& combination = ... ;
 1050
 1060const double chi2 = fun ( combination ) ;

AWM, C++ type LoKi::AParticles::WrongMass

Simple evaluator of the invariant mass of combination using some prescribed ("wrong") masses for daughter particles. It could be used e.g. for study of various mass-refelction of for generic topological selection using various mass-windows for different mass hypothesis

 1000const LHCb::Particle::ConstVector& comobination = ... ;
 1010
 1020// consider the mass hypotheses of 3 kaons:
 1030AFun m3K = AWM ( "K+" , "K+" , "K-" ) ;
 1040const double mass3kaon = m3k ( combination ) ;
 1050
 1060// consider the another mass hypotheses:
 1070AFun m2Kpi = AWM ( "K+" , "K+" , "pi+" ) ;
 1080
 1090const double mass2kaonAndPion = m2Kpi ( combination ) ;
 1100    

AWRONGMASS, C++ type LoKi::AParticles::WrongMass

Simple evaluator of the invariant mass of combination using some prescribed ("wrong") masses for daughter particles. It could be used e.g. for study of various mass-refelction of for generic topological selection using various mass-windows for different mass hypothesis

 1000const LHCb::Particle::ConstVector& comobination = ... ;
 1010
 1020// consider the mass hypotheses of 3 kaons:
 1030AFun m3K = AWRONGMASS ( "K+" , "K+" , "K-" ) ;
 1040const double mass3kaon = m3k ( combination ) ;
 1050
 1060// consider the another mass hypotheses:
 1070AFun m2Kpi = AWRONGMASS  ( "K+" , "K+" , "pi+" ) ;
 1080
 1090const double mass2kaonAndPion = m2Kpi ( combination ) ;
 1100    


-- Vanya Belyaev - 21 Jul 2007

Edit | Attach | Watch | Print version | History: r2 < r1 | Backlinks | Raw View | WYSIWYG | More topic actions
Topic revision: r2 - 2007-07-28 - VanyaBelyaev
 
    • Cern Search Icon Cern Search
    • TWiki Search Icon TWiki Search
    • Google Search Icon Google Search

    LHCb All webs login

This site is powered by the TWiki collaboration platform Powered by PerlCopyright & 2008-2019 by the contributing authors. All material on this collaboration platform is the property of the contributing authors.
Ideas, requests, problems regarding TWiki? Send feedback