trigger-sim C++ API Reference¶
-
class ClusterTriggerAlgorithm : public TriggerService¶
- #include <ClusterTriggerAlgorithm.h>
Public Functions
-
ClusterTriggerAlgorithm(double triggerWindow, unsigned int triggerThreshold, unsigned int coherenceLength, int domSet, I3MapKeyVectorIntConstPtr customDomSets)¶
-
inline ~ClusterTriggerAlgorithm()¶
-
virtual void Trigger()¶
-
int GetHash(int string, unsigned int position)¶
-
int GetString(int hash)¶
-
unsigned int GetPosition(int hash)¶
-
ClusterTriggerAlgorithm(double triggerWindow, unsigned int triggerThreshold, unsigned int coherenceLength, int domSet, I3MapKeyVectorIntConstPtr customDomSets)¶
-
class CylinderTriggerAlgorithm : public TriggerService¶
- #include <CylinderTriggerAlgorithm.h>
Public Functions
-
CylinderTriggerAlgorithm(double triggerWindow, unsigned int triggerThreshold, unsigned int simpleMultiplicity, I3GeometryConstPtr Geometry, double Radius, double Height, int domSet, I3MapKeyVectorIntConstPtr customDomSets)¶
-
inline ~CylinderTriggerAlgorithm()¶
-
virtual void Trigger()¶
-
CylinderTriggerAlgorithm(double triggerWindow, unsigned int triggerThreshold, unsigned int simpleMultiplicity, I3GeometryConstPtr Geometry, double Radius, double Height, int domSet, I3MapKeyVectorIntConstPtr customDomSets)¶
-
class FaintParticleTriggerAlgorithm : public TriggerService¶
- #include <FaintParticleTriggerAlgorithm.h>
The FaintParticleTriggerAlgorithm looks for faint signatures of particles dominantly producing SLC hits and receives also SLC hits as an input. Four cuts are calculated for each time window. All cuts are described on https://wiki.icecube.wisc.edu/index.php/Faint_Particle_Trigger. Cuts 1 and 4 are calculated in the FPTTimeWindow.h class. The first cut is a hit count and the cuts are set to remove detector noise with the lower threshold and too bright signatures with the upper threshold. For the second cut velocity consistent hit pair combinations (Doubles) are calcuated to further reduce the detector noise contribution.
For the third cut the directional consistency of Doubles is calculated, which is not expected for detector noise.For each Double the corresponding zenith and azimuth angles are calculated, binned and then the number of entries of the maximum bin for both histograms is taken and compared to a cut value. For the last cut the amount of SLC hits in the time window over the amount of all hits is calculated. The amount of double triggered events (standard triggers and faint particle trigger) is reduced by selecting events with a high SLC fraction.
Public Functions
-
FaintParticleTriggerAlgorithm(double time_window, double time_window_separation, double max_trigger_length, unsigned int hit_min, unsigned int hit_max, double double_velocity_min, double double_velocity_max, unsigned int double_min, unsigned int azimuth_histogram_min, unsigned int zenith_histogram_min, double histogram_binning, double slcfraction_min, I3GeometryConstPtr Geometry, int domSet, I3MapKeyVectorIntConstPtr customDomSets)¶
- Parameters:
time_window – The time window length.
time_window_separation – The separation between consecutive time windows.
max_trigger_length – If the max_trigger_length is exceeded the trigger is formed.
hit_min – The minimum number of hits required for triggering.
hit_max – The maximum number of hits required for triggering.
double_velocity_min – The minimum velocity for hit pairs to count as a Double.
double_velocity_max – The maximum velocity for hit pairs to count as a Double.
double_min – The minimum number of Doubles required for triggering.
azimuth_histogram_min – The minimum number of entries in the azimuth histogram. (DC version)
zenith_histogram_min – The minimum number of entries in the zenith histogram.(DC version)
histogram_binning – The binning value for the azimuth and zenith histograms.
slcfraction_min – The minimum slc fraction of the time window that is requried for triggering.
Geometry – Pointer to the I3Geometry
domSet – The DOMSet
customDomSets –
-
inline ~FaintParticleTriggerAlgorithm()¶
-
virtual void Trigger()¶
-
double getDistance(TriggerHit hit1, TriggerHit hit2, I3GeometryConstPtr Geometry)¶
Private Functions
- SET_LOGGER ("FaintParticleTriggerAlgorithm")
Private Members
-
double time_window_¶
-
double time_window_separation_¶
-
double max_trigger_length_¶
-
unsigned int hit_min_¶
-
unsigned int hit_max_¶
-
double double_velocity_min_¶
-
double double_velocity_max_¶
-
unsigned int double_min_¶
-
unsigned int azimuth_histogram_min_¶
-
unsigned int zenith_histogram_min_¶
-
double histogram_binning_¶
-
double slcfraction_min_¶
-
I3GeometryConstPtr geo_¶
-
FaintParticleTriggerAlgorithm(double time_window, double time_window_separation, double max_trigger_length, unsigned int hit_min, unsigned int hit_max, double double_velocity_min, double double_velocity_max, unsigned int double_min, unsigned int azimuth_histogram_min, unsigned int zenith_histogram_min, double histogram_binning, double slcfraction_min, I3GeometryConstPtr Geometry, int domSet, I3MapKeyVectorIntConstPtr customDomSets)¶
-
class FPTTimeWindow¶
- #include <FPTTimeWindow.h>
A simple timewindow class used by the triggers.
Public Functions
-
FPTTimeWindow(unsigned int hit_min, unsigned int hit_max, double slcfraction_min, double time_window, double time_window_separation)¶
- Parameters:
time_window – The time window length.
time_window_separation – The separation between consecutive time windows.
hit_min – The minimum number of hits required for the first cut.
hit_max – The maximum number of hits required for the first cut.
slcfraction_min – The minimum slc fraction of the time window that is required for the fourth cut.
-
~FPTTimeWindow() = default¶
Private Functions
-
FPTTimeWindow()¶
Default constructor is private until we define defaults
- SET_LOGGER ("FPTTimeWindow")
-
FPTTimeWindow(unsigned int hit_min, unsigned int hit_max, double slcfraction_min, double time_window, double time_window_separation)¶
-
class GlobalTriggerSim¶
- #include <GlobalTriggerSim.h>
Public Functions
-
inline GlobalTriggerSim(const I3DetectorStatus &d)¶
-
inline GlobalTriggerSim(I3DetectorStatusConstPtr d)¶
-
inline void SetDefaultReadoutConfig(I3TriggerReadoutConfig r)¶
-
void InsertThroughputTriggers(std::vector<I3Trigger> &src, I3TriggerPairVector &dest)¶
-
I3TriggerHierarchy::iterator FindOverlap(I3TriggerHierarchyConstPtr, const I3Trigger&)¶
Private Members
-
boost::shared_ptr<std::map<TriggerKey, I3TriggerStatus>> triggerStatusMap_¶
-
boost::shared_ptr<ReadoutWindowUtil> roUtil_¶
-
I3TriggerReadoutConfig defaultReadoutConfig_¶
-
inline GlobalTriggerSim(const I3DetectorStatus &d)¶
-
class I3GlobalTriggerSim : public I3Module¶
- #include <I3GlobalTriggerSim.h>
class: I3GlobalTriggerSim
Version
- Rcs
date:
- Rcs
Copyright (c) 2006 IceCube Collaboration
- Author
wikstrom
- Author
: Seon-Hee Seo
Public Functions
-
inline ~I3GlobalTriggerSim()¶
-
void Configure()¶
-
void DAQ(I3FramePtr frame)¶
-
inline void Finish()¶
Private Functions
-
I3GlobalTriggerSim()¶
-
I3GlobalTriggerSim(const I3GlobalTriggerSim&)¶
-
I3GlobalTriggerSim &operator=(const I3GlobalTriggerSim&)¶
-
inline void PushIf(bool triggerCondition, I3FramePtr frame)¶
- SET_LOGGER ("I3GlobalTriggerSim")
-
class I3Pruner : public I3ConditionalModule¶
- #include <I3Pruner.h>
IceTray module to remove launches outside the readout window.
-
class I3TriggerSimModule : public I3Module¶
- #include <I3TriggerSimModule.h>
Public Functions
-
inline ~I3TriggerSimModule()¶
-
void Configure()¶
-
void Geometry(I3FramePtr frame)¶
-
void DetectorStatus(I3FramePtr frame)¶
-
void DAQ(I3FramePtr frame)¶
-
void Finish()¶
Private Functions
- SET_LOGGER ("I3TriggerSimModule")
-
inline ~I3TriggerSimModule()¶
-
class ReadoutWindowUtil¶
- #include <ReadoutWindowUtil.h>
Utility class to provide methods to retrieve readout windows for various trigger conditions. It’s more complicated than you may think.
class: ReadoutWindowUtil
Version
- Rcs
date:
- Rcs
Copyright (c) 2006 IceCube Collaboration
- Author
toale
Public Functions
-
ReadoutWindowUtil(const I3DetectorStatus &detectorStatus)¶
-
~ReadoutWindowUtil()¶
-
std::pair<double, double> GetReadoutWindow(I3TriggerStatus::Subdetector subdetector, const I3Trigger &trigger)¶
Return one readout window corresponding to subdetector InIce or IceTop
Private Members
-
std::map<TriggerKey, I3TriggerStatus> triggerStatus_¶
-
class SimpleMajorityTriggerAlgorithm : public TriggerService¶
- #include <SimpleMajorityTriggerAlgorithm.h>
The SimpleMajorityTriggerAlgorithm class is a port the DAQ trigger algorithm of the same name. The algorithm has 2 configuration parameters:
1) TriggerWindow: The length of the sliding time window used to define an interesting time period.
2) TriggerThreshold: The minimum number of hits that must be present within a time period of TriggerWindow.
A trigger is satisfied if there are at least TriggerThreshold hits in a time window of length TriggerWindow.
The algorithm works in one step:
1) Find a time window: A sliding time window is applied to each hit and the number of hits falling in it is counted. As long as the threshold condition is met, the window will continue to slide. Once the number of hits drops below the threshold, a trigger is defined which includes all hits that were part of any valid window. A hit cannot be part of more than one active time window.
When a trigger is found, the hits that are included are any that were in the time window. The trigger start will be the earliest hit time in the time window and the trigger stop will be the latest hit time in the time window.
A time window with multiple independent clusters of hits will only produce a single trigger, but all hits in the time window will be used to define the length of the trigger.
Public Functions
-
SimpleMajorityTriggerAlgorithm(double triggerWindow, unsigned int triggerThreshold, int domSet, I3MapKeyVectorIntConstPtr customDomSets)¶
-
inline ~SimpleMajorityTriggerAlgorithm()¶
-
virtual void Trigger()¶
Private Functions
- SET_LOGGER ("SimpleMajorityTriggerAlgorithm")
-
SimpleMajorityTriggerAlgorithm(double triggerWindow, unsigned int triggerThreshold, int domSet, I3MapKeyVectorIntConstPtr customDomSets)¶
-
class SlowMonopoleTriggerAlgorithm : public TriggerService¶
- #include <SlowMonopoleTriggerAlgorithm.h>
Public Functions
-
SlowMonopoleTriggerAlgorithm(double t_proximity, double t_min, double t_max, boost::optional<double> deltad, boost::optional<double> alpha_min, boost::optional<bool> dc_algo, double relv, int min_tuples, double max_event_length, I3GeometryConstPtr Geometry, int domSet, I3MapKeyVectorIntConstPtr customDomSets)¶
-
inline ~SlowMonopoleTriggerAlgorithm()¶
-
virtual void Trigger()¶
Private Functions
-
int GetTrigStatusSetting(const I3TriggerStatus&, const std::string&)¶
-
void ConfigureFromDetStatus(const I3DetectorStatus &detstatus)¶
-
void RunTrigger(TriggerHitVector *one_hit_list__, TriggerHitVector *two_hit_list__, double *muon_time_window__, TriggerHit new_hit, const I3GeometryConstPtr &geo)¶
-
bool HLCPairCheck(TriggerHit hit1, TriggerHit hit2)¶
-
void CheckTriggerStatus(TriggerHitVector *two_hit_list__, const I3GeometryConstPtr &geo)¶
-
void CheckTriple(TriggerHit hit1, TriggerHit hit2, TriggerHit hit3, const I3GeometryConstPtr &geo)¶
-
double getDistance(TriggerHit hit1, TriggerHit hit2, const I3GeometryConstPtr &geo)¶
- SET_LOGGER ("SlowMonopoleTrigger")
Private Members
-
double t_proximity_¶
-
double t_min_¶
-
double t_max_¶
-
boost::optional<double> deltad_¶
-
boost::optional<double> alpha_min_¶
-
boost::optional<bool> dc_algo_¶
-
double relv_¶
-
int min_tuples_¶
-
double max_event_length_¶
-
I3GeometryConstPtr geo_¶
-
double cos_alpha_min_¶
-
TriggerContainerVector trigger_container_vector¶
-
SlowMonopoleTriggerAlgorithm(double t_proximity, double t_min, double t_max, boost::optional<double> deltad, boost::optional<double> alpha_min, boost::optional<bool> dc_algo, double relv, int min_tuples, double max_event_length, I3GeometryConstPtr Geometry, int domSet, I3MapKeyVectorIntConstPtr customDomSets)¶
-
class SlowMPHit¶
- #include <SlowMPHit.h>
Slow_MP_Hit
A simple hit class used by the trigger-sim. derived from trigger-sim/TriggerHit.h
-
class TimeWindow¶
- #include <TimeWindow.h>
A simple timewindow class used by the triggers.
Public Functions
-
TimeWindow(unsigned int threshold, double window)¶
Constructor requires a threshold and a time window (in ns)
-
~TimeWindow()¶
-
TriggerHitIterPairVectorPtr SlidingTimeWindows(TriggerHitVectorPtr hits)¶
Sliding time windows
Implementation of a sliding time window algorithm. Input is a vector of TriggerHit objects that should be time ordered: std::vector<TriggerHit> Output is a std::vector of pairs comprised of iterators pointing to the begining and end of each valid time window: std::vector<pair<TriggerHitVector::const_iterator,TriqggerHitVector::const_iterator> >
-
TriggerHitIterPairVectorPtr SlidingTimeWindows(TriggerHitVector &hits)¶
-
TriggerHitIterPairVectorPtr FixedTimeWindows(TriggerHitVectorPtr hits)¶
Fixed time windows
Implementation of a sliding time window algorithm. Input is a std::vector of TriggerHit objects that should be time ordered: std::vector<TriggerHit> Output is a std::vector of pairs comprised of iterators pointing to the begining and end of each valid time window: std::vector<pair<TriggerHitVector::const_iterator,TriggerHitVector::const_iterator> >
Private Functions
-
TimeWindow()¶
Default constructor is private until we define defaults
-
void DumpHits(TriggerHitList &hits, const std::string &head, const std::string &pad)¶
-
void CopyHits(TriggerHitList &sourceHits, TriggerHitList &destHits)¶
-
bool Overlap(TriggerHitList &hits1, TriggerHitList &hits2)¶
- SET_LOGGER ("TimeWindow")
-
TimeWindow(unsigned int threshold, double window)¶
-
class TriggerContainer¶
- #include <TriggerContainer.h>
Public Functions
-
inline TriggerContainer()¶
-
inline TriggerContainer(I3TriggerPtr aSlowmptrigger)¶
-
inline TriggerContainer(I3TriggerPtr aSlowmptrigger, int aNtuples)¶
-
inline I3TriggerPtr GetTrigger()¶
-
inline int GetNTuples()¶
-
inline void IncreaseNTuples()¶
-
inline I3ParticlePtr GetPos1(int tuple_nr)¶
-
inline I3ParticlePtr GetPos2(int tuple_nr)¶
-
inline I3ParticlePtr GetPos3(int tuple_nr)¶
-
inline void SetPos1(I3ParticlePtr tupel)¶
-
inline void SetPos2(I3ParticlePtr tupel)¶
-
inline void SetPos3(I3ParticlePtr tupel)¶
-
inline TriggerContainer()¶
-
class TriggerHit¶
- #include <TriggerHit.h>
A simple hit class used by the triggers.
Public Functions
-
inline TriggerHit()¶
-
inline TriggerHit(double aTime, unsigned int aPos, int aString, bool aLc)¶
-
inline bool operator<(const TriggerHit &rhs) const¶
-
inline bool operator==(const TriggerHit &rhs) const¶
-
inline TriggerHit()¶
-
class TriggerService¶
- #include <TriggerService.h>
Subclassed by ClusterTriggerAlgorithm, CylinderTriggerAlgorithm, FaintParticleTriggerAlgorithm, SimpleMajorityTriggerAlgorithm, SlowMonopoleTriggerAlgorithm
Public Functions
-
TriggerService(int domSet, I3MapKeyVectorIntConstPtr customDomSets)¶
-
inline virtual ~TriggerService()¶
-
void FillHits(I3DOMLaunchSeriesMapConstPtr launches, I3RecoPulseSeriesMapConstPtr pulses, bool useSLC)¶
-
virtual void Trigger() = 0¶
-
unsigned int GetNumberOfTriggers()¶
-
TriggerHitVectorPtr GetNextTrigger()¶
-
TriggerService(int domSet, I3MapKeyVectorIntConstPtr customDomSets)¶
-
namespace DetectorStatusUtils¶
Contains utility functions for getting information out of the DetectorStatus frame necessary to configure the modules.
Typedefs
-
typedef std::pair<TriggerKey, I3TriggerStatus> tk_ts_pair_t¶
-
typedef std::pair<TriggerKey, I3TriggerStatus> tk_ts_pair_t¶
-
namespace DOMSetFunctions¶
Contains functions useful for handling DOMsets.
Functions
-
I3MapKeyVectorIntPtr GetDefaultDOMSets()¶
Variables
-
const std::vector<unsigned> DOMSETS = boost::assign::list_of(2)(3)(4)(5)(6)(7)(8)(9)(10)(11)¶
The DomSets are defined by an XML file that lives with the pDAQ trigger definitions, here: http://code.icecube.wisc.edu/daq/projects/config/trunk/trigger/domset-definitions.xml As of August 2022, this includes: (0) and (1): AMANDA sets not used anymore (2) InIce strings 1-78 (3) IceTop IT81 (4), (5), and (6) Different versions of DeepCore (7) Scintillators (8) IceTop Infill (HG only) for the 2-station (Volume) trigger (9) IceAct (10) DMIce (11) InIce IC86 See DOMSetFunctions::InDOMSet_orig for details. Or the XML file. (Hopefully they match!)
-
I3MapKeyVectorIntPtr GetDefaultDOMSets()¶
-
namespace GTSUtils¶
Namespace that holds helper functions that the I3GlobalTriggerSim module uses.
class: GTSUtils
Version
- Rcs
date:
- Rcs
Copyright (c) 2006 IceCube Collaboration
- Author
olivas
Functions
-
I3TriggerStatus::Subdetector KeyToSubDetector(const TriggerKey &k)¶
-
void DumpChildren(const I3TriggerHierarchy&, I3TriggerHierarchy::iterator, std::stringstream&)¶
-
std::string Stringize(const I3TriggerHierarchy&)¶
-
namespace I3Units
-
namespace std
STL namespace.
-
namespace TimeShifterUtils¶
- file ClusterTriggerAlgorithm.cxx
-
#include <trigger-sim/algorithms/TimeWindow.h>#include <boost/foreach.hpp>#include <boost/assign/std/vector.hpp>
Copyright (C) 2006 the icecube collaboration $Id:
- Version
- Date
- Author
toale
- file ClusterTriggerAlgorithm.h
- #include “icetray/I3Logging.h”#include “trigger-sim/algorithms/TriggerService.h”#include “trigger-sim/algorithms/TriggerHit.h”
Copyright (C) 2006 the icecube collaboration $Id:
- Version
- Date
- Author
toale
- file CylinderTriggerAlgorithm.cxx
-
#include <boost/foreach.hpp>
Copyright (C) 2010 the icecube collaboration $Id:
- Version
- Date
- Author
danninger
- file CylinderTriggerAlgorithm.h
- #include “icetray/I3Logging.h”#include “trigger-sim/algorithms/TriggerService.h”#include “trigger-sim/algorithms/TriggerHit.h”#include <dataclasses/geometry/I3Geometry.h>
Copyright (C) 2010 the icecube collaboration $Id:
- Version
- Date
- Author
danninger
- file DetectorStatusUtils.h
- #include <dataclasses/status/I3DetectorStatus.h>#include <dataclasses/status/I3TriggerStatus.h>#include <dataclasses/TriggerKey.h>#include <boost/foreach.hpp>#include <boost/optional.hpp>#include <boost/parameter/keyword.hpp>#include <boost/parameter/name.hpp>#include <boost/parameter/preprocessor.hpp>#include <iostream>#include <sstream>
- file DOMSetFunctions.cxx
- #include “trigger-sim/utilities/DOMSetFunctions.h”#include <vector>#include <algorithm>#include <boost/foreach.hpp>#include <boost/assign/std/vector.hpp>#include <boost/assign/list_of.hpp>
Copyright (C) 2004 the icecube collaboration
- Rcs
DOMSetFunctions.cxx 2009/08/04 21: olivas Exp
- Version
- Rcs
1.1.1.1
- Date
- Rcs
2004/11/06 21:00:46
- Author
olivas
Functions
- const int LOWEST_DEEPCORE_STRING (79)
- const int DOMSET_2_HIGH_STRING (78)
- const int DOMSET_2_LOW_STRING (1)
- const unsigned DOMSET_2_HIGH_DOM (60)
- const unsigned DOMSET_2_LOW_DOM (1)
- const unsigned DOMSET_4_DC_HIGH_DOM (60)
- const unsigned DOMSET_4_DC_LOW_DOM (11)
- const unsigned DOMSET_4_II_HIGH_DOM (60)
- const unsigned DOMSET_4_II_LOW_DOM (41)
- const unsigned DOMSET_5_DC_HIGH_DOM (60)
- const unsigned DOMSET_5_DC_LOW_DOM (11)
- const unsigned DOMSET_5_II_HIGH_DOM (60)
- const unsigned DOMSET_5_II_LOW_DOM (39)
- const unsigned int DOMSET_6_DC_HIGH_DOM (60)
- const unsigned int DOMSET_6_DC_LOW_DOM (11)
- const unsigned int DOMSET_6_II_HIGH_DOM (60)
- const unsigned int DOMSET_6_II_LOW_DOM (39)
- file DOMSetFunctions.h
- #include <icetray/OMKey.h>#include <dataclasses/I3Map.h>#include <boost/assign/list_of.hpp>
Copyright (C) 2004 the icecube collaboration
- Rcs
DOMSetFunctions.h 2009/08/04 21: olivas Exp
- Version
- Rcs
1.1.1.1
- Date
- Rcs
2004/11/06 21:00:46
- Author
olivas
- file FaintParticleTriggerAlgorithm.cxx
-
#include <trigger-sim/algorithms/FPTTimeWindow.h>#include <boost/foreach.hpp>#include <boost/assign/std/vector.hpp>#include “trigger-sim/algorithms/TriggerHit.h”#include “dataclasses/geometry/I3Geometry.h”#include <dataclasses/I3Direction.h>#include <icetray/I3Units.h>
- Version
- Date
2023
- Author
tstuerwald
- file FaintParticleTriggerAlgorithm.h
- #include “icetray/I3Logging.h”#include “trigger-sim/algorithms/TriggerHit.h”#include “dataclasses/geometry/I3Geometry.h”#include “trigger-sim/algorithms/TriggerService.h”
- Version
- Date
- Author
toale
- file FPTTimeWindow.cxx
- #include “trigger-sim/algorithms/FPTTimeWindow.h”
- file FPTTimeWindow.h
- #include <string>#include “trigger-sim/algorithms/TriggerHit.h”#include “icetray/I3Logging.h”
- file GlobalTriggerSim.cxx
- #include “trigger-sim/algorithms/GlobalTriggerSim.h”#include “trigger-sim/utilities/GTSUtils.h”#include “dataclasses/physics/I3Trigger.h”#include “dataclasses/status/I3DetectorStatus.h”#include <dataclasses/I3Time.h>#include “trigger-sim/utilities/ReadoutWindowUtil.h”#include <boost/foreach.hpp>
Typedefs
-
typedef std::map<I3TriggerStatus::Subdetector, I3TriggerReadoutConfig> roconfigmap_t¶
-
typedef std::map<I3TriggerStatus::Subdetector, I3TriggerReadoutConfig> roconfigmap_t¶
- file GlobalTriggerSim.h
- #include “dataclasses/status/I3TriggerStatus.h”#include “dataclasses/status/I3DetectorStatus.h”#include “dataclasses/TriggerKey.h”#include “dataclasses/physics/I3TriggerHierarchy.h”#include <vector>#include “trigger-sim/utilities/ReadoutWindowUtil.h”
- file GTSUtils.cxx
- #include “trigger-sim/utilities/GTSUtils.h”#include “dataclasses/physics/I3Trigger.h”#include <boost/foreach.hpp>
class: GTSUtils
Version
- Rcs
Copyright (c) 2009 IceCube Collaboration
- Date
- Rcs
- Author
olivas
- file GTSUtils.h
- #include “dataclasses/status/I3TriggerStatus.h”#include “dataclasses/TriggerKey.h”#include “dataclasses/physics/I3TriggerHierarchy.h”#include <vector>
- file I3GlobalTriggerSim.cxx
- #include <limits>#include <boost/foreach.hpp>#include <trigger-sim/modules/I3GlobalTriggerSim.h>#include <trigger-sim/algorithms/GlobalTriggerSim.h>#include <icetray/I3TrayHeaders.h>#include <icetray/I3Module.h>#include <icetray/I3Units.h>#include <dataclasses/physics/I3EventHeader.h>#include <dataclasses/calibration/I3Calibration.h>#include <dataclasses/status/I3DetectorStatus.h>#include <dataclasses/status/I3TriggerStatus.h>#include <dataclasses/I3Time.h>#include <dataclasses/I3Double.h>#include <icetray/I3Bool.h>#include <dataclasses/physics/I3Trigger.h>#include <dataclasses/physics/I3TriggerHierarchy.h>#include <dataclasses/physics/I3DOMLaunch.h>#include <dataclasses/physics/I3RecoPulse.h>#include <phys-services/I3RandomService.h>
class: I3GlobalTriggerSim
Version
- Rcs
Copyright (c) 2009 IceCube Collaboration
class: I3GlobalTriggerSim
- Date
- Rcs
- Author
olivas
Version
- Rcs
Copyright (c) 2022 IceCube Collaboration
- Date
- Rcs
- Author
mlarson
Functions
-
I3_MODULE(I3GlobalTriggerSim)¶
- file I3GlobalTriggerSim.h
- #include <string>#include <boost/optional.hpp>#include <icetray/I3Module.h>#include <dataclasses/I3Time.h>
- file I3Pruner.cxx
- #include <iostream>#include “icetray/I3TrayHeaders.h”#include “dataclasses/geometry/I3Geometry.h”#include “dataclasses/physics/I3Trigger.h”#include “dataclasses/physics/I3TriggerHierarchy.h”#include “dataclasses/physics/I3DOMLaunch.h”#include “dataclasses/TriggerKey.h”#include “icetray/I3Units.h”#include “trigger-sim/modules/I3Pruner.h”#include “trigger-sim/utilities/ReadoutWindowUtil.h”#include <boost/foreach.hpp>
Copyright (C) 2009 the icecube collaboration
- file I3Pruner.h
- #include “icetray/I3ConditionalModule.h”
Copyright (C) 2004 the icecube collaboration
- Rcs
I3Pruner.h,v 1.1.1.1 2004/11/06 21:00:46 bouchta Exp
- Version
- Rcs
1.1.1.1
- Date
- Rcs
2004/11/06 21:00:46
- Author
olivas
- file I3TriggerSimModule.cxx
- #include <boost/foreach.hpp>#include <trigger-sim/modules/I3TriggerSimModule.h>
Functions
-
I3_MODULE(I3TriggerSimModule)¶
-
I3_MODULE(I3TriggerSimModule)¶
- file I3TriggerSimModule.h
- #include <icetray/I3Context.h>#include <icetray/I3Frame.h>#include <icetray/I3Module.h>#include <icetray/I3ConditionalModule.h>#include <dataclasses/I3Map.h>#include <dataclasses/physics/I3DOMLaunch.h>#include <dataclasses/physics/I3RecoPulse.h>#include <dataclasses/geometry/I3Geometry.h>#include <dataclasses/status/I3DetectorStatus.h>#include <dataclasses/TriggerKey.h>#include <dataclasses/physics/I3Trigger.h>#include <dataclasses/physics/I3TriggerHierarchy.h>#include <dataclasses/status/I3TriggerStatus.h>#include <trigger-sim/algorithms/TriggerHit.h>#include <trigger-sim/algorithms/TriggerService.h>
Copyright (C) 2022 the icecube collaboration $Id:
- Version
- Date
- Author
mlarson
- file ReadoutWindowUtil.cxx
- #include “trigger-sim/utilities/ReadoutWindowUtil.h”#include “dataclasses/TriggerKey.h”
- file ReadoutWindowUtil.h
- #include “dataclasses/physics/I3Trigger.h”#include “dataclasses/TriggerKey.h”#include “dataclasses/status/I3DetectorStatus.h”#include “dataclasses/status/I3TriggerStatus.h”
- file SimpleMajorityTriggerAlgorithm.cxx
-
#include <trigger-sim/algorithms/TimeWindow.h>#include <boost/foreach.hpp>#include <boost/assign/std/vector.hpp>
Copyright (C) 2006 the icecube collaboration $Id:
- Version
- Date
- Author
toale
- file SimpleMajorityTriggerAlgorithm.h
- #include “icetray/I3Logging.h”#include “trigger-sim/algorithms/TriggerService.h”#include “trigger-sim/algorithms/TriggerHit.h”
Copyright (C) 2006 the icecube collaboration $Id:
- Version
- Date
- Author
toale
- file SlowMonopoleTriggerAlgorithm.cxx
- #include <iostream>#include <math.h>#include <boost/foreach.hpp>#include “icetray/I3TrayHeaders.h”#include “icetray/I3Units.h”#include “icetray/I3Bool.h”#include “dataclasses/physics/I3TriggerHierarchy.h”#include “dataclasses/status/I3TriggerStatus.h”#include “dataclasses/status/I3DetectorStatus.h”#include “dataclasses/geometry/I3Geometry.h”#include “dataclasses/physics/I3DOMLaunch.h”#include “dataclasses/TriggerKey.h”#include “dataclasses/I3Constants.h”#include “dataclasses/I3Double.h”#include “dataclasses/I3Map.h”#include “dataclasses/I3Vector.h”#include “trigger-sim/utilities/DOMSetFunctions.h”#include “trigger-sim/algorithms/TriggerContainer.h”#include “trigger-sim/algorithms/TriggerHit.h”#include “trigger-sim/utilities/DetectorStatusUtils.h”
- file SlowMonopoleTriggerAlgorithm.h
- #include “dataclasses/physics/I3Trigger.h”#include “dataclasses/physics/I3Particle.h”#include <dataclasses/geometry/I3Geometry.h>#include “trigger-sim/algorithms/TriggerHit.h”#include “trigger-sim/algorithms/TriggerContainer.h”#include “trigger-sim/algorithms/TriggerService.h”
- file SlowMPHit.h
- #include <vector>#include <list>#include “dataclasses/I3Map.h”
Typedefs
-
typedef std::vector<SlowMPHitVector> SlowMPHitVectorVector¶
-
typedef I3Map<int, SlowMPHitVector> IntSlowMPHitVectorMap¶
-
typedef std::pair<SlowMPHitVector::const_iterator, SlowMPHitVector::const_iterator> SlowMPHitIterPair¶
-
typedef std::vector<SlowMPHitIterPair> SlowMPHitIterPairVector¶
Functions
-
I3_POINTER_TYPEDEFS(SlowMPHitVector)¶
-
I3_POINTER_TYPEDEFS(SlowMPHitVectorVector)¶
-
I3_POINTER_TYPEDEFS(IntSlowMPHitVectorMap)¶
-
I3_POINTER_TYPEDEFS(SlowMPHitIntMap)¶
-
I3_POINTER_TYPEDEFS(SlowMPHitIterPair)¶
-
I3_POINTER_TYPEDEFS(SlowMPHitIterPairVector)¶
-
I3_POINTER_TYPEDEFS(SlowMPHitList)¶
-
typedef std::vector<SlowMPHitVector> SlowMPHitVectorVector¶
- file TimeShifterUtils.cxx
- #include <algorithm>#include <vector>#include <boost/foreach.hpp>#include <trigger-sim/utilities/TimeShifterUtils.h>#include <icetray/I3Frame.h>#include <dataclasses/physics/I3MCHit.h>#include “dataclasses/physics/I3MCTree.h”#include “simclasses/I3MMCTrack.h”#include “dataclasses/physics/I3Trigger.h”#include “dataclasses/physics/I3TriggerHierarchy.h”#include “dataclasses/physics/I3DOMLaunch.h”#include “dataclasses/physics/I3RecoPulse.h”#include “dataclasses/I3Double.h”#include “dataclasses/physics/I3FlasherInfo.h”#include “simclasses/I3MCPE.h”#include “simclasses/I3MCPulse.h”
Functions
-
I3MCHitSeriesMapPtr Shift(I3MCHitSeriesMapConstPtr map, double dt)¶
-
I3MCTreePtr Shift(I3MCTreeConstPtr tree, double dt)¶
-
I3MMCTrackListPtr Shift(I3MMCTrackListConstPtr tracks, double dt)¶
-
I3DOMLaunchSeriesMapPtr Shift(I3DOMLaunchSeriesMapConstPtr launches, double dt)¶
-
I3RecoPulseSeriesMapPtr Shift(I3RecoPulseSeriesMapConstPtr pulses, double dt)¶
-
I3TriggerHierarchyPtr Shift(I3TriggerHierarchyConstPtr tree, double dt)¶
-
I3VectorI3TriggerPtr Shift(I3VectorI3TriggerConstPtr triggers, double dt)¶
-
I3DoublePtr Shift(I3DoubleConstPtr i3double, double dt)¶
-
I3ParticlePtr Shift(I3ParticleConstPtr i3particle, double dt)¶
-
I3FlasherInfoVectPtr Shift(I3FlasherInfoVectConstPtr oldflasher, double dt)¶
-
I3MCPESeriesMapPtr Shift(I3MCPESeriesMapConstPtr map, double dt)¶
-
I3MCPulseSeriesMapPtr Shift(I3MCPulseSeriesMapConstPtr map, double dt)¶
-
I3VectorI3ParticlePtr Shift(I3VectorI3ParticleConstPtr vect, double dt)¶
-
I3MCHitSeriesMapPtr Shift(I3MCHitSeriesMapConstPtr map, double dt)¶
- file TimeShifterUtils.h
- #include <icetray/I3PointerTypedefs.h>
Functions
-
I3_FORWARD_DECLARATION(I3Frame)¶
class: TimeShifterUtils
Version
- Rcs
date:
- Rcs
Copyright (c) 2006 IceCube Collaboration
- Author
olivas
-
I3_FORWARD_DECLARATION(I3Frame)¶
- file TimeWindow.cxx
- #include “trigger-sim/algorithms/TimeWindow.h”#include <boost/foreach.hpp>
- file TimeWindow.h
- #include <string>#include “trigger-sim/algorithms/TriggerHit.h”#include “icetray/I3Logging.h”
- file TriggerContainer.h
- #include <dataclasses/physics/I3Trigger.h>#include <dataclasses/physics/I3Particle.h>
Typedefs
-
typedef std::vector<TriggerContainer> TriggerContainerVector¶
-
typedef std::vector<TriggerContainer> TriggerContainerVector¶
- file TriggerHit.h
- #include <vector>#include <list>#include “dataclasses/I3Map.h”
Typedefs
-
typedef std::vector<TriggerHit> TriggerHitVector¶
-
typedef std::vector<TriggerHitVector> TriggerHitVectorVector¶
-
typedef I3Map<int, TriggerHitVector> IntTriggerHitVectorMap¶
-
typedef I3Map<TriggerHit, int> TriggerHitIntMap¶
-
typedef std::pair<TriggerHitVector::const_iterator, TriggerHitVector::const_iterator> TriggerHitIterPair¶
-
typedef std::vector<TriggerHitIterPair> TriggerHitIterPairVector¶
-
typedef std::list<TriggerHit> TriggerHitList¶
Functions
-
I3_POINTER_TYPEDEFS(TriggerHit)¶
-
I3_POINTER_TYPEDEFS(TriggerHitVector)¶
-
I3_POINTER_TYPEDEFS(TriggerHitVectorVector)¶
-
I3_POINTER_TYPEDEFS(IntTriggerHitVectorMap)¶
-
I3_POINTER_TYPEDEFS(TriggerHitIntMap)¶
-
I3_POINTER_TYPEDEFS(TriggerHitIterPair)¶
-
I3_POINTER_TYPEDEFS(TriggerHitIterPairVector)¶
-
I3_POINTER_TYPEDEFS(TriggerHitList)¶
-
typedef std::vector<TriggerHit> TriggerHitVector¶
- file TriggerService.cxx
- #include “trigger-sim/algorithms/TriggerService.h”#include <algorithm>#include <boost/foreach.hpp>
Copyright (C) 2022 the icecube collaboration $Id:
- Version
- Date
- Author
mlarson
- file TriggerService.h
- #include “icetray/I3Logging.h”#include “icetray/OMKey.h”#include “dataclasses/physics/I3DOMLaunch.h”#include “dataclasses/physics/I3RecoPulse.h”#include “dataclasses/I3Map.h”#include “trigger-sim/utilities/DOMSetFunctions.h”#include “trigger-sim/algorithms/TriggerHit.h”
Copyright (C) 2006 the icecube collaboration $Id:
- Version
- Date
- Author
mlarson
- dir algorithms
- dir algorithms
- dir icetray
- dir modules
- dir modules
- dir private
- dir public
- dir trigger-sim
- dir trigger-sim
- dir trigger-sim
- dir utilities
- dir utilities