CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
DBremWatcher Class Reference
Inheritance diagram for DBremWatcher:
SimProducer Observer< const BeginOfTrack *> Observer< const BeginOfEvent *> Observer< const BeginOfRun *> Observer< const EndOfEvent *> Observer< const EndOfTrack *> SimWatcher

Public Member Functions

 DBremWatcher (edm::ParameterSet const &p)
 
void produce (edm::Event &, const edm::EventSetup &) override
 
void update (const BeginOfTrack *trk) override
 This routine will be called when the appropriate signal arrives. More...
 
void update (const BeginOfEvent *event) override
 This routine will be called when the appropriate signal arrives. More...
 
void update (const EndOfEvent *event) override
 This routine will be called when the appropriate signal arrives. More...
 
void update (const BeginOfRun *run) override
 This routine will be called when the appropriate signal arrives. More...
 
void update (const EndOfTrack *trk) override
 This routine will be called when the appropriate signal arrives. More...
 
 ~DBremWatcher () override=default
 
- Public Member Functions inherited from SimProducer
const SimProduceroperator= (const SimProducer &)=delete
 
void registerProducts (edm::ProducesCollector producesCollector)
 
 SimProducer ()
 
 SimProducer (const SimProducer &)=delete
 
- Public Member Functions inherited from SimWatcher
virtual void beginRun (edm::EventSetup const &)
 
bool isMT () const
 
const SimWatcheroperator= (const SimWatcher &)=delete
 
virtual void registerConsumes (edm::ConsumesCollector)
 
 SimWatcher ()
 
 SimWatcher (const SimWatcher &)=delete
 
virtual ~SimWatcher ()
 
- Public Member Functions inherited from Observer< const BeginOfTrack *>
 Observer ()
 
void slotForUpdate (const BeginOfTrack * iT)
 
virtual ~Observer ()
 
- Public Member Functions inherited from Observer< const BeginOfEvent *>
 Observer ()
 
void slotForUpdate (const BeginOfEvent * iT)
 
virtual ~Observer ()
 
- Public Member Functions inherited from Observer< const BeginOfRun *>
 Observer ()
 
void slotForUpdate (const BeginOfRun * iT)
 
virtual ~Observer ()
 
- Public Member Functions inherited from Observer< const EndOfEvent *>
 Observer ()
 
void slotForUpdate (const EndOfEvent * iT)
 
virtual ~Observer ()
 
- Public Member Functions inherited from Observer< const EndOfTrack *>
 Observer ()
 
void slotForUpdate (const EndOfTrack * iT)
 
virtual ~Observer ()
 

Private Attributes

G4ThreeVector aPrimeTraj
 
double biasFactor
 
double f_energy
 
G4ThreeVector finaltraj
 
bool foundbrem
 
float m_weight
 
int MotherId
 
std::vector< int > pdgs_
 
G4ThreeVector VertexPos
 

Additional Inherited Members

- Protected Member Functions inherited from SimProducer
template<class T >
void produces ()
 
template<class T >
void produces (const std::string &instanceName)
 
- Protected Member Functions inherited from SimWatcher
void setMT (bool val)
 

Detailed Description

Definition at line 32 of file DBremWatcher.cc.

Constructor & Destructor Documentation

◆ DBremWatcher()

DBremWatcher::DBremWatcher ( edm::ParameterSet const &  p)

Definition at line 60 of file DBremWatcher.cc.

References aPrimeTraj, biasFactor, f_energy, finaltraj, foundbrem, edm::ParameterSet::getUntrackedParameter(), dqmdumpme::k, m_weight, AlCaHLTBitMon_ParallelJobs::p, pdgs_, and VertexPos.

