CMS 3D CMS Logo

List of all members | Public Member Functions | Protected Member Functions | Protected Attributes | Private Member Functions | Static Private Member Functions | Private Attributes
G4muDarkBremsstrahlungModel Class Reference

#include <G4muDarkBremsstrahlungModel.h>

Inheritance diagram for G4muDarkBremsstrahlungModel:

Public Member Functions

G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double tkin, G4double Z, G4double, G4double cut, G4double maxE=DBL_MAX) override
 
G4DataVector * ComputePartialSumSigma (const G4Material *material, G4double tkin, G4double cut)
 
 G4muDarkBremsstrahlungModel (const G4String &scalefile, const G4double biasFactor, const G4ParticleDefinition *p=nullptr, const G4String &nam="eDBrem")
 
 G4muDarkBremsstrahlungModel (const G4muDarkBremsstrahlungModel &)=delete
 
frame GetMadgraphData (double E0)
 
void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 
void LoadMG ()
 
void MakePlaceholders ()
 
G4muDarkBremsstrahlungModeloperator= (const G4muDarkBremsstrahlungModel &right)=delete
 
void SampleSecondaries (std::vector< G4DynamicParticle *> *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
 
void SetMethod (std::string)
 
 ~G4muDarkBremsstrahlungModel () override
 

Protected Member Functions

const G4Element * SelectRandomAtom (const G4MaterialCutsCouple *couple)
 

Protected Attributes

const G4double cxBias
 
G4ParticleChangeForLoss * fParticleChange
 
G4bool isMuon
 
G4double MA
 
const G4String & mgfile
 
G4double muonMass
 
const G4ParticleDefinition * particle
 
G4ParticleDefinition * theAPrime
 

Private Member Functions

void SetParticle (const G4ParticleDefinition *p)
 

Static Private Member Functions

static G4double chi (double t, void *pp)
 
static G4double DsigmaDx (double x, void *pp)
 

Private Attributes

std::vector< std::pair< double, int > > energies
 
G4double highKinEnergy
 
G4bool isInitialised
 
G4double lowKinEnergy
 
std::string method
 
G4bool mg_loaded
 
std::map< double, std::vector< frame > > mgdata
 
std::vector< G4DataVector * > partialSumSigma
 
G4double probsup
 

Detailed Description

Definition at line 36 of file G4muDarkBremsstrahlungModel.h.

Constructor & Destructor Documentation

◆ G4muDarkBremsstrahlungModel() [1/2]

G4muDarkBremsstrahlungModel::G4muDarkBremsstrahlungModel ( const G4String &  scalefile,
const G4double  biasFactor,
const G4ParticleDefinition *  p = nullptr,
const G4String &  nam = "eDBrem" 
)

Definition at line 21 of file G4muDarkBremsstrahlungModel.cc.

References G4APrime::APrime(), fParticleChange, highKinEnergy, lowKinEnergy, MA, muonMass, AlCaHLTBitMon_ParallelJobs::p, SetParticle(), and theAPrime.

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 }
G4ParticleChangeForLoss * fParticleChange
static G4APrime * APrime(double apmass=1000)
Definition: G4APrime.cc:43
void SetParticle(const G4ParticleDefinition *p)
const G4ParticleDefinition * particle

◆ ~G4muDarkBremsstrahlungModel()

G4muDarkBremsstrahlungModel::~G4muDarkBremsstrahlungModel ( )
override

Definition at line 45 of file G4muDarkBremsstrahlungModel.cc.

References mps_fire::i, dqmiodumpmetadata::n, and partialSumSigma.

45  {
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 }
std::vector< G4DataVector * > partialSumSigma

◆ G4muDarkBremsstrahlungModel() [2/2]

G4muDarkBremsstrahlungModel::G4muDarkBremsstrahlungModel ( const G4muDarkBremsstrahlungModel )
delete

Member Function Documentation

◆ chi()

