Back to home page

sPhenix code displayed by LXR

 
 

    


File indexing completed on 2025-08-05 08:18:07

0001 // Tell emacs that this is a C++ source
0002 //  -*- C++ -*-.
0003 #ifndef G4MAIN_FUN4ALLDSTPILEUPINPUTMANAGER_H
0004 #define G4MAIN_FUN4ALLDSTPILEUPINPUTMANAGER_H
0005 
0006 /*!
0007  * \file Fun4AllDstPileupInputManager.h
0008  * \author Hugo Pereira Da Costa <hugo.pereira-da-costa@cea.fr>
0009  */
0010 
0011 #include <fun4all/Fun4AllInputManager.h>
0012 #include <fun4all/Fun4AllReturnCodes.h>  // for SYNC_NOOBJECT, SYNC_OK
0013 
0014 #include <phool/PHCompositeNode.h>  // for PHCompositeNode
0015 #include <phool/PHNodeIOManager.h>  // for PHNodeIOManager
0016 #include <phool/sphenix_constants.h>
0017 
0018 #include <gsl/gsl_rng.h>
0019 
0020 #include <map>
0021 #include <memory>
0022 #include <string>
0023 #include <utility>  // for pair
0024 
0025 class SyncObject;
0026 
0027 /*!
0028  * dedicated input manager that merges single events into "merged" events, containing a trigger event
0029  * and a number of time-shifted pile-up events corresponding to a given pile-up rate
0030  */
0031 class Fun4AllDstPileupInputManager : public Fun4AllInputManager
0032 {
0033  public:
0034   Fun4AllDstPileupInputManager(const std::string &name = "DUMMY", const std::string &nodename = "DST", const std::string &topnodename = "TOP");
0035   int fileopen(const std::string &filenam) override;
0036   int fileclose() override;
0037   int run(const int nevents = 0) override;
0038   int BranchSelect(const std::string &branch, const int iflag) override;
0039   int setBranches() override;
0040   void Print(const std::string &what = "ALL") const override;
0041   int PushBackEvents(const int i) override;
0042 
0043   // Effectivly turn off the synchronization checking (copy from Fun4AllNoSyncDstInputManager)
0044   int SyncIt(const SyncObject * /*mastersync*/) override { return Fun4AllReturnCodes::SYNC_OK; }
0045   int GetSyncObject(SyncObject ** /*mastersync*/) override { return Fun4AllReturnCodes::SYNC_NOOBJECT; }
0046   int NoSyncPushBackEvents(const int nevt) override { return PushBackEvents(nevt); }
0047 
0048   /// collision rate in Hz
0049   void setCollisionRate(double Hz)
0050   {
0051     m_collision_rate = Hz;
0052   }
0053 
0054   /// time between bunch crossing in ns
0055   void setTimeBetweenCrossings(double nsec)
0056   {
0057     m_time_between_crossings = nsec;
0058   }
0059 
0060   //! set time window for pileup events (ns)
0061   void setPileupTimeWindow(double tmin, double tmax)
0062   {
0063     m_tmin = tmin;
0064     m_tmax = tmax;
0065   }
0066   //! for symmetric windows
0067   void setDetectorActiveCrossings(const std::string &name, const int nbcross);
0068 
0069   void setDetectorActiveCrossings(const std::string &name, const int min, const int max);
0070 
0071  private:
0072   //! loads one event on internal DST node
0073   int runOne(const int nevents = 0);
0074 
0075   //!@name event counters
0076   //@{
0077   bool m_ReadRunTTree = true;
0078   int m_ievent_total = 0;
0079   int m_ievent_thisfile = 0;
0080   //@}
0081 
0082   std::string m_fullfilename;
0083   std::string m_RunNode = "RUN";
0084   std::map<const std::string, int> m_branchread;
0085 
0086   //! dst node from TopNode
0087   PHCompositeNode *m_dstNode = nullptr;
0088 
0089   //! run node from TopNode
0090   PHCompositeNode *m_runNode = nullptr;
0091 
0092   //! internal dst node to copy background events
0093   std::unique_ptr<PHCompositeNode> m_dstNodeInternal;
0094 
0095   //!@name internal runnodes to perform the summation over multiple runs
0096   //@{
0097   std::unique_ptr<PHCompositeNode> m_runNodeCopy;
0098   std::unique_ptr<PHCompositeNode> m_runNodeSum;
0099   //@}
0100 
0101   //! input manager for active (trigger) events
0102   /*! corresponding nodes are copied directly to the topNode */
0103   std::unique_ptr<PHNodeIOManager> m_IManager;
0104 
0105   //! time between crossings. This is a RHIC constant (ns)
0106   double m_time_between_crossings = sphenix_constants::time_between_crossings;
0107 
0108   //! collision rate (Hz)
0109   double m_collision_rate = 5e4;
0110 
0111   //! min integration time for pileup in the TPC (ns)
0112   double m_tmin = -13500;
0113 
0114   //! max integration time for pileup in the TPC (ns)
0115   double m_tmax = 13500;
0116 
0117   //! random generator
0118   class Deleter
0119   {
0120    public:
0121     void operator()(gsl_rng *rng) const { gsl_rng_free(rng); }
0122   };
0123 
0124   std::unique_ptr<gsl_rng, Deleter> m_rng;
0125 
0126   std::map<std::string, std::pair<double, double>> m_DetectorTiming;
0127 };
0128 
0129 #endif /* __Fun4AllDstPileupInputManager_H__ */