60  {
61  edm::ParameterSet ps = p.getParameter<edm::ParameterSet>("DBremWatcher");
62  biasFactor = ps.getUntrackedParameter<double>("DBremBiasFactor", 1);
63  m_weight = 0;
64  foundbrem = false;
65  finaltraj = G4ThreeVector(0, 0, 0);
66  aPrimeTraj = G4ThreeVector(0, 0, 0);
67  VertexPos = G4ThreeVector(0, 0, 0);
68  f_energy = 0;
69  produces<float>("DBremEventWeight");
70  produces<float>("DBremLocationX");
71  produces<float>("DBremLocationY");
72  produces<float>("DBremLocationZ");
73  produces<float>("DBremAngle");
74  produces<float>("DBremInitialEnergy");
75  produces<float>("DBremFinalEnergy");
76  produces<float>("BiasFactor");
77  pdgs_ = ps.getUntrackedParameter<std::vector<int>>("PDGCodes");
78  edm::LogVerbatim("DBremWatcher") << "DBremWatcher:: Save Sim Track if PDG code "
79  << "is one from the list of " << pdgs_.size() << " items";
80  for (unsigned int k = 0; k < pdgs_.size(); ++k)
81  edm::LogVerbatim("DBremWatcher") << "[" << k << "] " << pdgs_[k];
82 }
Log< level::Info, true > LogVerbatim
double biasFactor
Definition: DBremWatcher.cc:52
double f_energy
Definition: DBremWatcher.cc:57
G4ThreeVector aPrimeTraj
Definition: DBremWatcher.cc:54
T getUntrackedParameter(std::string const &, T const &) const
G4ThreeVector VertexPos
Definition: DBremWatcher.cc:56
G4ThreeVector finaltraj
Definition: DBremWatcher.cc:55
std::vector< int > pdgs_
Definition: DBremWatcher.cc:49

◆ ~DBremWatcher()

DBremWatcher::~DBremWatcher ( )
overridedefault

Member Function Documentation

◆ produce()

void DBremWatcher::produce ( edm::Event fEvent,
const edm::EventSetup  
)
overridevirtual

Implements SimProducer.

Definition at line 145 of file DBremWatcher.cc.

References aPrimeTraj, biasFactor, f_energy, hcaldqm::fEvent, finaltraj, foundbrem, m_weight, eostools::move(), funct::pow(), mathSSE::sqrt(), and VertexPos.

145  {
146  if (foundbrem) {
147  std::unique_ptr<float> weight = std::make_unique<float>(m_weight);
148  fEvent.put(std::move(weight), "DBremEventWeight");
149  std::unique_ptr<float> vtxposx = std::make_unique<float>(VertexPos.x());
150  std::unique_ptr<float> vtxposy = std::make_unique<float>(VertexPos.y());
151  std::unique_ptr<float> vtxposz = std::make_unique<float>(VertexPos.z());
152  fEvent.put(std::move(vtxposx), "DBremLocationX");
153  fEvent.put(std::move(vtxposy), "DBremLocationY");
154  fEvent.put(std::move(vtxposz), "DBremLocationZ");
155  std::unique_ptr<float> finalE = std::make_unique<float>(f_energy / CLHEP::GeV);
156  fEvent.put(std::move(finalE), "DBremFinalEnergy");
157  float deflectionAngle = -1;
158  float initialEnergy = sqrt(pow(finaltraj.x() + aPrimeTraj.x(), 2) + pow(finaltraj.y() + aPrimeTraj.y(), 2) +
159  pow(finaltraj.z() + aPrimeTraj.z(), 2) + pow(0.1056 * CLHEP::GeV, 2));
160  G4ThreeVector mother(
161  finaltraj.x() + aPrimeTraj.x(), finaltraj.y() + aPrimeTraj.y(), finaltraj.z() + aPrimeTraj.z());
162  deflectionAngle = mother.angle(finaltraj);
163  std::unique_ptr<float> dAngle = std::make_unique<float>(deflectionAngle);
164  std::unique_ptr<float> initialE = std::make_unique<float>(initialEnergy / CLHEP::GeV);
165  fEvent.put(std::move(dAngle), "DBremAngle");
166  fEvent.put(std::move(initialE), "DBremInitialEnergy");
167  std::unique_ptr<float> bias = std::make_unique<float>(biasFactor);
168  fEvent.put(std::move(bias), "BiasFactor");
169  } else {
170  std::unique_ptr<float> weight = std::make_unique<float>(0.);
171  fEvent.put(std::move(weight), "DBremEventWeight");
172  }
173 }
double biasFactor
Definition: DBremWatcher.cc:52
double f_energy
Definition: DBremWatcher.cc:57
G4ThreeVector aPrimeTraj
Definition: DBremWatcher.cc:54
Definition: weight.py:1
G4ThreeVector VertexPos
Definition: DBremWatcher.cc:56
T sqrt(T t)
Definition: SSEVec.h:23
G4ThreeVector finaltraj
Definition: DBremWatcher.cc:55
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
def move(src, dest)
Definition: eostools.py:511

