GENIEGenerator
Loading...
Searching...
No Matches
HNLFluxCreator.h
Go to the documentation of this file.
1//----------------------------------------------------------------------------
2/*!
3
4 This is a module for GENIE to read in hadron flux ntuples and construct HNL fluxes
5 on the fly.
6
7 Core loop: + Open flux entry
8 + Get ancestry information
9 + Assume decay to HNL (discard SM nu info, is unneeded)
10 + Calculate HNL production mode based on parameter space read from config
11 + Calculate kinematics of HNL
12 + Return HNL as SimpleHNL object.
13
14\class genie::hnl::FluxCreator
15
16\brief Calculates HNL production kinematics & production vertex.
17 Is a concrete implementation of the FluxRecordVisitorI interface
18
19\author John Plows <komninos-john.plows@physics.ox.ac.uk>
20
21\created April 25th, 2022
22
23\cpright Copyright (c) 2003-2025, The GENIE Collaboration
24 For the full text of the license visit http://copyright.genie-mc.org
25
26 */
27//----------------------------------------------------------------------------
28
29#ifndef _HNL_FLUXCREATOR_H_
30#define _HNL_FLUXCREATOR_H_
31
32// -- C++ includes
33#include <array>
34#include <cassert>
35#include <iomanip> // for momentum balance stream
36#include <map>
37#include <sstream>
38
39// -- ROOT includes
40#include "TChain.h"
41#include "TDecayChannel.h"
42#include "TFile.h"
43#include "TGenPhaseSpace.h"
44#include "TH1.h"
45#include "TH3.h"
46#include "TLorentzVector.h"
47#include "TMath.h"
48#include "TObjArray.h"
49#include "TROOT.h"
50#include "TString.h"
51#include "TSystemDirectory.h"
52#include "TSystem.h"
53#include "TTree.h"
54#include "TVector3.h"
55
56#include <TGeoVolume.h>
57#include <TGeoManager.h>
58#include <TGeoShape.h>
59#include <TGeoBBox.h>
60
61// -- GENIE includes
75
77
84
85const double kRDET = 1.0; // calculate fluxes per m^2
86
87namespace genie{
88
89 namespace hnl{
90
91 class SimpleHNL;
92 class FluxContainer;
93
95
96 public:
97
99 FluxCreator(string name);
100 FluxCreator(string name, string config);
101 ~FluxCreator();
102
103 //-- implement the FluxRecordVisitorI interface
104 void ProcessEventRecord(GHepRecord * event_rec) const;
105
106 // overload the Algorithm::Configure() methods to load private data
107 // members from configuration options
108 void Configure(const Registry & config);
109 void Configure(string config);
110
111 // set the input path for a flux
112 void SetInputFluxPath( std::string finpath ) const;
113 // set path to geometry file
114 void SetGeomFile( string geomfile ) const;
115 // get N(flux input entries)
116 int GetNFluxEntries() const;
117 // set first entry for read-in from chain
118 void SetFirstFluxEntry( int iFirst ) const;
119
120 // get flux info
122
123 // return information about frames
124 std::vector< double > GetB2UTranslation() const { return fB2UTranslation; }
125 std::vector< double > GetB2URotation() const { return fB2URotation; }
126 std::vector< double > GetDetOffset() const { return fDetOffset; }
127 std::vector< double > GetDetRotation() const { return fDetRotation; }
128
129 private:
130
131 void LoadConfig(void);
132
133 // if using root geom, let this module know
134 void SetUsingRootGeom( bool IsUsingRootGeom ) const;
135 void ImportBoundingBox( TGeoBBox * box ) const;
136
137 void SetCurrentEntry( int iCurr ) const;
138
139 // workhorse methods
140 FluxContainer MakeTupleFluxEntry( int iEntry, std::string finpath ) const;
141 void FillNonsense( int iEntry, genie::hnl::FluxContainer & gnmf ) const;
142
143 // init
144 void OpenFluxInput( std::string finpath ) const;
145 void InitialiseTree() const;
146 void InitialiseMeta() const;
147
148 // returns HNL 4-momentum from random decay in same frame as p4par
149 TLorentzVector HNLEnergy( genie::hnl::HNLProd_t hnldm, TLorentzVector p4par ) const;
150 // gets random point in BBox and returns separation to it in BEAM FRAME
151 TVector3 PointToRandomPointInBBox( ) const;
152
153 void ReadBRs() const;
154 std::map< genie::hnl::HNLProd_t, double > GetProductionProbs( int parPDG ) const;
155
156 // Obtain detector dimensions + position
157 // RETHERE: BBox isn't good enough! But roll with it for now
158 void MakeBBox() const;
159 TVector3 ApplyUserRotation( TVector3 vec, bool doBackwards = false ) const;
160 TVector3 ApplyUserRotation( TVector3 vec, TVector3 oriVec, std::vector<double> rotVec, bool doBackwards = false ) const;
161
162 // calculate detector acceptance (== solid angle of proj of det onto unit-radius sphere / (4pi))
163 // NOTE THIS IS A LAB FRAME (==GEOMETRICAL) ACCEPTANCE!!!!
164 // detO == detector BBox centre wrt HNL prod vertex, L{x,y,z} BBox length on each axis. Both [m]
165 double CalculateDetectorAcceptanceSAA( TVector3 detO ) const;
166 // collimation effect calc, returns HNL_acc / geom_acc
167 double CalculateAcceptanceCorrection( TLorentzVector p4par, TLorentzVector p4HNL, double SMECM, double zm, double zp ) const;
168 static double labangle( double * x, double * par ); // function formula for correction
169 // get minimum and maximum deviation from parent momentum to hit detector, [deg]
170 double GetAngDeviation( TLorentzVector p4par, TVector3 detO, bool seekingMax ) const;
171 void GetAngDeviation( TLorentzVector p4par, TVector3 detO, double &zm, double &zp ) const;
172 // returns 1.0 / (area of flux calc)
174
175 // utility function -- is a copy of TGeoChecker::CheckPoint() but doesn't output to cout
176 std::string CheckGeomPoint( Double_t x, Double_t y, Double_t z ) const;
177
178 // current path to keep track of what is loaded
179 mutable std::string fCurrPath = ""; mutable bool fPathLoaded = false;
180 // and which entry we're on
181 mutable int iCurrEntry = 0;
182 // which one was first?
183 mutable int fFirstEntry = 0;
184 // out of how many?
185 mutable int fNEntries = 0;
186
187 // maps to keep P( production )
188 mutable std::map< genie::hnl::HNLProd_t, double > dynamicScores; // map in use
189 mutable std::map< genie::hnl::HNLProd_t, double > dynamicScores_pion;
190 mutable std::map< genie::hnl::HNLProd_t, double > dynamicScores_kaon;
191 mutable std::map< genie::hnl::HNLProd_t, double > dynamicScores_muon;
192 mutable std::map< genie::hnl::HNLProd_t, double > dynamicScores_neuk;
193
195
196 mutable bool isParentOnAxis = true;
197 mutable TGeoVolume * fTopVol = 0;
198 mutable string fGeomFile = "";
199 mutable bool fIsUsingRootGeom = false;
200
201 mutable TChain * ctree = 0, * cmeta = 0;
202
203 mutable double fMass; // HNL mass, GeV
204 mutable std::vector< double > fU4l2s; // couplings
205 mutable bool fIsMajorana;
206 //mutable int fType; // for hist fluxes. 0 ==> only particle, 1 ==> only anti, 2 ==> mix of both
207 //mutable double fMinE = 0.0, fMaxE = 100.0, fAngDev = 0.0; // for hist fluxes
208
209 mutable double fLx, fLy, fLz;
210 mutable double fLxR, fLyR, fLzR; // BBox side [m]
211
212 mutable std::vector< double > fB2UTranslation, fB2URotation;
213 mutable std::vector< double > fDetRotation; // rotation of detector wrt tgt hall
214 mutable std::vector< double > fDetOffset; // offset of det centre wrt geom file origin
215 mutable double fCx, fCy, fCz; // BBox centre wrt HNL prod [m]
216 mutable double fAx1, fAz, fAx2; // Euler angles, extrinsic x-z-x. Tgt-hall to beam
217 mutable double fBx1, fBz, fBx2; // Tgt-hall to detector frame
218
219 mutable double fDx, fDy, fDz; //HNL production point [m]
220
221 //std::vector< double > trVec, roVec;
222
223 mutable bool doPol = true, fixPol = false;
224 mutable double fLPx, fLPy, fLPz; // direction of co-produced lepton == polarisation vector
225 mutable double fLPE;
226 mutable std::vector< double > fFixedPolarisation;
227
228 mutable int fLepPdg; // pdg code of co-produced lepton
229 mutable int fNuPdg; // pdg code of SM neutrino from same decay type
230 mutable double parentMass, parentMomentum, parentEnergy; // GeV
231 mutable double fECM, fSMECM; // GeV
232 mutable double fZm, fZp; // deg
233
234 mutable int fProdChan, fNuProdChan;
235
237
238 static const int maxArray = 30, maxC = 100;
239
240 // tree variables. Add as per necessary.
241 mutable double potnum; ///< N POT for this SM-v
242 mutable int decay_ptype; ///< PDG code of parent
243 mutable double decay_vx, decay_vy, decay_vz; ///< coordinates of prod vtx [cm]
244 mutable double decay_pdpx, decay_pdpy, decay_pdpz; ///< final parent momentum [GeV]
245 mutable double decay_necm; ///< SM v CM energy [GeV]
246 mutable double decay_nimpwt; ///< Importance weight from beamsim
247
248 mutable int arSize, anArSize, trArSize;
249 mutable int djob;
250 mutable double ppvx, ppvy, ppvz;
253 mutable int decay_ppmedium;
255
257
258 mutable int ancestor_pdg[maxArray];
266
270
273
274 // meta variables. Add as necessary
275 mutable int job; ///< beamsim MC job number
276 mutable double pots; ///< how many pot in this job?
277
278 mutable int mArSize;
280 mutable char tgtcfg[maxC], horncfg[maxC], dkvolcfg[maxC];
281 mutable double beam0x, beam0y, beam0z;
282 mutable double beamhwidth, beamvwidth;
283 mutable double beamdxdz, beamdydz;
286
288
289 mutable double POTScaleWeight;
290 mutable std::vector<double> fScales;
291
292 mutable bool fDoingOldFluxCalc = false;
293 mutable bool fRerollPoints = false;
294 mutable double fRadius; // m
295 mutable bool fSupplyingBEAM = false;
296 mutable bool fIsConfigLoaded = false;
297
298 mutable bool fUsingDk2nu = true;
299 mutable string fFinPath, fProdHist;
300 mutable TH1D * fSpectrum = 0, * fIntegrals = 0;
301
302 }; // class FluxCreator
303
304 } // namespace hnl
305} // namespace genie
306
307#endif // #ifndef _HNL_FLUXCREATOR_H_
const double kRDET
Most commonly used PDG codes. A set of utility functions to handle PDG codes is provided in PDGUtils.
GENIE's GHEP MC event record.
Definition GHepRecord.h:45
A registry. Provides the container for algorithm configuration parameters.
Definition Registry.h:65
A GENIE flux container specific for HNL containers. Based on the dk2nu flux paradigm and genie::flux:...
double traj_trky[maxArray]
int decay_ptype
PDG code of parent.
void SetCurrentEntry(int iCurr) const
double GetAngDeviation(TLorentzVector p4par, TVector3 detO, bool seekingMax) const
double traj_trkpz[maxArray]
double decay_vz
coordinates of prod vtx [cm]
int ancestor_nucleus[maxArray]
void SetFirstFluxEntry(int iFirst) const
double ancestor_stoppz[maxArray]
double decay_pdpz
final parent momentum [GeV]
double ancestor_pprodpy[maxArray]
char ancestor_imat[maxArray *maxC]
double ancestor_poly[maxArray]
std::vector< double > fU4l2s
double ancestor_starty[maxArray]
std::vector< double > fB2UTranslation
static const int maxArray
TVector3 ApplyUserRotation(TVector3 vec, bool doBackwards=false) const
std::vector< double > fFixedPolarisation
double ancestor_polz[maxArray]
double ancestor_startt[maxArray]
char ancestor_proc[maxArray *maxC]
std::map< genie::hnl::HNLProd_t, double > dynamicScores_muon
std::map< genie::hnl::HNLProd_t, double > dynamicScores
std::vector< double > GetDetRotation() const
double nuray_wgt[maxArray]
void SetUsingRootGeom(bool IsUsingRootGeom) const
double ancestor_startx[maxArray]
std::map< genie::hnl::HNLProd_t, double > GetProductionProbs(int parPDG) const
std::vector< double > GetDetOffset() const
genie::hnl::FluxContainer fGnmf
int job
beamsim MC job number
double location_z[maxArray]
std::vector< double > fDetOffset
std::vector< double > fScales
double nuray_px[maxArray]
double traj_trkpx[maxArray]
std::string CheckGeomPoint(Double_t x, Double_t y, Double_t z) const
double CalculateAcceptanceCorrection(TLorentzVector p4par, TLorentzVector p4HNL, double SMECM, double zm, double zp) const
double nuray_py[maxArray]
FluxContainer MakeTupleFluxEntry(int iEntry, std::string finpath) const
void ProcessEventRecord(GHepRecord *event_rec) const
std::map< genie::hnl::HNLProd_t, double > dynamicScores_neuk
void OpenFluxInput(std::string finpath) const
double location_x[maxArray]
double ancestor_pprodpx[maxArray]
FluxContainer RetrieveFluxInfo() const
static double labangle(double *x, double *par)
double ancestor_startz[maxArray]
void Configure(const Registry &config)
double decay_necm
SM v CM energy [GeV].
std::map< genie::hnl::HNLProd_t, double > dynamicScores_pion
double ancestor_stoppx[maxArray]
double ancestor_pprodpz[maxArray]
double nuray_pz[maxArray]
double traj_trkx[maxArray]
std::vector< double > fB2URotation
double pots
how many pot in this job?
TLorentzVector HNLEnergy(genie::hnl::HNLProd_t hnldm, TLorentzVector p4par) const
double decay_nimpwt
Importance weight from beamsim.
void SetGeomFile(string geomfile) const
double ancestor_startpx[maxArray]
double ancestor_startpz[maxArray]
TVector3 PointToRandomPointInBBox() const
std::vector< double > fDetRotation
std::map< genie::hnl::HNLProd_t, double > dynamicScores_kaon
double traj_trkpy[maxArray]
void ImportBoundingBox(TGeoBBox *box) const
double traj_trkz[maxArray]
double ancestor_stoppy[maxArray]
char location_name[maxArray *maxC]
double ancestor_polx[maxArray]
double location_y[maxArray]
double potnum
N POT for this SM-v.
double ancestor_startpy[maxArray]
std::vector< double > GetB2URotation() const
double nuray_E[maxArray]
std::vector< double > GetB2UTranslation() const
double CalculateDetectorAcceptanceSAA(TVector3 detO) const
void SetInputFluxPath(std::string finpath) const
char ancestor_ivol[maxArray *maxC]
void FillNonsense(int iEntry, genie::hnl::FluxContainer &gnmf) const
Utilities for simulating the decay of Heavy Neutral Leptons.
enum genie::hnl::t_HNLProd HNLProd_t
THE MAIN GENIE PROJECT NAMESPACE
Definition AlgCmp.h:25