ATLAS Offline Software
Loading...
Searching...
No Matches
xAODTruthParticleSlimmerTau Class Reference

Algorithm to skim the xAOD truth particle container for tau filter. More...

#include <xAODTruthParticleSlimmerTau.h>

Inheritance diagram for xAODTruthParticleSlimmerTau:
Collaboration diagram for xAODTruthParticleSlimmerTau:

Public Member Functions

 xAODTruthParticleSlimmerTau (const std::string &name, ISvcLocator *svcLoc)
 Regular algorithm constructor.
virtual StatusCode initialize ()
 Function initialising the algorithm.
virtual StatusCode execute ()
 Function executing the algorithm.
CLHEP::HepLorentzVector sumDaughterNeutrinos (const xAOD::TruthParticle *tau)
virtual StatusCode sysInitialize () override
 Override sysInitialize.
virtual const DataObjIDColl & extraOutputDeps () const override
 Return the list of extra output dependencies.
ServiceHandle< StoreGateSvc > & evtStore ()
 The standard StoreGateSvc (event store) Returns (kind of) a pointer to the StoreGateSvc.
const ServiceHandle< StoreGateSvc > & detStore () const
 The standard StoreGateSvc/DetectorStore Returns (kind of) a pointer to the StoreGateSvc.
virtual StatusCode sysStart () override
 Handle START transition.
virtual std::vector< Gaudi::DataHandle * > inputHandles () const override
 Return this algorithm's input handles.
virtual std::vector< Gaudi::DataHandle * > outputHandles () const override
 Return this algorithm's output handles.
Gaudi::Details::PropertyBase & declareProperty (Gaudi::Property< T, V, H > &t)
void updateVHKA (Gaudi::Details::PropertyBase &)
MsgStream & msg () const
bool msgLvl (const MSG::Level lvl) const

Protected Member Functions

void renounceArray (SG::VarHandleKeyArray &handlesArray)
 remove all handles from I/O resolution
std::enable_if_t< std::is_void_v< std::result_of_t< decltype(&T::renounce)(T)> > &&!std::is_base_of_v< SG::VarHandleKeyArray, T > &&std::is_base_of_v< Gaudi::DataHandle, T >, void > renounce (T &h)
void extraDeps_update_handler (Gaudi::Details::PropertyBase &ExtraDeps)
 Add StoreName to extra input/output deps as needed.

Private Types

typedef ServiceHandle< StoreGateSvcStoreGateSvc_t

Private Member Functions

Gaudi::Details::PropertyBase & declareGaudiProperty (Gaudi::Property< T, V, H > &hndl, const SG::VarHandleKeyType &)
 specialization for handling Gaudi::Property<SG::VarHandleKey>

Private Attributes

SG::ReadHandleKey< xAOD::TruthParticleContainerm_xaodTruthParticleContainerName {this,"xAODTruthParticleContainerName","TruthParticles","Name of Truth Particle container"}
SG::WriteHandleKey< xAOD::TruthParticleContainerm_xaodTruthTauParticleContainerName {this, "xAODTruthTauParticleContainerName","TruthTaus","Name of Truth Taus contatiner from the slimmer"}
 The key for the output xAOD truth containers.
DoubleProperty m_tau_pt_selection {this, "tau_pt_selection", 0.001 * Gaudi::Units::GeV}
 Selection values for keeping taus and leptons.
DoubleProperty m_abseta_selection {this, "abseta_selection", 10.}
BooleanProperty m_forceRerun {this, "ForceRerun", false}
 a flag to force rerunning (useful for rerunning on ESDs)
PublicToolHandle< IMCTruthClassifierm_classifier {this, "MCTruthClassifier", "MCTruthClassifier/MCTruthClassifier"}
DataObjIDColl m_extendedExtraObjects
StoreGateSvc_t m_evtStore
 Pointer to StoreGate (event store by default)
StoreGateSvc_t m_detStore
 Pointer to StoreGate (detector store by default)
std::vector< SG::VarHandleKeyArray * > m_vhka
bool m_varHandleArraysDeclared

Detailed Description

Algorithm to skim the xAOD truth particle container for tau filter.

This algorithm is used to copy and skim the particles from the xAOD TruthParticles container, keeping just relevant taus from the event.
The design of this class heavily mirrors the DerivationFramework::TruthCollectionMaker.

Author
Jeff Dandoy Jeff..nosp@m.Dand.nosp@m.oy@ce.nosp@m.rn.c.nosp@m.h

Definition at line 28 of file xAODTruthParticleSlimmerTau.h.

Member Typedef Documentation

◆ StoreGateSvc_t

typedef ServiceHandle<StoreGateSvc> AthCommonDataStore< AthCommonMsg< Algorithm > >::StoreGateSvc_t
privateinherited

Definition at line 388 of file AthCommonDataStore.h.

Constructor & Destructor Documentation

◆ xAODTruthParticleSlimmerTau()

xAODTruthParticleSlimmerTau::xAODTruthParticleSlimmerTau ( const std::string & name,
ISvcLocator * svcLoc )

Regular algorithm constructor.

Definition at line 27 of file xAODTruthParticleSlimmerTau.cxx.

28 : AthAlgorithm(name, svcLoc)
29{
30}
AthAlgorithm()
Default constructor:

Member Function Documentation

◆ declareGaudiProperty()

Gaudi::Details::PropertyBase & AthCommonDataStore< AthCommonMsg< Algorithm > >::declareGaudiProperty ( Gaudi::Property< T, V, H > & hndl,
const SG::VarHandleKeyType &  )
inlineprivateinherited

specialization for handling Gaudi::Property<SG::VarHandleKey>

Definition at line 156 of file AthCommonDataStore.h.

158 {
160 hndl.value(),
161 hndl.documentation());
162
163 }
Gaudi::Details::PropertyBase & declareProperty(Gaudi::Property< T, V, H > &t)

◆ declareProperty()

Gaudi::Details::PropertyBase & AthCommonDataStore< AthCommonMsg< Algorithm > >::declareProperty ( Gaudi::Property< T, V, H > & t)
inlineinherited

Definition at line 145 of file AthCommonDataStore.h.

145 {
146 typedef typename SG::HandleClassifier<T>::type htype;
148 }
Gaudi::Details::PropertyBase & declareGaudiProperty(Gaudi::Property< T, V, H > &hndl, const SG::VarHandleKeyType &)
specialization for handling Gaudi::Property<SG::VarHandleKey>

◆ detStore()

const ServiceHandle< StoreGateSvc > & AthCommonDataStore< AthCommonMsg< Algorithm > >::detStore ( ) const
inlineinherited

The standard StoreGateSvc/DetectorStore Returns (kind of) a pointer to the StoreGateSvc.

Definition at line 95 of file AthCommonDataStore.h.

◆ evtStore()

ServiceHandle< StoreGateSvc > & AthCommonDataStore< AthCommonMsg< Algorithm > >::evtStore ( )
inlineinherited

The standard StoreGateSvc (event store) Returns (kind of) a pointer to the StoreGateSvc.

Definition at line 85 of file AthCommonDataStore.h.

◆ execute()

StatusCode xAODTruthParticleSlimmerTau::execute ( )
virtual

Function executing the algorithm.

Definition at line 61 of file xAODTruthParticleSlimmerTau.cxx.

62{
63
64 CLHEP::HepLorentzVector nutau;
65
66 // If the containers already exists then assume that nothing needs to be done
68 {
69 ATH_MSG_WARNING("xAOD Tau Truth Particles are already available in the event");
70 return StatusCode::SUCCESS;
71 }
72
73 // Create new output container
74 SG::WriteHandle<xAOD::TruthParticleContainer> xTruthTauParticleContainer(m_xaodTruthTauParticleContainerName);
75 ATH_CHECK(xTruthTauParticleContainer.record(std::make_unique<xAOD::TruthParticleContainer>(), std::make_unique<xAOD::TruthParticleAuxContainer>()));
76 ATH_MSG_INFO("Recorded TruthTauParticleContainer with key: " << m_xaodTruthTauParticleContainerName.key());
77
78 // Retrieve full TruthParticle container
79 SG::ReadHandle<xAOD::TruthParticleContainer> xTruthParticleContainer{m_xaodTruthParticleContainerName};
80 if ( !xTruthParticleContainer.isValid() )
81 {
82 ATH_MSG_ERROR("No TruthParticle collection with name " << m_xaodTruthParticleContainerName.key() << " found in StoreGate!");
83 return StatusCode::FAILURE;
84 }
85
86 // Set up decorators
87
88 const static SG::AuxElement::Decorator<unsigned int> originDecorator("classifierParticleOrigin");
89 const static SG::AuxElement::Decorator<unsigned int> typeDecorator("classifierParticleType");
90 const static SG::AuxElement::Decorator<unsigned int> outcomeDecorator("classifierParticleOutCome");
91 const static SG::AuxElement::Decorator<unsigned int> classificationDecorator("Classification");
92 const static SG::AuxElement::Decorator<int> parenthadronPIDDecorator("parentHadronID");
93
94 // sum of neutrinos 4-vector in tau decay products
95 const static SG::AuxElement::Decorator<CLHEP::HepLorentzVector> nuDecorator("nuVector");
96 const static SG::AuxElement::Decorator<int> tauTypeDecorator("tauType");
97
98 // Loop over full TruthParticle container
99
100 std::vector<int> uniqueID_list;
101 int zero_uniqueID=0;
102 int dup_uniqueID = 0;
103
104 unsigned int nParticles = xTruthParticleContainer->size();
105
106 for (unsigned int iPart = 0; iPart < nParticles; ++iPart)
107 {
108 ElementLink<xAOD::TruthParticleContainer> eltp(*xTruthParticleContainer, iPart);
109 const xAOD::TruthParticle *theParticle = (*xTruthParticleContainer)[iPart];
110
111 int my_uniqueID = HepMC::uniqueID(theParticle);
112 if ( my_uniqueID == HepMC::UNDEFINED_ID ) {
113 zero_uniqueID++;
114 continue;
115 }
116 bool found = false;
117 if (uniqueID_list.size() > 0){
118 found = (std::find(uniqueID_list.begin(), uniqueID_list.end(), my_uniqueID) != uniqueID_list.end());
119 if(found) {
120 dup_uniqueID++;
121 continue;}
122 }
123 uniqueID_list.push_back(my_uniqueID);
124
125 float this_abseta = theParticle->abseta();
126 float this_pt = theParticle->pt();
127
128
129 //Save Taus above 0.001 GeV, & with any eta (may be changed on JOs level eg. to dectector acceptance of eta 4.5)
130 // see GeneratorFilters/share/common/xAODTauFilter_Common.py
131 // we want to avoid status 3 taus
132 if (MC::isPhysical(theParticle) && MC::isTau(theParticle) && this_pt >= m_tau_pt_selection && this_abseta < m_abseta_selection)
133 {
134 xAOD::TruthParticle *xTruthParticle = new xAOD::TruthParticle();
135
136 xTruthTauParticleContainer->push_back(xTruthParticle);
137 // Fill with numerical content
138 *xTruthParticle=*theParticle;
139 const xAOD::TruthParticle *tau = theParticle;
140 nutau = sumDaughterNeutrinos(tau);
141
142 nuDecorator(*xTruthParticle) = nutau;
143
144 int tauType = 0;
145 for (size_t n = 0; n < tau->nChildren(); ++n)
146 {
147 if (tau->child(n)->absPdgId() == MC::NU_E)
148 tauType = 1; //Tau decays into an electron
149 else if (tau->child(n)->absPdgId() == MC::NU_MU)
150 tauType = 2; //Tau decays into a muon
151 else if (MC::isTau(tau->child(n)))
152 tauType = 11; //Tau radiates a particle and decays into another tau
153 }
154 tauTypeDecorator(*xTruthParticle) = tauType;
155
156unsigned int particleOutCome;
157unsigned int result;
158unsigned int particleType;
159unsigned int particleOrigin;
160int hadron_pdg;
161Common::classify(m_classifier,theParticle,particleOutCome,result,hadron_pdg,particleType,particleOrigin );
162 typeDecorator(*xTruthParticle) = particleType;
163 originDecorator(*xTruthParticle) = particleOrigin;
164 outcomeDecorator(*xTruthParticle) = particleOutCome;
165
166 classificationDecorator(*xTruthParticle) = result;
167 parenthadronPIDDecorator(*xTruthParticle) = hadron_pdg;
168 }
169
170 } //end of loop over particles
171 if(zero_uniqueID!=0 || dup_uniqueID != 0) ATH_MSG_INFO("Found " << zero_uniqueID << " uniqueID=0 particles and " << dup_uniqueID << " duplicated");
172 return StatusCode::SUCCESS;
173}
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_INFO(x)
#define ATH_MSG_WARNING(x)
SG::Decorator< T, ALLOC > Decorator
Definition AuxElement.h:575
virtual bool isValid() override final
Can the handle be successfully dereferenced?
SG::ReadHandleKey< xAOD::TruthParticleContainer > m_xaodTruthParticleContainerName
SG::WriteHandleKey< xAOD::TruthParticleContainer > m_xaodTruthTauParticleContainerName
The key for the output xAOD truth containers.
PublicToolHandle< IMCTruthClassifier > m_classifier
CLHEP::HepLorentzVector sumDaughterNeutrinos(const xAOD::TruthParticle *tau)
DoubleProperty m_tau_pt_selection
Selection values for keeping taus and leptons.
const TruthParticle_v1 * child(size_t i) const
Retrieve the i-th mother (TruthParticle) of this TruthParticle.
double abseta() const
The absolute pseudorapidity ( ) of the particle.
int absPdgId() const
Absolute PDG ID code (often useful)
virtual double pt() const override final
The transverse momentum ( ) of the particle.
size_t nChildren() const
Number of children of this particle.
bool contains(const std::string &s, const std::string &regx)
does a string contain the substring
Definition hcg.cxx:114
void classify(ToolHandle< IMCTruthClassifier > &m_classif, const xAOD::TruthParticle *theParticle, unsigned int &particleOutCome, unsigned int &result, int &hadron_pdg, unsigned int &particleType, unsigned int &particleOrigin)
Definition Common.cxx:96
int uniqueID(const T &p)
constexpr int UNDEFINED_ID
static const int NU_MU
static const int NU_E
bool isTau(const T &p)
bool isPhysical(const T &p)
Identify if the particle is physical, i.e. is stable or decayed.
TruthParticle_v1 TruthParticle
Typedef to implementation.