◆ update() [1/5]

void DBremWatcher::update ( const BeginOfTrack )
overridevirtual

This routine will be called when the appropriate signal arrives.

Implements Observer< const BeginOfTrack *>.

Definition at line 84 of file DBremWatcher.cc.

References aPrimeTraj, f_energy, finaltraj, spr::find(), foundbrem, LogDebug, m_weight, dileptonTrigSettings_cff::pdg, pdgs_, TrackInformation::setStoreTrack(), TrackInformation::storeTrack(), and VertexPos.

Referenced by progressbar.ProgressBar::__next__(), MatrixUtil.Matrix::__setitem__(), MatrixUtil.Steps::__setitem__(), progressbar.ProgressBar::finish(), and MatrixUtil.Steps::overwrite().

84  {
85  G4Track* theTrack = (G4Track*)((*trk)());
86  TrackInformation* trkInfo = (TrackInformation*)(theTrack->GetUserInformation());
87  if (trkInfo) {
88  int pdg = theTrack->GetDefinition()->GetPDGEncoding();
89  G4ThreeVector Vpos = theTrack->GetVertexPosition();
90  const G4VProcess* TrPro = theTrack->GetCreatorProcess();
91  if (TrPro != nullptr) {
92  if ((theTrack->GetCreatorProcess()->GetProcessName()) == "muDBrem") {
93  if (std::find(pdgs_.begin(), pdgs_.end(), pdg) == pdgs_.end()) {
94  //Found the deflected muon track
95  trkInfo->storeTrack();
96  if (!theTrack->IsGoodForTracking()) {
97  theTrack->SetGoodForTrackingFlag(true);
98  }
99  f_energy = theTrack->GetTotalEnergy();
100  foundbrem = true;
101  finaltraj = theTrack->GetMomentum();
102  } else {
103  m_weight = theTrack->GetWeight();
104  }
105  }
106  }
107  if (std::find(pdgs_.begin(), pdgs_.end(), pdg) != pdgs_.end()) {
108  //Found an A'
109  trkInfo->setStoreTrack();
110  VertexPos = Vpos;
111  aPrimeTraj = theTrack->GetMomentum();
112  LogDebug("DBremWatcher") << "Save SimTrack the Track " << theTrack->GetTrackID() << " Type "
113  << theTrack->GetDefinition()->GetParticleName() << " Momentum "
114  << theTrack->GetMomentum() / CLHEP::MeV << " MeV/c";
115  }
116  }
117 }
bool storeTrack() const
void setStoreTrack()
can only be set to true, cannot be reset to false!
double f_energy
Definition: DBremWatcher.cc:57
G4ThreeVector aPrimeTraj
Definition: DBremWatcher.cc:54
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
G4ThreeVector VertexPos
Definition: DBremWatcher.cc:56
G4ThreeVector finaltraj
Definition: DBremWatcher.cc:55
std::vector< int > pdgs_
Definition: DBremWatcher.cc:49
#define LogDebug(id)

◆ update() [2/5]

void DBremWatcher::update ( const BeginOfEvent )
overridevirtual

This routine will be called when the appropriate signal arrives.

Implements Observer< const BeginOfEvent *>.

Definition at line 121 of file DBremWatcher.cc.

References foundbrem, and unpackData-CaloStage2::pname.

Referenced by progressbar.ProgressBar::__next__(), MatrixUtil.Matrix::__setitem__(), MatrixUtil.Steps::__setitem__(), progressbar.ProgressBar::finish(), and MatrixUtil.Steps::overwrite().

121  {
122  G4String pname = "muDBrem";
123  G4ProcessTable* ptable = G4ProcessTable::GetProcessTable();
124  G4bool state = true;
125  ptable->SetProcessActivation(pname, state);
126  foundbrem = false;
127 }

◆ update() [3/5]

void DBremWatcher::update ( const EndOfEvent )
overridevirtual

This routine will be called when the appropriate signal arrives.

Implements Observer< const EndOfEvent *>.

Definition at line 129 of file DBremWatcher.cc.

Referenced by progressbar.ProgressBar::__next__(), MatrixUtil.Matrix::__setitem__(), MatrixUtil.Steps::__setitem__(), progressbar.ProgressBar::finish(), and MatrixUtil.Steps::overwrite().

129 {}

◆ update() [4/5]

void DBremWatcher::update ( const BeginOfRun )
overridevirtual

This routine will be called when the appropriate signal arrives.

Implements Observer< const BeginOfRun *>.

Definition at line 119 of file DBremWatcher.cc.

Referenced by progressbar.ProgressBar::__next__(), MatrixUtil.Matrix::__setitem__(), MatrixUtil.Steps::__setitem__(), progressbar.ProgressBar::finish(), and MatrixUtil.Steps::overwrite().

119 {}

◆ update() [5/5]

void DBremWatcher::update ( const EndOfTrack )
overridevirtual

This routine will be called when the appropriate signal arrives.

Implements Observer< const EndOfTrack *>.

Definition at line 131 of file DBremWatcher.cc.

References spr::find(), dileptonTrigSettings_cff::pdg, pdgs_, and TrackInformation::setStoreTrack().

Referenced by progressbar.ProgressBar::__next__(), MatrixUtil.Matrix::__setitem__(), MatrixUtil.Steps::__setitem__(), progressbar.ProgressBar::finish(), and MatrixUtil.Steps::overwrite().

131  {
132  G4Track* theTrack = (G4Track*)((*trk)());
133  TrackInformation* trkInfo = (TrackInformation*)(theTrack->GetUserInformation());
134  const G4VProcess* TrPro = theTrack->GetCreatorProcess();
135  if (trkInfo && TrPro != nullptr) {
136  int pdg = theTrack->GetDefinition()->GetPDGEncoding();
137 
138  if (std::find(pdgs_.begin(), pdgs_.end(), pdg) == pdgs_.end() &&
139  (theTrack->GetCreatorProcess()->GetProcessName()) == "muDBrem") {
140  trkInfo->setStoreTrack();
141  }
142  }
143 }
void setStoreTrack()
can only be set to true, cannot be reset to false!
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
std::vector< int > pdgs_
Definition: DBremWatcher.cc:49

Member Data Documentation

◆ aPrimeTraj

G4ThreeVector DBremWatcher::aPrimeTraj
private

Definition at line 54 of file DBremWatcher.cc.

Referenced by DBremWatcher(), produce(), and update().

◆ biasFactor

double DBremWatcher::biasFactor
private

Definition at line 52 of file DBremWatcher.cc.

Referenced by DBremWatcher(), and produce().

◆ f_energy

double DBremWatcher::f_energy
private

Definition at line 57 of file DBremWatcher.cc.

Referenced by DBremWatcher(), produce(), and update().

◆ finaltraj

G4ThreeVector DBremWatcher::finaltraj
private

Definition at line 55 of file DBremWatcher.cc.

Referenced by DBremWatcher(), produce(), and update().

◆ foundbrem

bool DBremWatcher::foundbrem
private

Definition at line 53 of file DBremWatcher.cc.

Referenced by DBremWatcher(), produce(), and update().

◆ m_weight

float DBremWatcher::m_weight
private

Definition at line 51 of file DBremWatcher.cc.

Referenced by DBremWatcher(), produce(), and update().

◆ MotherId

int DBremWatcher::MotherId
private

Definition at line 50 of file DBremWatcher.cc.

◆ pdgs_

std::vector<int> DBremWatcher::pdgs_
private

Definition at line 49 of file DBremWatcher.cc.

Referenced by DBremWatcher(), and update().

◆ VertexPos

G4ThreeVector DBremWatcher::VertexPos
private

Definition at line 56 of file DBremWatcher.cc.

Referenced by DBremWatcher(), produce(), and update().