CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
Adjuster.h
Go to the documentation of this file.
1 #ifndef SimGeneral_MixingModule_Adjuster_h
2 #define SimGeneral_MixingModule_Adjuster_h
3 
15 
16 #include <memory>
17 #include <vector>
18 #include <string>
19 #include <iostream>
20 
21 class FastTrackerRecHit;
22 
23 namespace edm {
24  class ModuleCallingContext;
25 
26  class AdjusterBase {
27  public:
28  virtual ~AdjusterBase() {}
29  virtual void doOffset(int bunchspace,
30  int bcr,
31  const edm::EventPrincipal&,
32  ModuleCallingContext const*,
33  unsigned int EventNr,
34  int vertexOffset) = 0;
35  virtual bool checkSignal(edm::Event const& event) = 0;
36  };
37 
38  template <typename T>
39  class Adjuster : public AdjusterBase {
40  public:
41  Adjuster(InputTag const& tag, edm::ConsumesCollector&& iC, bool wrap);
42 
43  ~Adjuster() override {}
44 
45  void doOffset(int bunchspace,
46  int bcr,
47  const edm::EventPrincipal&,
48  ModuleCallingContext const*,
49  unsigned int EventNr,
50  int vertexOffset) override;
51 
52  bool checkSignal(edm::Event const& event) override {
53  bool got = false;
54  edm::Handle<T> result_t;
55  got = event.getByToken(token_, result_t);
56  return got;
57  }
58 
59  private:
61  bool WrapT_ = false;
63  };
64 
65  //==============================================================================
66  // implementations
67  //==============================================================================
68 
69  namespace detail {
70  void doTheOffset(int bunchspace,
71  int bcr,
72  std::vector<SimTrack>& product,
73  unsigned int eventNr,
74  int vertexOffset,
75  bool wraptimes);
76  void doTheOffset(int bunchspace,
77  int bcr,
78  std::vector<SimVertex>& product,
79  unsigned int eventNr,
80  int vertexOffset,
81  bool wraptimes);
82  void doTheOffset(int bunchspace,
83  int bcr,
84  std::vector<PCaloHit>& product,
85  unsigned int eventNr,
86  int vertexOffset,
87  bool wraptimes);
88  void doTheOffset(
89  int bunchspace, int bcr, std::vector<PSimHit>& product, unsigned int eventNr, int vertexOffset, bool wraptimes);
90  void doTheOffset(int bunchspace,
91  int bcr,
92  TrackingRecHitCollection& product,
93  unsigned int eventNr,
94  int vertexOffset,
95  bool wraptimes);
96  } // namespace detail
97 
98  template <typename T>
100  int bcr,
101  const EventPrincipal& ep,
102  ModuleCallingContext const* mcc,
103  unsigned int eventNr,
104  int vertexOffset) {
105  std::shared_ptr<Wrapper<T> const> shPtr = getProductByTag<T>(ep, tag_, mcc);
106  if (shPtr) {
107  T& product = const_cast<T&>(*shPtr->product());
108  detail::doTheOffset(bunchspace, bcr, product, eventNr, vertexOffset, WrapT_);
109  }
110  }
111 
112  template <typename T>
113  Adjuster<T>::Adjuster(InputTag const& tag, ConsumesCollector&& iC, bool wrapLongTimes)
114  : tag_(tag), token_(iC.consumes<T>(tag)) {
115  if (wrapLongTimes) {
116  std::string Musearch = tag_.instance();
117  if (Musearch.find("Muon") == 0)
118  WrapT_ = true; // wrap time for neutrons in Muon system subdetectors
119  }
120  }
121 } // namespace edm
122 
123 #endif
void doOffset(int bunchspace, int bcr, const edm::EventPrincipal &, ModuleCallingContext const *, unsigned int EventNr, int vertexOffset) override
Definition: Adjuster.h:99
tuple bunchspace
in terms of 25 nsec
~Adjuster() override
Definition: Adjuster.h:43
virtual void doOffset(int bunchspace, int bcr, const edm::EventPrincipal &, ModuleCallingContext const *, unsigned int EventNr, int vertexOffset)=0
virtual ~AdjusterBase()
Definition: Adjuster.h:28
void doTheOffset(int bunchSpace, int bcr, std::vector< SimTrack > &simtracks, unsigned int evtNr, int vertexOffset, bool wrap)
Definition: Adjuster.cc:7
Adjuster(InputTag const &tag, edm::ConsumesCollector &&iC, bool wrap)
Definition: Adjuster.h:113
bool WrapT_
Definition: Adjuster.h:61
bool checkSignal(edm::Event const &event) override
Definition: Adjuster.h:52
EDGetTokenT< T > token_
Definition: Adjuster.h:62
virtual bool checkSignal(edm::Event const &event)=0
auto wrap(F iFunc) -> decltype(iFunc())
long double T
std::string const & instance() const
Definition: InputTag.h:37
InputTag tag_
Definition: Adjuster.h:60