◆ extraDeps_update_handler()

void AthCommonDataStore< AthCommonMsg< Algorithm > >::extraDeps_update_handler ( Gaudi::Details::PropertyBase & ExtraDeps)
protectedinherited

Add StoreName to extra input/output deps as needed.

use the logic of the VarHandleKey to parse the DataObjID keys supplied via the ExtraInputs and ExtraOuputs Properties to add the StoreName if it's not explicitly given

◆ extraOutputDeps()

const DataObjIDColl & AthAlgorithm::extraOutputDeps ( ) const
overridevirtualinherited

Return the list of extra output dependencies.

This list is extended to include symlinks implied by inheritance relations.

Definition at line 50 of file AthAlgorithm.cxx.

51{
52 // If we didn't find any symlinks to add, just return the collection
53 // from the base class. Otherwise, return the extended collection.
54 if (!m_extendedExtraObjects.empty()) {
56 }
57 return Algorithm::extraOutputDeps();
58}
DataObjIDColl m_extendedExtraObjects

◆ initialize()

StatusCode xAODTruthParticleSlimmerTau::initialize ( )
virtual

Function initialising the algorithm.

Definition at line 32 of file xAODTruthParticleSlimmerTau.cxx.

33{
34 ATH_CHECK(m_classifier.retrieve());
36 ATH_MSG_INFO("xAOD input TruthParticleContainer name = " << m_xaodTruthParticleContainerName.key());
38 ATH_MSG_INFO("xAOD output TruthTauParticleContainer name = " << m_xaodTruthTauParticleContainerName.key());
39 return StatusCode::SUCCESS;
40}

◆ inputHandles()

virtual std::vector< Gaudi::DataHandle * > AthCommonDataStore< AthCommonMsg< Algorithm > >::inputHandles ( ) const
overridevirtualinherited

Return this algorithm's input handles.

We override this to include handle instances from key arrays if they have not yet been declared. See comments on updateVHKA.

◆ msg()

MsgStream & AthCommonMsg< Algorithm >::msg ( ) const
inlineinherited

Definition at line 24 of file AthCommonMsg.h.

24 {
25 return this->msgStream();
26 }

◆ msgLvl()

bool AthCommonMsg< Algorithm >::msgLvl ( const MSG::Level lvl) const
inlineinherited

Definition at line 30 of file AthCommonMsg.h.

30 {
31 return this->msgLevel(lvl);
32 }