G4double G4muDarkBremsstrahlungModel::chi ( double  t,
void *  pp 
)
staticprivate

Definition at line 197 of file G4muDarkBremsstrahlungModel.cc.

References a, ztail::d, submitPVValidationJobs::params, funct::pow(), createTree::pp, and submitPVValidationJobs::t.

Referenced by ComputeCrossSectionPerAtom().

197  {
199  /* Reminder II:
200  * params->AA;
201  * params->ZZ;
202  * params->MMA;
203  * params->MMu;
204  * params->EE0;
205  */
206  G4double Mel = 5.1E-04;
207  G4double MUp = 2.79;
208  G4double Mpr = 0.938;
209  G4double d = 0.164 / pow((params->AA), 2. / 3.);
210  G4double ap = 773.0 / Mel / pow((params->ZZ), 2. / 3.);
211  G4double a = 111.0 / Mel / pow((params->ZZ), 1. / 3.);
212  G4double G2el = (params->ZZ) * (params->ZZ) * a * a * a * a * t * t / (1.0 + a * a * t) / (1.0 + a * a * t) /
213  (1.0 + t / d) / (1.0 + t / d);
214  G4double G2in = (params->ZZ) * ap * ap * ap * ap * t * t / (1.0 + ap * ap * t) / (1.0 + ap * ap * t) /
215  (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) /
216  (1.0 + t / 0.71) / (1.0 + t / 0.71) / (1.0 + t / 0.71) *
217  (1.0 + t * (MUp * MUp - 1.0) / 4.0 / Mpr / Mpr) * (1.0 + t * (MUp * MUp - 1.0) / 4.0 / Mpr / Mpr);
218  G4double G2 = G2el + G2in;
219  G4double ttmin = (params->MMA) * (params->MMA) * (params->MMA) * (params->MMA) / 4.0 / (params->EE0) / (params->EE0);
220  // G4double ttmin = lowerLimit(x,theta,p);
221  G4double Under = G2 * (t - ttmin) / t / t;
222 
223  return Under;
224 }
d
Definition: ztail.py:151
double a
Definition: hdecay.h:119
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29

◆ ComputeCrossSectionPerAtom()

G4double G4muDarkBremsstrahlungModel::ComputeCrossSectionPerAtom ( const G4ParticleDefinition *  ,
G4double  tkin,
G4double  Z,
G4double  A,
G4double  cut,
G4double  maxE = DBL_MAX 
)
override

Definition at line 226 of file G4muDarkBremsstrahlungModel.cc.

References A, alpha, chi(), cross(), PA_MinBiasSkim_cff::cut, cxBias, DsigmaDx(), MillePedeFileConverter_cfg::e, submitPVResolutionJobs::err, relativeConstraints::error, F(), cmssw_cycle_finder::G, MA, muonMass, mps_fire::result, tmax, muonTiming_cfi::tmin, w(), TrackerOfflineValidation_Dqm_cff::xmax, TrackerOfflineValidation_Dqm_cff::xmin, and BeamSpotPI::Z.

Referenced by ComputePartialSumSigma().

229 {
230  G4double cross = 0.0;
231  if (E0 < keV || E0 < cut) {
232  return cross;
233  }
234 
235  E0 = E0 / CLHEP::GeV; //Change energy to GeV.
236  if (E0 < 2. * MA)
237  return 0.;
238 
239  //begin: chi-formfactor calculation
240  gsl_integration_workspace* w = gsl_integration_workspace_alloc(1000);
241 
242  G4double result, error;
243  G4double tmin = MA * MA * MA * MA / (4. * E0 * E0);
244  G4double tmax = MA * MA;
245 
246  gsl_function F;
247  ParamsForChi alpha = {1.0, 1.0, 1.0, 1.0, 1.0};
249  F.params = &alpha;
250  alpha.AA = A;
251  alpha.ZZ = Z;
252  alpha.MMA = MA;
253  alpha.MMu = muonMass;
254  alpha.EE0 = E0;
255 
256  //Integrate over chi.
257  gsl_integration_qags(&F, tmin, tmax, 0, 1e-7, 1000, w, &result, &error);
258 
259  G4double ChiRes = result;
260  gsl_integration_workspace_free(w);
261 
262  //Integrate over x. Can use log approximation instead, which falls off at high A' mass.
263  gsl_integration_workspace* dxspace = gsl_integration_workspace_alloc(1000);
264  gsl_function G;
265  G.function = &DsigmaDx;
266  G.params = &alpha;
267  G4double xmin = 0;
268  G4double xmax = 1;
269  if ((muonMass / E0) > (MA / E0))
270  xmax = 1 - muonMass / E0;
271  else
272  xmax = 1 - MA / E0;
273  G4double res, err;
274 
275  gsl_integration_qags(&G, xmin, xmax, 0, 1e-7, 1000, dxspace, &res, &err);
276 
277  G4double DsDx = res;
278  gsl_integration_workspace_free(dxspace);
279 
280  G4double GeVtoPb = 3.894E08;
281  G4double alphaEW = 1.0 / 137.0;
282  G4double epsilBench = 1;
283 
284  cross = GeVtoPb * 4. * alphaEW * alphaEW * alphaEW * epsilBench * epsilBench * ChiRes * DsDx * CLHEP::picobarn;
285  if (cross < 0.) {
286  cross = 0.;
287  }
288  E0 = E0 * CLHEP::GeV;
289 
290  cross = cross * cxBias;
291  return cross;
292 }
float alpha
Definition: AMPTWrapper.h:105
Basic3DVector cross(const Basic3DVector &v) const
Vector product, or "cross" product, with a vector of same type.
T w() const
Definition: Electron.h:6
static G4double chi(double t, void *pp)
static const double tmax[3]
static G4double DsigmaDx(double x, void *pp)
Definition: APVGainStruct.h:7
static uInt32 F(BLOWFISH_CTX *ctx, uInt32 x)
Definition: blowfish.cc:163

◆ ComputePartialSumSigma()

G4DataVector * G4muDarkBremsstrahlungModel::ComputePartialSumSigma ( const G4Material *  material,
G4double  tkin,
G4double  cut 
)

Definition at line 294 of file G4muDarkBremsstrahlungModel.cc.

References ComputeCrossSectionPerAtom(), cross(), PA_MinBiasSkim_cff::cut, mps_fire::i, and particle.

Referenced by Initialise().

301 {
302  G4int nElements = material->GetNumberOfElements();
303  const G4ElementVector* theElementVector = material->GetElementVector();
304  const G4double* theAtomNumDensityVector = material->GetAtomicNumDensityVector();
305  G4DataVector* dv = new G4DataVector();
306  G4double cross = 0.0;
307 
308  for (G4int i = 0; i < nElements; i++) {
309  cross += theAtomNumDensityVector[i] *
311  particle, kineticEnergy, (*theElementVector)[i]->GetZ(), (*theElementVector)[i]->GetN(), cut);
312  dv->push_back(cross);
313  }
314  return dv;
315 }
Basic3DVector cross(const Basic3DVector &v) const
Vector product, or "cross" product, with a vector of same type.
G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double tkin, G4double Z, G4double, G4double cut, G4double maxE=DBL_MAX) override
const G4ParticleDefinition * particle

◆ DsigmaDx()

G4double G4muDarkBremsstrahlungModel::DsigmaDx ( double  x,
void *  pp 
)
staticprivate

Definition at line 186 of file G4muDarkBremsstrahlungModel.cc.

References HLT_2022v14_cff::beta, makePileupJSON::denom, EgammaValidation_cff::num, submitPVValidationJobs::params, createTree::pp, mathSSE::sqrt(), and x.

Referenced by ComputeCrossSectionPerAtom().

