SKKinematicsGenerator.cxx
Go to the documentation of this file.
1 //____________________________________________________________________________
2 /*
3  Copyright (c) 2003-2019, The GENIE Collaboration
4  For the full text of the license visit http://copyright.genie-mc.org
5  or see $GENIE/LICENSE
6 
7  Authors: Chris Marshall <marshall \at pas.rochester.edu>
8  University of Rochester
9 
10  Martti Nirkko
11  University of Berne
12 
13  For the class documentation see the corresponding header file.
14 
15 */
16 //____________________________________________________________________________
17 
18 #include <cstdlib>
19 
20 #include <TMath.h>
21 
39 
40 using namespace genie;
41 using namespace genie::constants;
42 using namespace genie::controls;
43 using namespace genie::utils;
44 
45 //___________________________________________________________________________
47 KineGeneratorWithCache("genie::SKKinematicsGenerator")
48 {
49  //fEnvelope = 0;
50 }
51 //___________________________________________________________________________
53 KineGeneratorWithCache("genie::SKKinematicsGenerator", config)
54 {
55  //fEnvelope = 0;
56 }
57 //___________________________________________________________________________
59 {
60  //if(fEnvelope) delete fEnvelope;
61 }
62 //___________________________________________________________________________
64 {
65  if(fGenerateUniformly) {
66  LOG("SKKinematics", pNOTICE)
67  << "Generating kinematics uniformly over the allowed phase space";
68  }
69 
70  //-- Access cross section algorithm for running thread
72  const EventGeneratorI * evg = rtinfo->RunningThread();
73  fXSecModel = evg->CrossSectionAlg();
75 }
76 //___________________________________________________________________________
78 {
79  // Get the Primary Interacton object
80  Interaction * interaction = evrec->Summary();
81  interaction->SetBit(kISkipProcessChk);
82  interaction->SetBit(kISkipKinematicChk);
83 
84  // Initialise a random number generator
86 
87  //-- For the subsequent kinematic selection with the rejection method:
88  // Calculate the max differential cross section or retrieve it from the
89  // cache. Throw an exception and quit the evg thread if a non-positive
90  // value is found.
91  // If the kinematics are generated uniformly over the allowed phase
92  // space the max xsec is irrelevant
93  double xsec_max = (fGenerateUniformly) ? -1 : this->MaxXSec(evrec);
94 
95  // Determine lepton and kaon masses
96  int leppdg = interaction->FSPrimLeptonPdg();
97  const TLorentzVector pnuc4 = interaction->InitState().Tgt().HitNucP4(); // 4-momentum of struck nucleon in lab frame
98  TVector3 beta = pnuc4.BoostVector();
99  TLorentzVector P4_nu = *(interaction->InitStatePtr()->GetProbeP4(kRfHitNucRest)); // struck nucleon rest frame
100 
101  double enu = P4_nu.E(); // in nucleon rest frame
102  int kaon_pdgc = interaction->ExclTag().StrangeHadronPdg();
103  double mk = PDGLibrary::Instance()->Find(kaon_pdgc)->Mass();
104  double ml = PDGLibrary::Instance()->Find(leppdg)->Mass();
105 
106  // Maximum possible kinetic energy
107  const double Tkmax = enu - mk - ml;
108  const double Tlmax = enu - mk - ml;
109 
110  // Tkmax <= 0 means we are below threshold for sure
111  if( Tkmax <= 0.0 ) {
112  LOG("SKKinematics", pWARN) << "No available phase space";
113  evrec->EventFlags()->SetBitNumber(kKineGenErr, true);
115  exception.SetReason("No available phase space");
116  exception.SwitchOnFastForward();
117  throw exception;
118  }
119 
120  const double Tkmin = 0.0;
121  const double Tlmin = 0.0;
122  // for performance, use log( 1 - cos(theta_lepton) ) instead of cos(theta_lepton) because it is sharply peaked near 1.0
123  const double xmin = fMinLog1MinusCosTheta; // set in LoadConfig
124  const double xmax = 0.69314718056; // log(2) is physical boundary
125  const double phikqmin = 0.0;
126  const double phikqmax = 2.0 * kPi;
127  const double dtk = Tkmax - Tkmin;
128  const double dtl = Tlmax - Tlmin;
129  const double dx = xmax - xmin;
130  const double dphikq = phikqmax - phikqmin;
131 
132  //------ Try to select a valid tk, tl, costhetal, phikq quadruplet
133 
134  unsigned int iter = 0;
135  bool accept = false;
136  double xsec = -1; // diff XS
137  double tk = -1; // kaon kinetic energy
138  double tl = -1; // lepton kinetic energy
139  double costhetal = -1; // cosine of lepton angle
140  double phikq = -1; // azimuthal angle between kaon and q vector
141 
142  while(1) {
143  iter++;
144  if(iter > kRjMaxIterations) {
145  LOG("SKKinematics", pWARN)
146  << "*** Could not select a valid (tk, tl, costhetal) triplet after "
147  << iter << " iterations";
148  evrec->EventFlags()->SetBitNumber(kKineGenErr, true);
150  exception.SetReason("Couldn't select kinematics");
151  exception.SwitchOnFastForward();
152  throw exception;
153  }
154 
155  if(fGenerateUniformly) {
156  tk = Tkmin + dtk * rnd->RndKine().Rndm();
157  tl = Tlmin + dtl * rnd->RndKine().Rndm();
158  double x = xmin + dx * rnd->RndKine().Rndm(); // log(1-costheta)
159  costhetal = 1.0 - TMath::Exp(x);
160  phikq = phikqmin + dphikq * rnd->RndKine().Rndm();
161  } else {
162  tk = Tkmin + dtk * rnd->RndKine().Rndm();
163  tl = Tlmin + dtl * rnd->RndKine().Rndm();
164  double x = xmin + dx * rnd->RndKine().Rndm(); // log(1-costheta)
165  costhetal = 1.0 - TMath::Exp(x);
166  phikq = phikqmin + dphikq * rnd->RndKine().Rndm();
167  }
168 
169  LOG("SKKinematics", pDEBUG) << "Trying: Tk = " << tk << ", Tl = " << tl << ", cosThetal = " << costhetal << ", phikq = " << phikq;
170 
171  // nucleon rest frame! these need to be boosted to the lab frame before they become actual particles
172  interaction->KinePtr()->SetKV(kKVTk, tk);
173  interaction->KinePtr()->SetKV(kKVTl, tl);
174  interaction->KinePtr()->SetKV(kKVctl, costhetal);
175  interaction->KinePtr()->SetKV(kKVphikq, phikq);
176 
177  // lorentz invariant stuff, but do all the calculations in the nucleon rest frame
178  double el = tl + ml;
179  double pl = TMath::Sqrt(el*el - ml*ml);
180  double M = interaction->InitState().Tgt().Mass();
181  TVector3 lepton_3vector = TVector3(0,0,0);
182  lepton_3vector.SetMagThetaPhi(pl,TMath::ACos(costhetal),0.0);
183  TLorentzVector P4_lep( lepton_3vector, tl+ml );
184  TLorentzVector q = P4_nu - P4_lep;
185  double Q2 = -q.Mag2();
186  double xbj = Q2/(2*M*q.E());
187  double y = q.E()/P4_nu.E();
188  double W2 = (pnuc4+q).Mag2();
189 
190 
191  // computing cross section for the current kinematics
192  xsec = fXSecModel->XSec(interaction, kPSTkTlctl);
193 
194  //-- decide whether to accept the current kinematics
195  if(!fGenerateUniformly) {
196  // Jacobian is 1-costheta for x = log(1-costheta)
197  double max = xsec_max;
198  double t = max * rnd->RndKine().Rndm();
199  double J = TMath::Abs(1. - costhetal);
200 
201  this->AssertXSecLimits(interaction, xsec*J, max);
202 
203  if( xsec*J > xsec_max ) { // freak out if this happens
204  LOG("SKKinematics", pWARN)
205  << "!!!!!!XSEC ABOVE MAX!!!!! xsec= " << xsec << ", J= " << J << ", xsec*J = " << xsec*J << " max= " << xsec_max;
206  }
207 
208  accept = (t< J*xsec);
209  }
210  else {
211  accept = (xsec>0);
212  }
213 
214  //-- If the generated kinematics are accepted, finish-up module's job
215  if(accept) {
216 
217  // calculate the stuff
218 
219  // for uniform kinematics, compute an event weight as
220  // wght = (phase space volume)*(differential xsec)/(event total xsec)
221  if(fGenerateUniformly) {
222  double wght = 1.0; // change this
223  wght *= evrec->Weight();
224  LOG("SKKinematics", pNOTICE) << "Current event wght = " << wght;
225  evrec->SetWeight(wght);
226  }
227  LOG("SKKinematics", pWARN) << "\nLepton energy (rest frame) = " << el << " kaon = " << tl + mk;
228 
229  // reset bits
230  interaction->ResetBit(kISkipProcessChk);
231  interaction->ResetBit(kISkipKinematicChk);
232 
233  interaction->KinePtr()->SetKV(kKVSelTk, tk); // nucleon rest frame
234  interaction->KinePtr()->SetKV(kKVSelTl, tl); // nucleon rest frame
235  interaction->KinePtr()->SetKV(kKVSelctl, costhetal); // nucleon rest frame
236  interaction->KinePtr()->SetKV(kKVSelphikq, phikq); // nucleon rest frame
237  interaction->KinePtr()->SetQ2(Q2, true);
238  interaction->KinePtr()->SetW(TMath::Sqrt(W2), true);
239  interaction->KinePtr()->Setx( xbj, true );
240  interaction->KinePtr()->Sety( y, true );
241  interaction->KinePtr()->ClearRunningValues();
242 
243  // set the cross section for the selected kinematics
244  evrec->SetDiffXSec(xsec*(1.0-costhetal),kPSTkTlctl); // phase space is really log(1-costheta)
245 
246  return;
247  }
248  }// iterations
249 }
250 //___________________________________________________________________________
252 {
253 // Computes the maximum differential cross section in the requested phase
254 // space. This method overloads KineGeneratorWithCache::ComputeMaxXSec
255 // method and the value is cached at a circular cache branch for retrieval
256 // during subsequent event generation.
257 
258 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
259  SLOG("SKKinematics", pDEBUG)
260  << "Scanning the allowed phase space {K} for the max(dxsec/d{K})";
261 #endif
262 
263  double max_xsec = 0;
264  double max_tk = -1;
265  double max_tl = -1;
266  double max_ctl = -1;
267 
268  const int Ntk = 100;
269  const int Ntl = 100;
270  const int Nctl = 100;
271  // don't do phi_kq -- the maximum will always occur at phi_kq = pi
272 
273  int leppdg = in->FSPrimLeptonPdg();
274  double enu = in->InitState().ProbeE(kRfHitNucRest); // Enu in nucleon rest frame
275  int kaon_pdgc = in->ExclTag().StrangeHadronPdg();
276  double mk = PDGLibrary::Instance()->Find(kaon_pdgc)->Mass();
277  double ml = PDGLibrary::Instance()->Find(leppdg)->Mass();
278 
279  const double Tkmax = enu - mk - ml;
280  const double Tlmax = enu - mk - ml;
281  const double Tkmin = 0.0;
282  const double Tlmin = 0.0;
283  // cross section is sharply peaked at forward lepton
284  // faster to scan over log(1 - cos(theta)) = x
285  const double xmin = fMinLog1MinusCosTheta; // set in LoadConfig
286  const double xmax = 0.69314718056; // physical limit -- this is log(2)
287  const double dtk = (Tkmax - Tkmin)/Ntk;
288  const double dtl = (Tlmax - Tlmin)/Ntl;
289  const double dx = (xmax - xmin)/Nctl;
290 
291  // XS is 0 when the kinetic energies are == 0, so start at 1
292  for(int i = 1; i < Ntk; i++) {
293  double tk = Tkmin + dtk*i;
294  for(int j = 1; j < Ntl; j++) {
295  double tl = Tlmin + dtl*j;
296  for(int k = 0; k < Nctl; k++) {
297  double log_1_minus_cosine_theta_lepton = xmin + dx*k;
298  // XS takes cos(theta_lepton), we are scanning over log(1-that) to save time, convert to cos(theta_lepton) here
299  double ctl = 1.0 - TMath::Exp(log_1_minus_cosine_theta_lepton);
300 
301  // The cross section is 4D, but the maximum differential cross section always occurs at phi_kq = pi (or -pi)
302  // this is because there is more phase space for the kaon when it is opposite the lepton
303  // to save time in this loop, just set phi_kq to pi
304  double phikq = kPi;
305 
306  in->KinePtr()->SetKV(kKVTk, tk);
307  in->KinePtr()->SetKV(kKVTl, tl);
308  in->KinePtr()->SetKV(kKVctl, ctl);
309  in->KinePtr()->SetKV(kKVphikq, phikq);
310 
311  double xsec = fXSecModel->XSec(in, kPSTkTlctl);
312 
313  // xsec returned by model is d4sigma/(dtk dtl dcosthetal dphikq)
314  // convert lepton theta to log(1-costheta) by multiplying by jacobian 1 - costheta
315  xsec *= (1.0 - ctl);
316 
317  if( xsec > max_xsec ) {
318  max_xsec = xsec;
319  max_tk = tk;
320  max_tl = tl;
321  max_ctl = ctl;
322  }
323  }//ctl
324  }//tl
325  }//tk
326 
327  LOG("SKKinematics", pINFO) << "Max XSec is " << max_xsec << " for enu = " << enu << " tk = " << max_tk << " tl = " << max_tl << " cosine theta = " << max_ctl;
328 
329  // Apply safety factor, since value retrieved from the cache might
330  // correspond to a slightly different energy.
331  max_xsec *= fSafetyFactor;
332 
333 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
334  SLOG("SKKinematics", pDEBUG) << in->AsString();
335  SLOG("SKKinematics", pDEBUG) << "Max xsec in phase space = " << max_xsec;
336  SLOG("SKKinematics", pDEBUG) << "Computed using alg = " << fXSecModel->Id();
337 #endif
338 
339 
340 
341  return max_xsec;
342 }
343 
344 //___________________________________________________________________________
345 double SKKinematicsGenerator::Energy(const Interaction * interaction) const
346 {
347 // Override the base class Energy() method to cache the max xsec for the
348 // neutrino energy in the LAB rather than in the hit nucleon rest frame.
349 
350  const InitialState & init_state = interaction->InitState();
351  double E = init_state.ProbeE(kRfLab);
352  return E;
353 }
354 //___________________________________________________________________________
356 {
357  Algorithm::Configure(config);
358  this->LoadConfig();
359 }
360 //____________________________________________________________________________
362 {
363  Algorithm::Configure(config);
364  this->LoadConfig();
365 }
366 //____________________________________________________________________________
368 {
369  // max xsec safety factor (for rejection method) and min cached energy
370  this->GetParamDef("MaxXSec-SafetyFactor", fSafetyFactor, 1.5);
371  this->GetParamDef("Cache-MinEnergy", fEMin, 0.6);
372 
373  // Generate kinematics uniformly over allowed phase space and compute
374  // an event weight?
375  this->GetParamDef("UniformOverPhaseSpace", fGenerateUniformly, false);
376 
377  // Maximum allowed fractional cross section deviation from maxim cross
378  // section used in rejection method
379  this->GetParamDef("MaxXSec-DiffTolerance", fMaxXSecDiffTolerance, 999999.);
381 
382  //
383  this->GetParamDef("MaxXSec-MinLog1MinusCosTheta", fMinLog1MinusCosTheta, -20.0);
384 }
385 //____________________________________________________________________________
const double kPi
virtual double MaxXSec(GHepRecord *evrec) const
Basic constants.
virtual void SetWeight(double wght)
Definition: GHepRecord.h:131
std::map< std::string, double > xmax
bool fGenerateUniformly
uniform over allowed phase space + event weight?
double J(double q0, double q3, double Enu, double ml)
Definition: MECUtils.cxx:141
THE MAIN GENIE PROJECT NAMESPACE
Definition: GeneratorBase.h:8
void Configure(const Registry &config)
double Q2(const Interaction *const i)
Definition: KineUtils.cxx:991
virtual Interaction * Summary(void) const
Definition: GHepRecord.cxx:135
static RandomGen * Instance()
Access instance.
Definition: RandomGen.cxx:79
void SetQ2(double Q2, bool selected=false)
Definition: Kinematics.cxx:265
Kinematics * KinePtr(void) const
Definition: Interaction.h:76
double fSafetyFactor
maxxsec -> maxxsec * safety_factor
virtual double XSec(const Interaction *i, KinePhaseSpace_t k=kPSfE) const =0
Compute the cross section for the input interaction.
virtual void AssertXSecLimits(const Interaction *in, double xsec, double xsec_max) const
double fMaxXSecDiffTolerance
max{100*(xsec-maxxsec)/.5*(xsec+maxxsec)} if xsec>maxxsec
Defines the EventGeneratorI interface.
virtual double Weight(void) const
Definition: GHepRecord.h:125
::xsd::cxx::tree::exception< char > exception
Definition: Database.h:225
Definition: config.py:1
Double_t beta
double Mass(void) const
Definition: Target.cxx:241
Abstract class. Provides a data caching mechanism for for concrete implementations of the EventRecord...
A singleton holding random number generator classes. All random number generation in GENIE should tak...
Definition: RandomGen.h:30
string AsString(void) const
int FSPrimLeptonPdg(void) const
final state primary lepton pdg
Summary information for an interaction.
Definition: Interaction.h:56
const TLorentzVector & HitNucP4(void) const
Definition: Target.h:92
double enu
Definition: runWimpSim.h:113
double ComputeMaxXSec(const Interaction *in) const
An exception thrown by EventRecordVisitorI when the normal processing sequence has to be disrupted (f...
double dx[NP][NC]
int StrangeHadronPdg(void) const
Definition: XclsTag.h:53
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:97
Float_t E
Definition: plot.C:20
virtual void Configure(const Registry &config)
Definition: Algorithm.cxx:70
void ProcessEventRecord(GHepRecord *event_rec) const
const double j
Definition: BetheBloch.cxx:29
#define pINFO
Definition: Messenger.h:63
TRandom3 & RndKine(void) const
rnd number generator used by kinematics generators
Definition: RandomGen.h:51
Misc GENIE control constants.
Double_t xsec[nknots]
Definition: testXsec.C:47
#define pWARN
Definition: Messenger.h:61
const UInt_t kISkipKinematicChk
if set, skip kinematic validity checks
Definition: Interaction.h:48
virtual const XSecAlgorithmI * CrossSectionAlg(void) const =0
void Setx(double x, bool selected=false)
Definition: Kinematics.cxx:241
void SetKV(KineVar_t kv, double value)
Definition: Kinematics.cxx:345
float Mag2() const
static RunningThreadInfo * Instance(void)
virtual const AlgId & Id(void) const
Get algorithm ID.
Definition: Algorithm.h:98
void SetW(double W, bool selected=false)
Definition: Kinematics.cxx:289
ifstream in
Definition: comparison.C:7
static PDGLibrary * Instance(void)
Definition: PDGLibrary.cxx:43
virtual TBits * EventFlags(void) const
Definition: GHepRecord.h:118
Var Sqrt(const Var &v)
Use to take sqrt of a var.
Definition: Var.cxx:326
A registry. Provides the container for algorithm configuration parameters.
Definition: Registry.h:66
void Sety(double y, bool selected=false)
Definition: Kinematics.cxx:253
const XclsTag & ExclTag(void) const
Definition: Interaction.h:72
TParticlePDG * Find(int pdgc)
Definition: PDGLibrary.cxx:61
double fEMin
min E for which maxxsec is cached - forcing explicit calc.
assert(nhit_max >=nhit_nbins)
InitialState * InitStatePtr(void) const
Definition: Interaction.h:74
const InitialState & InitState(void) const
Definition: Interaction.h:69
double Energy(const Interaction *in) const
#define pNOTICE
Definition: Messenger.h:62
void ClearRunningValues(void)
Definition: Kinematics.cxx:357
bool GetParamDef(const RgKey &name, T &p, const T &def) const
const Target & Tgt(void) const
Definition: InitialState.h:67
static const unsigned int kRjMaxIterations
Definition: Controls.h:27
void CalculateKin_AtharSingleKaon(GHepRecord *event_rec) const
T max(sqlite3 *const db, std::string const &table_name, std::string const &column_name)
Definition: statistics.h:68
const EventGeneratorI * RunningThread(void)
#define SLOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a short string (using the FUNCTION and...
Definition: Messenger.h:85
double ProbeE(RefFrame_t rf) const
GENIE&#39;s GHEP MC event record.
Definition: GHepRecord.h:46
Most commonly used PDG codes. A set of utility functions to handle PDG codes is provided in PDGUtils...
Keep info on the event generation thread currently on charge. This is used so that event generation m...
Root of GENIE utility namespaces.
TLorentzVector * GetProbeP4(RefFrame_t rf=kRfHitNucRest) const
const UInt_t kISkipProcessChk
if set, skip process validity checks
Definition: Interaction.h:47
virtual void SetDiffXSec(double xsec, KinePhaseSpace_t ps)
Definition: GHepRecord.h:134
Initial State information.
Definition: InitialState.h:49
#define pDEBUG
Definition: Messenger.h:64