◆ outputHandles()

virtual std::vector< Gaudi::DataHandle * > AthCommonDataStore< AthCommonMsg< Algorithm > >::outputHandles ( ) const
overridevirtualinherited

Return this algorithm's output handles.

We override this to include handle instances from key arrays if they have not yet been declared. See comments on updateVHKA.

◆ renounce()

std::enable_if_t< std::is_void_v< std::result_of_t< decltype(&T::renounce)(T)> > &&!std::is_base_of_v< SG::VarHandleKeyArray, T > &&std::is_base_of_v< Gaudi::DataHandle, T >, void > AthCommonDataStore< AthCommonMsg< Algorithm > >::renounce ( T & h)
inlineprotectedinherited

Definition at line 380 of file AthCommonDataStore.h.

381 {
382 h.renounce();
384 }
std::enable_if_t< std::is_void_v< std::result_of_t< decltype(&T::renounce)(T)> > &&!std::is_base_of_v< SG::VarHandleKeyArray, T > &&std::is_base_of_v< Gaudi::DataHandle, T >, void > renounce(T &h)

◆ renounceArray()

void AthCommonDataStore< AthCommonMsg< Algorithm > >::renounceArray ( SG::VarHandleKeyArray & handlesArray)
inlineprotectedinherited

remove all handles from I/O resolution

Definition at line 364 of file AthCommonDataStore.h.

364 {
366 }

◆ sumDaughterNeutrinos()

CLHEP::HepLorentzVector xAODTruthParticleSlimmerTau::sumDaughterNeutrinos ( const xAOD::TruthParticle * tau)

Definition at line 42 of file xAODTruthParticleSlimmerTau.cxx.

43{
44 CLHEP::HepLorentzVector nu(0, 0, 0, 0);
45 if (MC::isSMNeutrino(part) && MC::isPhysical(part))
46 {
47 nu.setPx(part->px());
48 nu.setPy(part->py());
49 nu.setPz(part->pz());
50 nu.setE(part->e());
51 }
52 if (!part->hasDecayVtx())
53 return nu;
54
55 for (size_t n = 0; n < part->nChildren(); ++n)
56 nu += sumDaughterNeutrinos(part->child(n));
57
58 return nu;
59}
bool isSMNeutrino(const T &p)

◆ sysInitialize()

StatusCode AthAlgorithm::sysInitialize ( )
overridevirtualinherited

Override sysInitialize.

Override sysInitialize from the base class.

Loop through all output handles, and if they're WriteCondHandles, automatically register them and this Algorithm with the CondSvc

Scan through all outputHandles, and if they're WriteCondHandles, register them with the CondSvc

Reimplemented from AthCommonDataStore< AthCommonMsg< Algorithm > >.

Reimplemented in AthAnalysisAlgorithm, AthFilterAlgorithm, AthHistogramAlgorithm, and PyAthena::Alg.

Definition at line 66 of file AthAlgorithm.cxx.

66 {
68
69 if (sc.isFailure()) {
70 return sc;
71 }
72 ServiceHandle<ICondSvc> cs("CondSvc",name());
73 for (auto h : outputHandles()) {
74 if (h->isCondition() && h->mode() == Gaudi::DataHandle::Writer) {
75 // do this inside the loop so we don't create the CondSvc until needed
76 if ( cs.retrieve().isFailure() ) {
77 ATH_MSG_WARNING("no CondSvc found: won't autoreg WriteCondHandles");
78 return StatusCode::SUCCESS;
79 }
80 if (cs->regHandle(this,*h).isFailure()) {
81 sc = StatusCode::FAILURE;
82 ATH_MSG_ERROR("unable to register WriteCondHandle " << h->fullKey()
83 << " with CondSvc");
84 }
85 }
86 }
87 return sc;
88}
static Double_t sc
virtual StatusCode sysInitialize() override
Override sysInitialize.
AthCommonDataStore(const std::string &name, T... args)
virtual std::vector< Gaudi::DataHandle * > outputHandles() const override
::StatusCode StatusCode
StatusCode definition for legacy code.

◆ sysStart()

virtual StatusCode AthCommonDataStore< AthCommonMsg< Algorithm > >::sysStart ( )
overridevirtualinherited

Handle START transition.

