All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Private Member Functions | Private Attributes | List of all members
evgen::ldm::ReThrowRayTraceBox Class Reference

ReThrowRayTraceBox class definiton. More...

Inheritance diagram for evgen::ldm::ReThrowRayTraceBox:
evgen::ldm::IRayTrace evgen::ldm::IMeVPrtlStage

Public Member Functions

 ReThrowRayTraceBox (fhicl::ParameterSet const &pset)
 Constructor. More...
 
 ~ReThrowRayTraceBox ()
 Destructor. More...
 
void configure (const fhicl::ParameterSet &) override
 Interface for configuring the particular algorithm tool. More...
 
bool IntersectDetector (MeVPrtlFlux &flux, std::array< TVector3, 2 > &intersection, double &weight) override
 
TLorentzVector ThrowMeVPrtlMomentum (const MeVPrtlFlux &flux)
 
double MaxWeight () override
 
- Public Member Functions inherited from evgen::ldm::IRayTrace
virtual ~IRayTrace () noexcept=default
 Virtual Destructor. More...
 
- Public Member Functions inherited from evgen::ldm::IMeVPrtlStage
virtual ~IMeVPrtlStage () noexcept
 Virtual Destructor. More...
 
 IMeVPrtlStage (const char *name)
 
TVector3 RandomUnitVector ()
 
double GetRandom ()
 
const char * Name ()
 

Private Member Functions

void CalculateMaxWeight ()
 

Private Attributes

geo::BoxBoundedGeo fBox
 
unsigned fNThrows
 
double fReferenceLabSolidAngle
 
double fReferencePrtlMass
 
int fReferenceScndPDG
 
double fReferenceKaonEnergy
 
double fMaxWeight
 

Additional Inherited Members

- Protected Attributes inherited from evgen::ldm::IMeVPrtlStage
CLHEP::HepRandomEngine * fEngine
 
const char * fName
 

Detailed Description

ReThrowRayTraceBox class definiton.

Definition at line 39 of file ReThrowRayTraceBox_tool.cc.

Constructor & Destructor Documentation

evgen::ldm::ReThrowRayTraceBox::ReThrowRayTraceBox ( fhicl::ParameterSet const &  pset)

Constructor.

Definition at line 77 of file ReThrowRayTraceBox_tool.cc.

77  :
78  IMeVPrtlStage("ReThrowRayTraceBox")
79 {
80  this->configure(pset);
81 }
void configure(const fhicl::ParameterSet &) override
Interface for configuring the particular algorithm tool.
IMeVPrtlStage(const char *name)
Definition: IMeVPrtlStage.h:46
evgen::ldm::ReThrowRayTraceBox::~ReThrowRayTraceBox ( )

Destructor.

Definition at line 85 of file ReThrowRayTraceBox_tool.cc.

86 {
87 }

Member Function Documentation

void evgen::ldm::ReThrowRayTraceBox::CalculateMaxWeight ( )
private

Definition at line 120 of file ReThrowRayTraceBox_tool.cc.

120  {
121  double secondary_mass = 0.;
122  switch (fReferenceScndPDG) {
123  case 11:
124  case -11:
125  secondary_mass = Constants::Instance().elec_mass;
126  break;
127  case 13:
128  case -13:
129  secondary_mass = Constants::Instance().muon_mass;
130  break;
131  case 211:
132  case -211:
133  secondary_mass = Constants::Instance().piplus_mass;
134  break;
135  default:
136  std::cerr << "RETHROWRAYTACE -- bad secondary pdg: " << fReferenceScndPDG << std::endl;
137  std::cout << "RETHROWRAYTACE -- bad secondary pdg: " << fReferenceScndPDG << std::endl;
138  break;
139  }
140 
141  double p = twobody_momentum(Constants::Instance().kplus_mass, secondary_mass, fReferencePrtlMass);
142  double E = sqrt(p*p + fReferencePrtlMass*fReferencePrtlMass);
143 
145  double kbeta = sqrt(1 - 1. / (kgamma * kgamma));
146 
147  double scale = kgamma * kgamma * (p + kbeta*E) * (p + kbeta*E) / (p*p);
148 
150 }
BEGIN_PROLOG could also be cerr
pdgs p
Definition: selectors.fcl:22
double twobody_momentum(double parent_mass, double childA_mass, double childB_mass)
Definition: Constants.cpp:73
process_name E
static const Constants & Instance()
Definition: Constants.h:68
BEGIN_PROLOG could also be cout
void evgen::ldm::ReThrowRayTraceBox::configure ( const fhicl::ParameterSet &  )
overridevirtual

Interface for configuring the particular algorithm tool.

Parameters
ParameterSetThe input set of parameters for configuration

Implements evgen::ldm::IMeVPrtlStage.

Definition at line 90 of file ReThrowRayTraceBox_tool.cc.

