CMS 3D CMS Logo

G4muDarkBremsstrahlungModel.cc
Go to the documentation of this file.
4 
5 //Geant 4
6 #include "G4ProcessTable.hh"
7 #include "G4ProcTblElement.hh"
8 #include "G4MuonMinus.hh"
9 #include "G4MuonPlus.hh"
10 #include "G4ProductionCutsTable.hh"
11 #include "G4SystemOfUnits.hh"
12 //Root
13 #include "TFile.h"
14 #include "TTree.h"
15 //gsl
16 #include <gsl/gsl_math.h>
17 #include <gsl/gsl_integration.h>
18 
19 using namespace std;
20 
22  const G4double biasFactor,
23  const G4ParticleDefinition* p,
24  const G4String& nam)
25  : G4VEmModel(nam),
26  mgfile(scalefile),
27  cxBias(biasFactor),
28  particle(nullptr),
29  isMuon(true),
30  probsup(1.0),
31  isInitialised(false),
32  method("forward_only"),
33  mg_loaded(false) {
34  if (p) {
35  SetParticle(p);
36  } //Verify that the particle is a muon.
38  MA = G4APrime::APrime()->GetPDGMass() / CLHEP::GeV; //Get the A' mass.
39  muonMass = G4MuonMinus::MuonMinusDefinition()->GetPDGMass() / CLHEP::GeV; //Get the muon mass
40  highKinEnergy = HighEnergyLimit();
41  lowKinEnergy = LowEnergyLimit();
42  fParticleChange = nullptr;
43 }
44 
46  size_t n = partialSumSigma.size();
47  if (n > 0) {
48  for (size_t i = 0; i < n; i++) {
49  delete partialSumSigma[i];
50  }
51  }
52 }
53 
54 void G4muDarkBremsstrahlungModel::SetParticle(const G4ParticleDefinition* p) {
55  particle = p;
56 
57  if ((p == G4MuonPlus::MuonPlus()) || (p == G4MuonMinus::MuonMinus())) {
58  isMuon = true;
59  } else {
60  isMuon = false;
61  }
62 }
63 
64 void G4muDarkBremsstrahlungModel::Initialise(const G4ParticleDefinition* p, const G4DataVector& cuts) {
65  if (p) {
66  SetParticle(p);
67  }
68 
69  highKinEnergy = HighEnergyLimit();
70  lowKinEnergy = LowEnergyLimit();
71  const G4ProductionCutsTable* theCoupleTable = G4ProductionCutsTable::GetProductionCutsTable();
72 
73  if (theCoupleTable) {
74  G4int numOfCouples = theCoupleTable->GetTableSize();
75  G4int nn = partialSumSigma.size();
76  G4int nc = cuts.size();
77  if (nn > 0) {
78  for (G4int ii = 0; ii < nn; ii++) {
79  G4DataVector* a = partialSumSigma[ii];
80  if (a)
81  delete a;
82  }
83  partialSumSigma.clear();
84  }
85  if (numOfCouples > 0) {
86  for (G4int i = 0; i < numOfCouples; i++) {
87  G4double cute = DBL_MAX;
88  if (i < nc)
89  cute = cuts[i];
90  const G4MaterialCutsCouple* couple = theCoupleTable->GetMaterialCutsCouple(i);
91  const G4Material* material = couple->GetMaterial();
92  G4DataVector* dv = ComputePartialSumSigma(material, 0.5 * highKinEnergy, std::min(cute, 0.25 * highKinEnergy));
93  partialSumSigma.push_back(dv);
94  }
95  }
96  }
97 
98  if (isInitialised)
99  return;
100  fParticleChange = GetParticleChangeForLoss();
101  isInitialised = true;
102 }
103 
105 //Parses a Madgraph root file to extract the kinetic energy fraction and pt of the outgoing electron in each event. Loads the two numbers from every event into a map of vectors of pairs (mgdata). Map is keyed by energy, vector pairs are energy fraction + pt. Also creates an list of energies and placeholders (energies), so that different energies can be looped separately.
106 {
107  TFile* f = TFile::Open(mgfile.c_str());
108  TTree* tree = (TTree*)f->Get("Events");
109  TLorentzVector* mvec = new TLorentzVector();
110  TLorentzVector* avec = new TLorentzVector();
111  TLorentzVector* nvec = new TLorentzVector();
112  tree->SetBranchAddress("IncidentParticle", &mvec);
113  tree->SetBranchAddress("APrime", &avec);
114  tree->SetBranchAddress("Nucleus", &nvec);
115  int entries = tree->GetEntries();
116  int start = int(G4UniformRand() * entries);
117  for (int i = start; i < int(start + entries / 10.); i++) {
118  if (i < entries) {
119  tree->GetEntry(i);
120  } else {
121  tree->GetEntry(i - entries);
122  }
123  frame evnt;
124  evnt.fEl = (TLorentzVector*)mvec->Clone();
125  evnt.cm = (TLorentzVector*)avec->Clone();
126  *evnt.cm = *evnt.cm + *evnt.fEl;
127  TLorentzVector* ebeam = (TLorentzVector*)nvec->Clone();
128  *ebeam = *ebeam + *evnt.cm;
129  evnt.E = round(ebeam->Z() * 10.) / 10.;
130  if (mgdata.count(evnt.E) == 0) {
131  mgdata[evnt.E];
132  }
133  mgdata[evnt.E].push_back(evnt);
134  }
135  f->Close();
136 }
137 
139  //Need to do this to set up random sampling of mg distributions
140  for (const auto& iter : mgdata) {
141  energies.push_back(std::make_pair(iter.first, iter.second.size()));
142  }
143 
144  for (uint64_t i = 0; i < energies.size(); i++) {
145  energies[i].second = int(G4UniformRand() * mgdata[energies[i].first].size());
146  }
147 }
148 
150 //Gets the energy fraction and Pt from the imported LHE data. E0 should be in GeV, returns the total energy and Pt in GeV. Scales from the closest imported beam energy above the given value (scales down to avoid biasing issues).
151 {
152  frame cmdata;
153  uint64_t i = 0;
154  bool pass = false;
155  //G4double Mel = 5.1E-04;
156  //Cycle through imported beam energies until the closest one above is found, or the max is reached.
157  while (!pass) {
158  i++;
159  if ((E0 <= energies[i].first) || (i >= energies.size())) {
160  pass = true;
161  }
162  }
163 
164  if (i == energies.size()) {
165  i--;
166  } //Decrement if the loop hit the maximum, to prevent a segfault.
167  //energies is a vector of pairs. energies[i].first holds the imported beam energy,
168  //energies[i].second holds the place as we loop through different energy sampling files.
169  //Need to loop around if we hit the end, when the size of mgdata is smaller than our index
170  //on energies[i].second.
171  if (energies[i].second >= double(mgdata[energies[i].first].size())) {
172  energies[i].second = 0;
173  }
174 
175  //Get the lorentz vectors from the index given by the placeholder.
176  cmdata = mgdata[energies[i].first].at(energies[i].second);
177 
178  //Increment the placeholder.
179  energies[i].second = energies[i].second + 1;
180 
181  return cmdata;
182 }
183 
184 G4double G4muDarkBremsstrahlungModel::DsigmaDx(double x, void* pp) {
186 
187  G4double beta = sqrt(1 - (params->MMA) * (params->MMA) / (params->EE0) / (params->EE0));
188  G4double num = 1. - x + x * x / 3.;
189  G4double denom = (params->MMA) * (params->MMA) * (1. - x) / x + (params->MMu) * (params->MMu) * x;
190  G4double DsDx = beta * num / denom;
191 
192  return DsDx;
193 }
194 
195 G4double G4muDarkBremsstrahlungModel::chi(double t, void* pp) {
197  /* Reminder II:
198  * params->AA;
199  * params->ZZ;
200  * params->MMA;
201  * params->MMu;
202  * params->EE0;
203  */
204  G4double Mel = 5.1E-04;
205  G4double MUp = 2.79;
206  G4double Mpr = 0.938;
207  G4double d = 0.164 / pow((params->AA), 2. / 3.);
208  G4double ap = 773.0 / Mel / pow((params->ZZ), 2. / 3.);
209  G4double a = 111.0 / Mel / pow((params->ZZ), 1. / 3.);
210  G4double G2el = (params->ZZ) * (params->ZZ) * a * a * a * a * t * t / (1.0 + a * a * t) / (1.0 + a * a * t) /
211  (1.0 + t / d) / (1.0 + t / d);
212  G4double G2in = (params->ZZ) * ap * ap * ap * ap * t * t / (1.0 + ap * ap * t) / (1.0 + ap * ap * t) /
213  (1.0 + t / 0.71) / (1.0 + t / 0.71) / (1.0 + t / 0.71) / (1.0 + t / 0.71) / (1.0 + t / 0.71) /
214  (1.0 + t / 0.71) / (1.0 + t / 0.71) / (1.0 + t / 0.71) *
215  (1.0 + t * (MUp * MUp - 1.0) / 4.0 / Mpr / Mpr) * (1.0 + t * (MUp * MUp - 1.0) / 4.0 / Mpr / Mpr);
216  G4double G2 = G2el + G2in;
217  G4double ttmin = (params->MMA) * (params->MMA) * (params->MMA) * (params->MMA) / 4.0 / (params->EE0) / (params->EE0);
218  // G4double ttmin = lowerLimit(x,theta,p);
219  G4double Under = G2 * (t - ttmin) / t / t;
220 
221  return Under;
222 }
223 
225  const G4ParticleDefinition*, G4double E0, G4double Z, G4double A, G4double cut, G4double)
226 // Calculates the cross section per atom in GEANT4 internal units. Uses WW approximation to find the total cross section, performing numerical integrals over x and theta.
227 {
228  G4double cross = 0.0;
229  if (E0 < keV || E0 < cut) {
230  return cross;
231  }
232 
233  E0 = E0 / CLHEP::GeV; //Change energy to GeV.
234  if (E0 < 2. * MA)
235  return 0.;
236 
237  //begin: chi-formfactor calculation
238  gsl_integration_workspace* w = gsl_integration_workspace_alloc(1000);
239 
240  G4double result, error;
241  G4double tmin = MA * MA * MA * MA / (4. * E0 * E0);
242  G4double tmax = MA * MA;
243 
244  gsl_function F;
245  ParamsForChi alpha = {1.0, 1.0, 1.0, 1.0, 1.0};
247  F.params = &alpha;
248  alpha.AA = A;
249  alpha.ZZ = Z;
250  alpha.MMA = MA;
251  alpha.MMu = muonMass;
252  alpha.EE0 = E0;
253 
254  //Integrate over chi.
255  gsl_integration_qags(&F, tmin, tmax, 0, 1e-7, 1000, w, &result, &error);
256 
257  G4double ChiRes = result;
258  gsl_integration_workspace_free(w);
259 
260  //Integrate over x. Can use log approximation instead, which falls off at high A' mass.
261  gsl_integration_workspace* dxspace = gsl_integration_workspace_alloc(1000);
262  gsl_function G;
263  G.function = &DsigmaDx;
264  G.params = &alpha;
265  G4double xmin = 0;
266  G4double xmax = muonMass > MA ? 1. - muonMass / E0 : 1. - MA / E0;
267  G4double res, err;
268 
269  gsl_integration_qags(&G, xmin, xmax, 0, 1e-7, 1000, dxspace, &res, &err);
270 
271  G4double DsDx = res;
272  gsl_integration_workspace_free(dxspace);
273 
274  G4double GeVtoPb = 3.894E08;
275  G4double alphaEW = 1.0 / 137.0;
276  G4double epsilBench = 1;
277 
278  cross = GeVtoPb * 4. * alphaEW * alphaEW * alphaEW * epsilBench * epsilBench * ChiRes * DsDx * CLHEP::picobarn;
279  if (cross < 0.) {
280  cross = 0.;
281  }
282 
283  cross = cross * cxBias;
284  return cross;
285 }
286 
287 G4DataVector* G4muDarkBremsstrahlungModel::ComputePartialSumSigma(const G4Material* material,
288  G4double kineticEnergy,
289  G4double cut)
290 
291 // Build the table of cross section per element.
292 //The table is built for MATERIALS.
293 // This table is used by DoIt to select randomly an element in the material.
294 {
295  G4int nElements = material->GetNumberOfElements();
296  const G4ElementVector* theElementVector = material->GetElementVector();
297  const G4double* theAtomNumDensityVector = material->GetAtomicNumDensityVector();
298  G4DataVector* dv = new G4DataVector();
299  G4double cross = 0.0;
300 
301  for (G4int i = 0; i < nElements; i++) {
302  cross += theAtomNumDensityVector[i] *
304  particle, kineticEnergy, (*theElementVector)[i]->GetZ(), (*theElementVector)[i]->GetN(), cut);
305  dv->push_back(cross);
306  }
307  return dv;
308 }
309 
310 void G4muDarkBremsstrahlungModel::SampleSecondaries(std::vector<G4DynamicParticle*>* vdp,
311  const G4MaterialCutsCouple* couple,
312  const G4DynamicParticle* dp,
313  G4double tmin,
314  G4double maxEnergy)
315 //Simulates the emission of a dark photon + electron. Gets an energy fraction and Pt from madgraph files. Scales the energy so that the fraction of kinectic energy is constant, keeps the Pt constant. If the Pt is larger than the new energy, that event is skipped, and a new one is taken from the file.
316 {
317  //Deactivate the process after one dark brem. Needs to be reactivated in the end of event action. If this is in the stepping action instead, more than one brem can occur within each step.
318  G4bool state = false;
319  G4String pname = "muDBrem";
320  G4ProcessTable* ptable = G4ProcessTable::GetProcessTable();
321  ptable->SetProcessActivation(pname, state);
322  if (mg_loaded == false) {
323  LoadMG();
324  MakePlaceholders(); //Setup the placeholder offsets for getting data.
325  mg_loaded = true;
326  }
327  G4double E0 = dp->GetTotalEnergy();
328  G4double tmax = min(maxEnergy, E0);
329  if (tmin >= tmax) {
330  return;
331  } // limits of the energy sampling
332 
333  E0 = E0 / CLHEP::GeV; //Convert the energy to GeV, the units used in the sampling files.
334 
335  frame data = GetMadgraphData(E0);
336  double EAcc, Pt, P, PhiAcc;
337  if (method == "forward_only") {
338  EAcc = (data.fEl->E() - muonMass) / (data.E - muonMass - MA) * (E0 - muonMass - MA);
339  EAcc = muonMass + EAcc;
340  Pt = data.fEl->Pt();
341  P = sqrt(EAcc * EAcc - muonMass * muonMass);
342  PhiAcc = data.fEl->Phi();
343  int i = 0;
344  while (Pt * Pt + muonMass * muonMass > EAcc * EAcc) //Skip events until the Pt is less than the energy.
345  {
346  i++;
347  data = GetMadgraphData(E0);
348  EAcc = (data.fEl->E() - muonMass) / (data.E - muonMass - MA) * (E0 - muonMass - MA);
349  EAcc = muonMass + EAcc;
350  Pt = data.fEl->Pt();
351  P = sqrt(EAcc * EAcc - muonMass * muonMass);
352  PhiAcc = data.fEl->Phi();
353 
354  if (i > 10000) {
355  break;
356  }
357  }
358  } else if (method == "cm_scaling") {
359  TLorentzVector* el = new TLorentzVector(data.fEl->X(), data.fEl->Y(), data.fEl->Z(), data.fEl->E());
360  double ediff = data.E - E0;
361  TLorentzVector* newcm = new TLorentzVector(data.cm->X(), data.cm->Y(), data.cm->Z() - ediff, data.cm->E() - ediff);
362  el->Boost(-1. * data.cm->BoostVector());
363  el->Boost(newcm->BoostVector());
364  double newE = (data.fEl->E() - muonMass) / (data.E - muonMass - MA) * (E0 - muonMass - MA);
365  el->SetE(newE);
366  EAcc = el->E();
367  Pt = el->Pt();
368  P = el->P();
369  PhiAcc = el->Phi();
370  } else {
371  EAcc = E0;
372  P = dp->GetTotalMomentum();
373  Pt = sqrt(dp->Get4Momentum().px() * dp->Get4Momentum().px() + dp->Get4Momentum().py() * dp->Get4Momentum().py());
374  PhiAcc = dp->Get4Momentum().phi();
375  }
376 
377  EAcc = EAcc * CLHEP::GeV; //Change the energy back to MeV, the internal GEANT unit.
378 
379  G4double muon_mass_mev = muonMass * CLHEP::GeV;
380  G4double momentum = sqrt(EAcc * EAcc - muon_mass_mev * muon_mass_mev); //Muon momentum in MeV.
381  G4ThreeVector newDirection;
382  double ThetaAcc = std::asin(Pt / P);
383  newDirection.set(std::sin(ThetaAcc) * std::cos(PhiAcc), std::sin(ThetaAcc) * std::sin(PhiAcc), std::cos(ThetaAcc));
384  newDirection.rotateUz(dp->GetMomentumDirection());
385  newDirection.setMag(momentum);
386  // create g4dynamicparticle object for the dark photon.
387  G4ThreeVector direction = (dp->GetMomentum() - newDirection);
388  G4DynamicParticle* dphoton = new G4DynamicParticle(theAPrime, direction);
389  vdp->push_back(dphoton);
390 
391  // energy of primary
392  G4double finalKE = EAcc - muon_mass_mev;
393 
394  // stop tracking and create new secondary instead of primary
395  bool makeSecondary = true;
396  if (makeSecondary) {
397  fParticleChange->ProposeTrackStatus(fStopAndKill);
398  fParticleChange->SetProposedKineticEnergy(0.0);
399  G4DynamicParticle* mu =
400  new G4DynamicParticle(const_cast<G4ParticleDefinition*>(particle), newDirection.unit(), finalKE);
401  vdp->push_back(mu);
402  // continue tracking
403  } else {
404  fParticleChange->SetProposedMomentumDirection(newDirection.unit());
405  fParticleChange->SetProposedKineticEnergy(finalKE);
406  }
407 }
408 
409 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
410 
411 const G4Element* G4muDarkBremsstrahlungModel::SelectRandomAtom(const G4MaterialCutsCouple* couple) {
412  // select randomly 1 element within the material
413 
414  const G4Material* material = couple->GetMaterial();
415  G4int nElements = material->GetNumberOfElements();
416  const G4ElementVector* theElementVector = material->GetElementVector();
417 
418  const G4Element* elm = nullptr;
419 
420  if (1 < nElements) {
421  --nElements;
422  G4DataVector* dv = partialSumSigma[couple->GetIndex()];
423  G4double rval = G4UniformRand() * ((*dv)[nElements]);
424 
425  elm = (*theElementVector)[nElements];
426  for (G4int i = 0; i < nElements; ++i) {
427  if (rval <= (*dv)[i]) {
428  elm = (*theElementVector)[i];
429  break;
430  }
431  }
432  } else {
433  elm = (*theElementVector)[0];
434  }
435 
436  SetCurrentElement(elm);
437  return elm;
438 }
439 
441  method = method_in;
442  return;
443 }
size
Write out results.
std::map< double, std::vector< frame > > mgdata
Definition: start.py:1
Basic3DVector cross(const Basic3DVector &v) const
Vector product, or "cross" product, with a vector of same type.
std::vector< std::pair< double, int > > energies
bool isMuon(const Candidate &part)
Definition: pdgIdUtils.h:9
T w() const
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double tkin, G4double Z, G4double, G4double cut, G4double maxE=DBL_MAX) override
G4ParticleChangeForLoss * fParticleChange
constexpr int pow(int x)
Definition: conifer.h:24
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
Class creating the A&#39; particle in Geant.
Definition: Electron.h:6
U second(std::pair< T, U > const &p)
T sqrt(T t)
Definition: SSEVec.h:19
Class provided to simulate the dark brem cross section and interaction.
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
double f[11][100]
d
Definition: ztail.py:151
static G4double chi(double t, void *pp)
ii
Definition: cuy.py:589
G4muDarkBremsstrahlungModel(const G4String &scalefile, const G4double biasFactor, const G4ParticleDefinition *p=nullptr, const G4String &nam="eDBrem")
void SampleSecondaries(std::vector< G4DynamicParticle *> *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
static const double tmax[3]
unsigned long long uint64_t
Definition: Time.h:13
TLorentzVector * cm
static G4double DsigmaDx(double x, void *pp)
Class providing the Dark Bremsstrahlung process class.
std::pair< OmniClusterRef, TrackingParticleRef > P
TLorentzVector * fEl
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:80
std::vector< G4DataVector * > partialSumSigma
double a
Definition: hdecay.h:121
static G4APrime * APrime(double apmass=1000)
Definition: G4APrime.cc:43
float x
void SetParticle(const G4ParticleDefinition *p)
Definition: APVGainStruct.h:7
Definition: tree.py:1
static uInt32 F(BLOWFISH_CTX *ctx, uInt32 x)
Definition: blowfish.cc:163
G4DataVector * ComputePartialSumSigma(const G4Material *material, G4double tkin, G4double cut)
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *couple)
const G4ParticleDefinition * particle