186  {
188 
189  G4double beta = sqrt(1 - (params->MMA) * (params->MMA) / (params->EE0) / (params->EE0));
190  G4double num = 1. - x + x * x / 3.;
191  G4double denom = (params->MMA) * (params->MMA) * (1. - x) / x + (params->MMu) * (params->MMu) * x;
192  G4double DsDx = beta * num / denom;
193 
194  return DsDx;
195 }
T sqrt(T t)
Definition: SSEVec.h:19

◆ GetMadgraphData()

frame G4muDarkBremsstrahlungModel::GetMadgraphData ( double  E0)

Definition at line 149 of file G4muDarkBremsstrahlungModel.cc.

References energies, first, mps_fire::i, mgdata, edm::second(), and findQualityFiles::size.

Referenced by SampleSecondaries().

151 {
152  double samplingE = energies[0].first;
153  frame cmdata;
154  uint64_t i = 0;
155  bool pass = false;
156  //G4double Mel = 5.1E-04;
157  //Cycle through imported beam energies until the closest one above is found, or the max is reached.
158  while (!pass) {
159  i++;
160  samplingE = energies[i].first;
161  if ((E0 <= samplingE) || (i >= energies.size())) {
162  pass = true;
163  }
164  }
165 
166  if (i == energies.size()) {
167  i--;
168  } //Decrement if the loop hit the maximum, to prevent a segfault.
169  //energies is a vector of pairs. energies[i].first holds the imported beam energy,
170  //energies[i].second holds the place as we loop through different energy sampling files.
171  //Need to loop around if we hit the end, when the size of mgdata is smaller than our index
172  //on energies[i].second.
173  if (energies[i].second >= double(mgdata[energies[i].first].size())) {
174  energies[i].second = 0;
175  }
176 
177  //Get the lorentz vectors from the index given by the placeholder.
178  cmdata = mgdata[energies[i].first].at(energies[i].second);
179 
180  //Increment the placeholder.
181  energies[i].second = energies[i].second + 1;
182 
183  return cmdata;
184 }
size
Write out results.
std::map< double, std::vector< frame > > mgdata
std::vector< std::pair< double, int > > energies
U second(std::pair< T, U > const &p)
unsigned long long uint64_t
Definition: Time.h:13

◆ Initialise()

void G4muDarkBremsstrahlungModel::Initialise ( const G4ParticleDefinition *  p,
const G4DataVector &  cuts 
)
override

Definition at line 64 of file G4muDarkBremsstrahlungModel.cc.

References a, ComputePartialSumSigma(), cuts, fParticleChange, highKinEnergy, mps_fire::i, cuy::ii, isInitialised, lowKinEnergy, SiStripPI::min, groupFilesInBlocks::nn, AlCaHLTBitMon_ParallelJobs::p, partialSumSigma, and SetParticle().

64  {
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 }
G4ParticleChangeForLoss * fParticleChange
TkSoA const *__restrict__ CAHitNtupletGeneratorKernelsGPU::QualityCuts cuts
ii
Definition: cuy.py:589
std::vector< G4DataVector * > partialSumSigma
double a
Definition: hdecay.h:119
void SetParticle(const G4ParticleDefinition *p)
G4DataVector * ComputePartialSumSigma(const G4Material *material, G4double tkin, G4double cut)

◆ LoadMG()

void G4muDarkBremsstrahlungModel::LoadMG ( )

Definition at line 104 of file G4muDarkBremsstrahlungModel.cc.

References frame::cm, frame::E, f, frame::fEl, mps_fire::i, createfilelist::int, mgdata, and mgfile.

Referenced by SampleSecondaries().

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 }
std::map< double, std::vector< frame > > mgdata
Definition: start.py:1
double f[11][100]
TLorentzVector * cm
TLorentzVector * fEl
Definition: tree.py:1

◆ MakePlaceholders()

void G4muDarkBremsstrahlungModel::MakePlaceholders ( )

Definition at line 138 of file G4muDarkBremsstrahlungModel.cc.

References energies, first, mps_fire::i, createfilelist::int, mgdata, and findQualityFiles::size.

Referenced by SampleSecondaries().

138  {
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 }
size
Write out results.
std::map< double, std::vector< frame > > mgdata
std::vector< std::pair< double, int > > energies
unsigned long long uint64_t
Definition: Time.h:13

◆ operator=()

G4muDarkBremsstrahlungModel& G4muDarkBremsstrahlungModel::operator= ( const G4muDarkBremsstrahlungModel right)
delete

◆ SampleSecondaries()

void G4muDarkBremsstrahlungModel::SampleSecondaries ( std::vector< G4DynamicParticle *> *  vdp,
const G4MaterialCutsCouple *  couple,
const G4DynamicParticle *  dp,
G4double  tmin,
G4double  maxEnergy 
)
override

Definition at line 317 of file G4muDarkBremsstrahlungModel.cc.

References funct::cos(), data, Calorimetry_cff::dp, fParticleChange, GetMadgraphData(), mps_fire::i, LoadMG(), MA, MakePlaceholders(), particleFlowClusterECALTimeSelected_cfi::maxEnergy, method, mg_loaded, SiStripPI::min, amptDefaultParameters_cff::mu, muonMass, particle, unpackData-CaloStage2::pname, funct::sin(), mathSSE::sqrt(), theAPrime, tmax, and muonTiming_cfi::tmin.

323 {
324  //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.
325  G4bool state = false;
326  G4String pname = "muDBrem";
327  G4ProcessTable* ptable = G4ProcessTable::GetProcessTable();
328  ptable->SetProcessActivation(pname, state);
329  if (mg_loaded == false) {
330  LoadMG();
331  MakePlaceholders(); //Setup the placeholder offsets for getting data.
332  mg_loaded = true;
333  }
334  G4double E0 = dp->GetTotalEnergy();
335  G4double tmax = min(maxEnergy, E0);
336  if (tmin >= tmax) {
337  return;
338  } // limits of the energy sampling
339 
340  E0 = E0 / CLHEP::GeV; //Convert the energy to GeV, the units used in the sampling files.
341 
342  frame data = GetMadgraphData(E0);
343  double EAcc, Pt, P, PhiAcc;
344  if (method == "forward_only") {
345  EAcc = (data.fEl->E() - muonMass) / (data.E - muonMass - MA) * (E0 - muonMass - MA);
346  EAcc = muonMass + EAcc;
347  Pt = data.fEl->Pt();
348  P = sqrt(EAcc * EAcc - muonMass * muonMass);
349  PhiAcc = data.fEl->Phi();
350  int i = 0;
351  while (Pt * Pt + muonMass * muonMass > EAcc * EAcc) //Skip events until the Pt is less than the energy.
352  {
353  i++;
354  data = GetMadgraphData(E0);
355  EAcc = (data.fEl->E() - muonMass) / (data.E - muonMass - MA) * (E0 - muonMass - MA);
356  EAcc = muonMass + EAcc;
357  Pt = data.fEl->Pt();
358  P = sqrt(EAcc * EAcc - muonMass * muonMass);
359  PhiAcc = data.fEl->Phi();
360 
361  if (i > 10000) {
362  break;
363  }
364  }
365  } else if (method == "cm_scaling") {
366  TLorentzVector* el = new TLorentzVector(data.fEl->X(), data.fEl->Y(), data.fEl->Z(), data.fEl->E());
367  double ediff = data.E - E0;
368  TLorentzVector* newcm = new TLorentzVector(data.cm->X(), data.cm->Y(), data.cm->Z() - ediff, data.cm->E() - ediff);
369  el->Boost(-1. * data.cm->BoostVector());
370  el->Boost(newcm->BoostVector());
371  double newE = (data.fEl->E() - muonMass) / (data.E - muonMass - MA) * (E0 - muonMass - MA);
372  el->SetE(newE);
373  EAcc = el->E();
374  Pt = el->Pt();
375  P = el->P();
376  PhiAcc = el->Phi();
377  } else {
378  EAcc = E0;
379  P = dp->GetTotalMomentum();
380  Pt = sqrt(dp->Get4Momentum().px() * dp->Get4Momentum().px() + dp->Get4Momentum().py() * dp->Get4Momentum().py());
381  PhiAcc = dp->Get4Momentum().phi();
382  }
383 
384  EAcc = EAcc * CLHEP::GeV; //Change the energy back to MeV, the internal GEANT unit.
385 
386  G4double muon_mass_mev = muonMass * CLHEP::GeV;
387  G4double momentum = sqrt(EAcc * EAcc - muon_mass_mev * muon_mass_mev); //Muon momentum in MeV.
388  G4ThreeVector newDirection;
389  double ThetaAcc = std::asin(Pt / P);
390  newDirection.set(std::sin(ThetaAcc) * std::cos(PhiAcc), std::sin(ThetaAcc) * std::sin(PhiAcc), std::cos(ThetaAcc));
391  newDirection.rotateUz(dp->GetMomentumDirection());
392  newDirection.setMag(momentum);
393  // create g4dynamicparticle object for the dark photon.
394  G4ThreeVector direction = (dp->GetMomentum() - newDirection);
395  G4DynamicParticle* dphoton = new G4DynamicParticle(theAPrime, direction);
396  vdp->push_back(dphoton);
397 
398  // energy of primary
399  G4double finalKE = EAcc - muon_mass_mev;
400 
401  // stop tracking and create new secondary instead of primary
402  bool makeSecondary = true;
403  if (makeSecondary) {
404  fParticleChange->ProposeTrackStatus(fStopAndKill);
405  fParticleChange->SetProposedKineticEnergy(0.0);
406  G4DynamicParticle* mu =
407  new G4DynamicParticle(const_cast<G4ParticleDefinition*>(particle), newDirection.unit(), finalKE);
408  vdp->push_back(mu);
409  // continue tracking
410  } else {
411  fParticleChange->SetProposedMomentumDirection(newDirection.unit());
412  fParticleChange->SetProposedKineticEnergy(finalKE);
413  }
414 }
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
G4ParticleChangeForLoss * fParticleChange
T sqrt(T t)
Definition: SSEVec.h:19
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
static const double tmax[3]
std::pair< OmniClusterRef, TrackingParticleRef > P
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:79
const G4ParticleDefinition * particle

◆ SelectRandomAtom()

const G4Element * G4muDarkBremsstrahlungModel::SelectRandomAtom ( const G4MaterialCutsCouple *  couple)
protected

Definition at line 418 of file G4muDarkBremsstrahlungModel.cc.

References mps_fire::i, and partialSumSigma.

418  {
419  // select randomly 1 element within the material
420 
421  const G4Material* material = couple->GetMaterial();
422  G4int nElements = material->GetNumberOfElements();
423  const G4ElementVector* theElementVector = material->GetElementVector();
424 
425  const G4Element* elm = nullptr;
426 
427  if (1 < nElements) {
428  --nElements;
429  G4DataVector* dv = partialSumSigma[couple->GetIndex()];
430  G4double rval = G4UniformRand() * ((*dv)[nElements]);
431 
432  elm = (*theElementVector)[nElements];
433  for (G4int i = 0; i < nElements; ++i) {
434  if (rval <= (*dv)[i]) {
435  elm = (*theElementVector)[i];
436  break;
437  }
438  }
439  } else {
440  elm = (*theElementVector)[0];
441  }
442 
443  SetCurrentElement(elm);
444  return elm;
445 }
std::vector< G4DataVector * > partialSumSigma

◆ SetMethod()

void G4muDarkBremsstrahlungModel::SetMethod ( std::string  method_in)

Definition at line 447 of file G4muDarkBremsstrahlungModel.cc.

References method.

447  {
448  method = method_in;
449  return;
450 }

◆ SetParticle()

void G4muDarkBremsstrahlungModel::SetParticle ( const G4ParticleDefinition *  p)
private

Definition at line 54 of file G4muDarkBremsstrahlungModel.cc.

References isMuon, AlCaHLTBitMon_ParallelJobs::p, and particle.

Referenced by G4muDarkBremsstrahlungModel(), and Initialise().

54  {
55  particle = p;
56 
57  if ((p == G4MuonPlus::MuonPlus()) || (p == G4MuonMinus::MuonMinus())) {
58  isMuon = true;
59  } else {
60  isMuon = false;
61  }
62 }
const G4ParticleDefinition * particle

Member Data Documentation

◆ cxBias

const G4double G4muDarkBremsstrahlungModel::cxBias
protected

Definition at line 80 of file G4muDarkBremsstrahlungModel.h.

Referenced by ComputeCrossSectionPerAtom().

◆ energies

std::vector<std::pair<double, int> > G4muDarkBremsstrahlungModel::energies
private

Definition at line 96 of file G4muDarkBremsstrahlungModel.h.

Referenced by GetMadgraphData(), and MakePlaceholders().

◆ fParticleChange

G4ParticleChangeForLoss* G4muDarkBremsstrahlungModel::fParticleChange
protected

◆ highKinEnergy

G4double G4muDarkBremsstrahlungModel::highKinEnergy
private

Definition at line 89 of file G4muDarkBremsstrahlungModel.h.

Referenced by G4muDarkBremsstrahlungModel(), and Initialise().

◆ isInitialised

G4bool G4muDarkBremsstrahlungModel::isInitialised
private

Definition at line 92 of file G4muDarkBremsstrahlungModel.h.

Referenced by Initialise().

◆ isMuon

G4bool G4muDarkBremsstrahlungModel::isMuon
protected

Definition at line 86 of file G4muDarkBremsstrahlungModel.h.

Referenced by SetParticle().

◆ lowKinEnergy

G4double G4muDarkBremsstrahlungModel::lowKinEnergy
private

Definition at line 90 of file G4muDarkBremsstrahlungModel.h.

Referenced by G4muDarkBremsstrahlungModel(), and Initialise().

◆ MA

G4double G4muDarkBremsstrahlungModel::MA
protected

◆ method

std::string G4muDarkBremsstrahlungModel::method
private

Definition at line 93 of file G4muDarkBremsstrahlungModel.h.

Referenced by SampleSecondaries(), and SetMethod().

◆ mg_loaded

G4bool G4muDarkBremsstrahlungModel::mg_loaded
private

Definition at line 94 of file G4muDarkBremsstrahlungModel.h.

Referenced by SampleSecondaries().

◆ mgdata

std::map<double, std::vector<frame> > G4muDarkBremsstrahlungModel::mgdata
private

Definition at line 95 of file G4muDarkBremsstrahlungModel.h.

Referenced by GetMadgraphData(), LoadMG(), and MakePlaceholders().

◆ mgfile

const G4String& G4muDarkBremsstrahlungModel::mgfile
protected

Definition at line 79 of file G4muDarkBremsstrahlungModel.h.

Referenced by LoadMG().

◆ muonMass

G4double G4muDarkBremsstrahlungModel::muonMass
protected

◆ partialSumSigma

std::vector<G4DataVector*> G4muDarkBremsstrahlungModel::partialSumSigma
private

◆ particle

const G4ParticleDefinition* G4muDarkBremsstrahlungModel::particle
protected

◆ probsup

G4double G4muDarkBremsstrahlungModel::probsup
private

Definition at line 91 of file G4muDarkBremsstrahlungModel.h.

◆ theAPrime

G4ParticleDefinition* G4muDarkBremsstrahlungModel::theAPrime
protected

Definition at line 82 of file G4muDarkBremsstrahlungModel.h.

Referenced by G4muDarkBremsstrahlungModel(), and SampleSecondaries().