91 {
92  if (pset.has_key("Box")) {
93  std::array<double, 6> box_config = pset.get<std::array<double, 6>>("Box");
94  // xmin, xmax, ymin, ymax, zmin, zmax
95  fBox = geo::BoxBoundedGeo(box_config[0], box_config[1], box_config[2], box_config[3], box_config[4], box_config[5]);
96  }
97  else {
98  const geo::GeometryCore *geometry = lar::providerFrom<geo::Geometry>();
99  fBox = geometry->DetectorEnclosureBox(pset.get<std::string>("Volume"));
100  }
101 
102  fNThrows = pset.get<unsigned>("NThrows", 10000);
103 
104  std::cout << "Detector Box." << std::endl;
105  std::cout << "X " << fBox.MinX() << " " << fBox.MaxX() << std::endl;
106  std::cout << "Y " << fBox.MinY() << " " << fBox.MaxY() << std::endl;
107  std::cout << "Z " << fBox.MinZ() << " " << fBox.MaxZ() << std::endl;
108 
109  fReferenceLabSolidAngle = pset.get<double>("ReferenceLabSolidAngle");
110  fReferencePrtlMass = pset.get<double>("ReferencePrtlMass");
111  fReferenceScndPDG = pset.get<int>("ReferenceScndPDG");
112  fReferenceKaonEnergy = pset.get<double>("ReferenceKaonEnergy");
113 
115 
116  // TEST: compare results of this weight to calcENuWeight
117 
118 }
const geo::GeometryCore * geometry
geo::BoxBoundedGeo DetectorEnclosureBox(std::string const &name="volDetEnclosure") const
Description of geometry of one entire detector.
A base class aware of world box coordinatesAn object describing a simple shape can inherit from this ...
Definition: BoxBoundedGeo.h:33
BEGIN_PROLOG could also be cout
bool evgen::ldm::ReThrowRayTraceBox::IntersectDetector ( MeVPrtlFlux flux,
std::array< TVector3, 2 > &  intersection,
double &  weight 
)
overridevirtual

Implements evgen::ldm::IRayTrace.

Definition at line 171 of file ReThrowRayTraceBox_tool.cc.

171  {
172  // try out the mevprtl direction a bunch of times
173  std::vector<std::vector<TVector3>> allIntersections;
174  std::vector<TLorentzVector> allHMom;
175 
176  for (unsigned i = 0; i <= fNThrows; i++) {
177  TLorentzVector mevprtl_mom = ThrowMeVPrtlMomentum(flux);
178  std::vector<TVector3> box_intersections = fBox.GetIntersections(flux.pos.Vect(), mevprtl_mom.Vect().Unit());
179 
180  // Does this ray intersect the box?
181  if (box_intersections.size() != 2) {
182  continue;
183  }
184 
185  TVector3 A = box_intersections[0];
186 
187  // if the ray points the wrong way, it doesn't intersect
188  if (mevprtl_mom.Vect().Unit().Dot((A - flux.pos.Vect()).Unit()) < 0.) {
189  std::cout << "RAYTRACE: MeVPrtl points wrong way" << std::endl;
190  std::cout << "Pos: " << flux.pos.X() << " " << flux.pos.Y() << " " << flux.pos.Z() << std::endl;
191  std::cout << "A: " << A.X() << " " << A.Y() << " " << A.Z() << std::endl;
192  std::cout << "P: " << mevprtl_mom.Vect().Unit().X() << " " << mevprtl_mom.Vect().Unit().Y() << " " << mevprtl_mom.Vect().Unit().Z() << std::endl;
193  continue;
194  }
195 
196  // if we're here, we have a valid ray
197  allIntersections.push_back(box_intersections);
198  allHMom.push_back(mevprtl_mom);
199  }
200 
201  std::cout << "Prtl intersected (" << allIntersections.size() << " / " << fNThrows << ") times.\n";
202 
203  // did we get a hit?
204  if (allIntersections.size() == 0) {
205  return false;
206  }
207 
208  unsigned ind = CLHEP::RandFlat::shootInt(fEngine, 0, allIntersections.size()-1); // inclusive?
209 
210  TLorentzVector mevprtl_mom = allHMom[ind];
211  std::vector<TVector3> box_intersections = allIntersections[ind];
212 
213  TVector3 A = box_intersections[0];
214  TVector3 B = box_intersections[1];
215 
216  // make sure that the flux start lies outside the detector
217  if ((flux.pos.Vect() - A).Mag() < (A-B).Mag() && (flux.pos.Vect() - B).Mag() < (A-B).Mag()) {
218  throw cet::exception("ReThrowRayTraceBox Exception", "Input mevprtl flux starts inside detector volume: "
219  "MeVPrtl start At: (" + std::to_string(flux.pos.X()) + ", " + std::to_string(flux.pos.Y()) + ", " + std::to_string(flux.pos.Z()) + "). "
220  "Intersection A At: " + std::to_string(A.X()) + ", " + std::to_string(A.Y()) + ", " + std::to_string(A.Z()) + "). "
221  "Intersection B At: " + std::to_string(B.X()) + ", " + std::to_string(B.Y()) + ", " + std::to_string(B.Z()) + ").\n"
222  );
223  }
224 
225  // set things
226  weight = (double)allIntersections.size() / fNThrows;
227  flux.mom = mevprtl_mom;
228  // transform to beam-coord frame
229  flux.mom_beamcoord = mevprtl_mom;
230  flux.mom_beamcoord.Boost(-flux.kmom.BoostVector()); // Boost to kaon rest frame
231  flux.mom_beamcoord.Boost(flux.kmom_beamcoord.BoostVector()); // And to beam coordinate frame
232  // Two boosts make a rotation!
233 
234  // Set the secondary 4-momentum
235  flux.sec = flux.kmom - flux.mom;
236  flux.sec_beamcoord = flux.kmom_beamcoord - flux.mom_beamcoord;
237 
238  if ((flux.pos.Vect() - A).Mag() < (flux.pos.Vect() - B).Mag()) {
239  intersection = {A, B}; // A is entry, B is exit
240  }
241  else {
242  intersection = {B, A}; // reversed
243  }
244 
245  std::cout << "Kaon 4P: " << flux.kmom.E() << " " << flux.kmom.Px() << " " << flux.kmom.Py() << " " << flux.kmom.Pz() << std::endl;
246  std::cout << "Selected Prtl 4P: " << flux.mom.E() << " " << flux.mom.Px() << " " << flux.mom.Py() << " " << flux.mom.Pz() << std::endl;
247  std::cout << "Selected Scdy 4P: " << flux.sec.E() << " " << flux.sec.Px() << " " << flux.sec.Py() << " " << flux.sec.Pz() << std::endl;
248 
249  return true;
250 }
CLHEP::HepRandomEngine * fEngine
Definition: IMeVPrtlStage.h:79
TLorentzVector ThrowMeVPrtlMomentum(const MeVPrtlFlux &flux)
std::string to_string(WindowPattern const &pattern)
std::vector< TVector3 > GetIntersections(TVector3 const &TrajectoryStart, TVector3 const &TrajectoryDirect) const
Calculates the entry and exit points of a trajectory on the box surface.
float A
Definition: dedx.py:137
BEGIN_PROLOG could also be cout
double evgen::ldm::ReThrowRayTraceBox::MaxWeight ( )
inlineoverridevirtual