We override this in order to make sure that conditions handle keys can cache a pointer to the conditions container.

◆ updateVHKA()

void AthCommonDataStore< AthCommonMsg< Algorithm > >::updateVHKA ( Gaudi::Details::PropertyBase & )
inlineinherited

Definition at line 308 of file AthCommonDataStore.h.

308 {
309 // debug() << "updateVHKA for property " << p.name() << " " << p.toString()
310 // << " size: " << m_vhka.size() << endmsg;
311 for (auto &a : m_vhka) {
313 for (auto k : keys) {
314 k->setOwner(this);
315 }
316 }
317 }
std::vector< SG::VarHandleKeyArray * > m_vhka

Member Data Documentation

◆ m_abseta_selection

DoubleProperty xAODTruthParticleSlimmerTau::m_abseta_selection {this, "abseta_selection", 10.}
private

Definition at line 50 of file xAODTruthParticleSlimmerTau.h.

50{this, "abseta_selection", 10.};

◆ m_classifier

PublicToolHandle<IMCTruthClassifier> xAODTruthParticleSlimmerTau::m_classifier {this, "MCTruthClassifier", "MCTruthClassifier/MCTruthClassifier"}
private

Definition at line 55 of file xAODTruthParticleSlimmerTau.h.

55{this, "MCTruthClassifier", "MCTruthClassifier/MCTruthClassifier"};

◆ m_detStore

StoreGateSvc_t AthCommonDataStore< AthCommonMsg< Algorithm > >::m_detStore
privateinherited

Pointer to StoreGate (detector store by default)

Definition at line 393 of file AthCommonDataStore.h.

◆ m_evtStore

StoreGateSvc_t AthCommonDataStore< AthCommonMsg< Algorithm > >::m_evtStore
privateinherited

Pointer to StoreGate (event store by default)

Definition at line 390 of file AthCommonDataStore.h.

◆ m_extendedExtraObjects

DataObjIDColl AthAlgorithm::m_extendedExtraObjects
privateinherited

Definition at line 79 of file AthAlgorithm.h.

◆ m_forceRerun

BooleanProperty xAODTruthParticleSlimmerTau::m_forceRerun {this, "ForceRerun", false}
private

a flag to force rerunning (useful for rerunning on ESDs)

Definition at line 53 of file xAODTruthParticleSlimmerTau.h.

53{this, "ForceRerun", false};

◆ m_tau_pt_selection

DoubleProperty xAODTruthParticleSlimmerTau::m_tau_pt_selection {this, "tau_pt_selection", 0.001 * Gaudi::Units::GeV}
private

Selection values for keeping taus and leptons.

Definition at line 49 of file xAODTruthParticleSlimmerTau.h.

49{this, "tau_pt_selection", 0.001 * Gaudi::Units::GeV}; //in GeV

◆ m_varHandleArraysDeclared

bool AthCommonDataStore< AthCommonMsg< Algorithm > >::m_varHandleArraysDeclared
privateinherited

Definition at line 399 of file AthCommonDataStore.h.

◆ m_vhka

std::vector<SG::VarHandleKeyArray*> AthCommonDataStore< AthCommonMsg< Algorithm > >::m_vhka
privateinherited

Definition at line 398 of file AthCommonDataStore.h.

◆ m_xaodTruthParticleContainerName

SG::ReadHandleKey<xAOD::TruthParticleContainer> xAODTruthParticleSlimmerTau::m_xaodTruthParticleContainerName {this,"xAODTruthParticleContainerName","TruthParticles","Name of Truth Particle container"}
private

Definition at line 42 of file xAODTruthParticleSlimmerTau.h.

43{this,"xAODTruthParticleContainerName","TruthParticles","Name of Truth Particle container"};

◆ m_xaodTruthTauParticleContainerName

SG::WriteHandleKey<xAOD::TruthParticleContainer> xAODTruthParticleSlimmerTau::m_xaodTruthTauParticleContainerName {this, "xAODTruthTauParticleContainerName","TruthTaus","Name of Truth Taus contatiner from the slimmer"}
private

The key for the output xAOD truth containers.

Definition at line 45 of file xAODTruthParticleSlimmerTau.h.

46{this, "xAODTruthTauParticleContainerName","TruthTaus","Name of Truth Taus contatiner from the slimmer"};

The documentation for this class was generated from the following files: