CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
DBremWatcher Class Reference

#include <DBremWatcher.h>

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
 
- 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
 
float 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 17 of file DBremWatcher.h.

Constructor & Destructor Documentation

◆ DBremWatcher()

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

Definition at line 28 of file DBremWatcher.cc.

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

28  {
29  edm::ParameterSet ps = p.getParameter<edm::ParameterSet>("DBremWatcher");
30  biasFactor = ps.getUntrackedParameter<double>("DBremBiasFactor", 1);
31  m_weight = 0;
32  foundbrem = false;
33  finaltraj = G4ThreeVector(0, 0, 0);
34  aPrimeTraj = G4ThreeVector(0, 0, 0);
35  VertexPos = G4ThreeVector(0, 0, 0);
36  f_energy = 0;
37  produces<float>("DBremEventWeight");
38  produces<float>("DBremLocationX");
39  produces<float>("DBremLocationY");
40  produces<float>("DBremLocationZ");
41  //produces<std::string>("DBremMaterial");
42  produces<float>("DBremAngle");
43  produces<float>("DBremInitialEnergy");
44  produces<float>("DBremFinalEnergy");
45  produces<float>("BiasFactor");
46  pdgs_ = ps.getUntrackedParameter<std::vector<int>>("PDGCodes");
47  edm::LogInfo("DBremWatcher") << "DBremWatcher:: Save Sim Track if PDG code "
48  << "is one from the list of " << pdgs_.size() << " items";
49  for (unsigned int k = 0; k < pdgs_.size(); ++k)
50  edm::LogInfo("DBremWatcher") << "[" << k << "] " << pdgs_[k];
51 }
double biasFactor
Definition: DBremWatcher.h:37
G4ThreeVector aPrimeTraj
Definition: DBremWatcher.h:39
T getUntrackedParameter(std::string const &, T const &) const
G4ThreeVector VertexPos
Definition: DBremWatcher.h:41
float m_weight
Definition: DBremWatcher.h:36
G4ThreeVector finaltraj
Definition: DBremWatcher.h:40
Log< level::Info, false > LogInfo
float f_energy
Definition: DBremWatcher.h:42
std::vector< int > pdgs_
Definition: DBremWatcher.h:34

◆ ~DBremWatcher()

DBremWatcher::~DBremWatcher ( )
override

Definition at line 53 of file DBremWatcher.cc.

53 {}

Member Function Documentation

◆ produce()

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

Implements SimProducer.

Definition at line 116 of file DBremWatcher.cc.

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

116  {
117  if (foundbrem) {
118  std::unique_ptr<float> weight = std::make_unique<float>(m_weight);
119  fEvent.put(std::move(weight), "DBremEventWeight");
120  std::unique_ptr<float> vtxposx = std::make_unique<float>(VertexPos.x());
121  std::unique_ptr<float> vtxposy = std::make_unique<float>(VertexPos.y());
122  std::unique_ptr<float> vtxposz = std::make_unique<float>(VertexPos.z());
123  fEvent.put(std::move(vtxposx), "DBremLocationX");
124  fEvent.put(std::move(vtxposy), "DBremLocationY");
125  fEvent.put(std::move(vtxposz), "DBremLocationZ");
126  std::unique_ptr<float> finalE = std::make_unique<float>(f_energy / GeV);
127  fEvent.put(std::move(finalE), "DBremFinalEnergy");
128  float deflectionAngle = -1;
129  float initialEnergy = sqrt(pow(finaltraj.x() + aPrimeTraj.x(), 2) + pow(finaltraj.y() + aPrimeTraj.y(), 2) +
130  pow(finaltraj.z() + aPrimeTraj.z(), 2) + pow(0.1056 * GeV, 2));
131  G4ThreeVector mother(
132  finaltraj.x() + aPrimeTraj.x(), finaltraj.y() + aPrimeTraj.y(), finaltraj.z() + aPrimeTraj.z());
133  deflectionAngle = mother.angle(finaltraj);
134  std::unique_ptr<float> dAngle = std::make_unique<float>(deflectionAngle);
135  std::unique_ptr<float> initialE = std::make_unique<float>(initialEnergy / GeV);
136  fEvent.put(std::move(dAngle), "DBremAngle");
137  fEvent.put(std::move(initialE), "DBremInitialEnergy");
138  std::unique_ptr<float> bias = std::make_unique<float>(biasFactor);
139  fEvent.put(std::move(bias), "BiasFactor");
140  } else {
141  std::unique_ptr<float> weight = std::make_unique<float>(0.);
142  fEvent.put(std::move(weight), "DBremEventWeight");
143  }
144 }
double biasFactor
Definition: DBremWatcher.h:37
G4ThreeVector aPrimeTraj
Definition: DBremWatcher.h:39
Definition: weight.py:1
G4ThreeVector VertexPos
Definition: DBremWatcher.h:41
float m_weight
Definition: DBremWatcher.h:36
T sqrt(T t)
Definition: SSEVec.h:19
G4ThreeVector finaltraj
Definition: DBremWatcher.h:40
float f_energy
Definition: DBremWatcher.h:42
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 55 of file DBremWatcher.cc.

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

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

55  {
56  G4Track* theTrack = (G4Track*)((*trk)());
57  TrackInformation* trkInfo = (TrackInformation*)(theTrack->GetUserInformation());
58  if (trkInfo) {
59  int pdg = theTrack->GetDefinition()->GetPDGEncoding();
60  G4ThreeVector Vpos = theTrack->GetVertexPosition();
61  const G4VProcess* TrPro = theTrack->GetCreatorProcess();
62  if (TrPro != nullptr) {
63  if ((theTrack->GetCreatorProcess()->GetProcessName()) == "muDBrem") {
64  if (std::find(pdgs_.begin(), pdgs_.end(), pdg) == pdgs_.end()) {
65  //Found the deflected muon track
66  trkInfo->storeTrack(true);
67  if (!theTrack->IsGoodForTracking()) {
68  theTrack->SetGoodForTrackingFlag(true);
69  }
70  f_energy = theTrack->GetTotalEnergy();
71  foundbrem = true;
72  finaltraj = theTrack->GetMomentum();
73  } else {
74  m_weight = theTrack->GetWeight();
75  }
76  }
77  }
78  if (std::find(pdgs_.begin(), pdgs_.end(), pdg) != pdgs_.end()) {
79  //Found an A'
80  trkInfo->storeTrack(true);
81  VertexPos = Vpos;
82  aPrimeTraj = theTrack->GetMomentum();
83  LogDebug("DBremWatcher") << "Save SimTrack the Track " << theTrack->GetTrackID() << " Type "
84  << theTrack->GetDefinition()->GetParticleName() << " Momentum "
85  << theTrack->GetMomentum() / MeV << " MeV/c";
86  }
87  }
88 }
bool storeTrack() const
G4ThreeVector aPrimeTraj
Definition: DBremWatcher.h:39
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.h:41
float m_weight
Definition: DBremWatcher.h:36
G4ThreeVector finaltraj
Definition: DBremWatcher.h:40
float f_energy
Definition: DBremWatcher.h:42
std::vector< int > pdgs_
Definition: DBremWatcher.h:34
#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 92 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().

92  {
93  G4String pname = "muDBrem";
94  G4ProcessTable* ptable = G4ProcessTable::GetProcessTable();
95  G4bool state = true;
96  ptable->SetProcessActivation(pname, state);
97  foundbrem = false;
98 }

◆ 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 100 of file DBremWatcher.cc.

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

100 {}

◆ 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 90 of file DBremWatcher.cc.

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

90 {}

◆ 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 102 of file DBremWatcher.cc.

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

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

102  {
103  G4Track* theTrack = (G4Track*)((*trk)());
104  TrackInformation* trkInfo = (TrackInformation*)(theTrack->GetUserInformation());
105  const G4VProcess* TrPro = theTrack->GetCreatorProcess();
106  if (trkInfo && TrPro != nullptr) {
107  int pdg = theTrack->GetDefinition()->GetPDGEncoding();
108 
109  if (std::find(pdgs_.begin(), pdgs_.end(), pdg) == pdgs_.end() &&
110  (theTrack->GetCreatorProcess()->GetProcessName()) == "muDBrem") {
111  trkInfo->storeTrack(true);
112  }
113  }
114 }
bool storeTrack() const
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.h:34

Member Data Documentation

◆ aPrimeTraj

G4ThreeVector DBremWatcher::aPrimeTraj
private

Definition at line 39 of file DBremWatcher.h.

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

◆ biasFactor

double DBremWatcher::biasFactor
private

Definition at line 37 of file DBremWatcher.h.

Referenced by DBremWatcher(), and produce().

◆ f_energy

float DBremWatcher::f_energy
private

Definition at line 42 of file DBremWatcher.h.

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

◆ finaltraj

G4ThreeVector DBremWatcher::finaltraj
private

Definition at line 40 of file DBremWatcher.h.

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

◆ foundbrem

bool DBremWatcher::foundbrem
private

Definition at line 38 of file DBremWatcher.h.

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

◆ m_weight

float DBremWatcher::m_weight
private

Definition at line 36 of file DBremWatcher.h.

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

◆ MotherId

int DBremWatcher::MotherId
private

Definition at line 35 of file DBremWatcher.h.

◆ pdgs_

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

Definition at line 34 of file DBremWatcher.h.

Referenced by DBremWatcher(), and update().

◆ VertexPos

G4ThreeVector DBremWatcher::VertexPos
private

Definition at line 41 of file DBremWatcher.h.

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