Implements evgen::ldm::IMeVPrtlStage.

Definition at line 59 of file ReThrowRayTraceBox_tool.cc.

59  {
60  return fMaxWeight;
61  }
TLorentzVector evgen::ldm::ReThrowRayTraceBox::ThrowMeVPrtlMomentum ( const MeVPrtlFlux flux)

Definition at line 152 of file ReThrowRayTraceBox_tool.cc.

152  {
153  // pick a direction for the rest-frame
154  TVector3 dir = RandomUnitVector();
155 
156  // make the mevprtl momentum this
157  //
158  // Move the mevprtl momentum back to the kaon rest frame
159  TLorentzVector mevprtl_mom = flux.mom;
160  mevprtl_mom.Boost(-flux.kmom.BoostVector());
161 
162  mevprtl_mom = TLorentzVector(mevprtl_mom.P() * dir, mevprtl_mom.E());
163 
164  // boost back
165  mevprtl_mom.Boost(flux.kmom.BoostVector());
166 
167  return mevprtl_mom;
168 
169 }
tuple dir
Definition: dropbox.py:28

Member Data Documentation

geo::BoxBoundedGeo evgen::ldm::ReThrowRayTraceBox::fBox
private

Definition at line 64 of file ReThrowRayTraceBox_tool.cc.

double evgen::ldm::ReThrowRayTraceBox::fMaxWeight
private

Definition at line 74 of file ReThrowRayTraceBox_tool.cc.

unsigned evgen::ldm::ReThrowRayTraceBox::fNThrows
private

Definition at line 65 of file ReThrowRayTraceBox_tool.cc.

double evgen::ldm::ReThrowRayTraceBox::fReferenceKaonEnergy
private

Definition at line 70 of file ReThrowRayTraceBox_tool.cc.

double evgen::ldm::ReThrowRayTraceBox::fReferenceLabSolidAngle
private

Definition at line 67 of file ReThrowRayTraceBox_tool.cc.

double evgen::ldm::ReThrowRayTraceBox::fReferencePrtlMass
private

Definition at line 68 of file ReThrowRayTraceBox_tool.cc.

int evgen::ldm::ReThrowRayTraceBox::fReferenceScndPDG
private

Definition at line 69 of file ReThrowRayTraceBox_tool.cc.


The documentation for this class was generated from the following file: