ATLAS Offline Software
Loading...
Searching...
No Matches
TRTTransitionRadiation.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2026 CERN for the benefit of the ATLAS collaboration
3*/
4
5//#define ARTRU // Choice of TR generator
6
7// class header
9
10// package includes
12#include "TRRegionXMLHandler.h"
13
14// Athena includes
15#include "CLHEP/Random/RandGaussZiggurat.h"
16#include "GeoModelInterfaces/StoredMaterialManager.h" //Material Manager
17#include "GeoMaterial2G4/Geo2G4MaterialFactory.h" //Converting GeoMaterial -> G4Material
19#include "StoreGate/StoreGateSvc.h" //Detector Store
20
25
26// Geant4 includes
27#include "G4DynamicParticle.hh"
28#include "G4Gamma.hh"
29#include "G4Material.hh"
30#include "G4Element.hh"
31#include "G4EmProcessSubType.hh" // for fTransitionRadiation
32#include "G4ProcessType.hh" // for fElectromagnetic
33
34#include "Randomize.hh"
35
36// CLHEP includes
37#include "CLHEP/Random/RandPoisson.h"
38
39// STL includes
40#include <cmath>
41#include <limits>
42
44//
45// Constructor, destructor
46
47TRTTransitionRadiation::TRTTransitionRadiation( const G4String& processName, const std::string & xmlfilename) :
48 G4VDiscreteProcess(processName,fElectromagnetic),
49 AthMessaging("TRTTransitionRadiation"),
50 m_XMLhandler(nullptr),m_xmlfilename(xmlfilename),
52 m_WplasmaFoil(0.0),m_GammaMin(0.0),m_EkinMin(0.0),m_Ey(nullptr),m_Sr(nullptr),
53 m_om(nullptr),m_Omg(nullptr),m_sigmaGas(nullptr),m_sigmaFoil(nullptr)
54{
55 m_radiators.clear();
56 m_XMLhandler = new TRRegionXMLHandler( this );
57
58 SetProcessSubType( fTransitionRadiation );
59
60 ATH_MSG_DEBUG ( "##### Constructor TRTTransitionRadiation" );
61
62 Initialize();
63
64 ATH_MSG_DEBUG ( "##### Constructor TRTTransitionRadiation done" );
65
66}
67
69
71 delete [] m_Ey;
72 delete [] m_Sr;
73 delete [] m_Omg;
74 delete [] m_om;
75 delete [] m_sigmaGas;
76 delete [] m_sigmaFoil;
77}
78
81
82 ATH_MSG_DEBUG(" New Radiator parameters being defined for TR process");
83 ATH_MSG_DEBUG(" Volume " << p.GetLogicalVolume()->GetName());
84 ATH_MSG_DEBUG(" FoilThickness " << p.GetFoilThickness());
85 ATH_MSG_DEBUG(" GasThickness " << p.GetGasThickness());
86
87 m_radiators.push_back(p);
88
89 return;
90}
91
92
94 ATH_MSG_DEBUG ( "Constructor TRTTransitionRadiation::Initialize()" );
95 m_NumBins = 100; // Bins for TR generation (Nevski has 25 bins)
96 m_MinEnergyTR = 2.0*CLHEP::keV; // Minimum TR energy
97 m_MaxEnergyTR = 50.0*CLHEP::keV; // Maximum TR energy
98
99 m_GammaMin = 500.; // Min. value of Lorentz gamma for TR generation
100 m_EkinMin = (m_GammaMin-1.)*CLHEP::electron_mass_c2;
101
102 //We get the materials from GeoMaterial manager and convert them to G4Materials
103 //
104 //NB: Maybe we should control this from TRT_GeoModel instead?
105 //
106 // On a related note we should get rid of the xml files too -
107 // they are a loophole in the versioning.
108
109 G4Material* g4mat_Gas{nullptr};
110 G4Material* g4mat_FoilMaterial{nullptr};
111 std::string errorMessage;
112
113 // Get material information from storegate.
114 ISvcLocator *svcLocator = Gaudi::svcLocator(); // from Bootstrap
115 SmartIF<StoreGateSvc> detStore{svcLocator->service( "DetectorStore")};
116 if (!detStore) {
117 errorMessage = "Can not access Detector Store";
118 ATH_MSG_FATAL (errorMessage);
119 throw std::runtime_error(errorMessage);
120 }
121
122 StoredMaterialManager* materialManager = detStore->tryRetrieve<StoredMaterialManager>("MATERIALS");
123 if(materialManager) {
124 Geo2G4MaterialFactory geo2g4_material_fact;//Note - this is a very lightweight class!
125 g4mat_Gas = geo2g4_material_fact.Build(materialManager->getMaterial("trt::CO2"));
126 g4mat_FoilMaterial = geo2g4_material_fact.Build(materialManager->getMaterial("trt::CH2"));
127
128 //Note: One might wonder if we now own g4mat_Gas and
129 //g4mat_FoilMaterial, but since the elementfactory used inside
130 //Geo2G4MaterialFactory is static, we must NOT delete them!!
131 }
132 else {
133 ATH_MSG_INFO("GeoModel Material is not available. Construct TR materials using IRDBAccessSvc interface");
134 SmartIF<IGeoDbTagSvc> geoDbTagSvc{svcLocator->service("GeoDbTagSvc")};
135 if (!geoDbTagSvc) {
136 errorMessage = "Can not access GeoDbTagSvc";
137 ATH_MSG_FATAL (errorMessage);
138 throw std::runtime_error(errorMessage);
139
140 return;
141 };
142 SmartIF<IRDBAccessSvc> pAccessSvc{svcLocator->service(geoDbTagSvc->getParamSvcName())};
143 if ( !pAccessSvc) {
144 errorMessage = "Can not access " + geoDbTagSvc->getParamSvcName();
145 ATH_MSG_FATAL (errorMessage);
146 throw std::runtime_error(errorMessage);
147 }
148
149 std::map<std::string,int> materialComponentsMap; // Number of elements for each material
150 std::map<std::string,G4Material*> materialMap;
151 IRDBRecordset_ptr materialsRec = pAccessSvc->getRecordsetPtr("TRMaterials","","");
152 // Step #1. Count elements per material
153 for (const IRDBRecord_ptr& material : *materialsRec) {
154 std::string key = material->getString("NAME");
155 auto mapIt = materialComponentsMap.find(key);
156 if (mapIt == materialComponentsMap.end()) {
157 materialComponentsMap.emplace(key,1);
158 }
159 else {
160 materialComponentsMap.at(key) = mapIt->second + 1;
161 }
162 }
163 // Step #2. Build materials
164 for (const IRDBRecord_ptr& material :*materialsRec) {
165 std::string key = material->getString("NAME");
166 G4Material* g4Material{nullptr};
167 if(materialMap.find(key)==materialMap.end()) {
168 auto itNElements = materialComponentsMap.find(key);
169 if (itNElements==materialComponentsMap.end()) continue;
170 g4Material = new G4Material(key
171 ,material->getDouble("DENSITY")*(CLHEP::gram / CLHEP::cm3)
172 ,itNElements->second);
173 materialMap.emplace(key,g4Material);
174 }
175 else {
176 g4Material = materialMap.at(key);
177 }
178
179 G4Element* g4Element = G4Element::GetElement(material->getString("ELEMENT_NAME"));
180 if(!g4Element) {
181 errorMessage = "Wrong name for element " + material->getString("ELEMENT_NAME")
182 + " found in the description of material " + key;
183 ATH_MSG_FATAL (errorMessage);
184 throw std::runtime_error(errorMessage);
185 }
186 g4Material->AddElement(g4Element,(G4int)material->getInt("N"));
187 }
188 g4mat_Gas = materialMap.at("CO2");
189 g4mat_FoilMaterial = materialMap.at("CH2");
190 }
191
192 G4double PlasmaCof(4.0*M_PI*CLHEP::fine_structure_const*CLHEP::hbarc*CLHEP::hbarc*CLHEP::hbarc/CLHEP::electron_mass_c2);
193
194 // Plate and gas indicies
195 m_WplasmaFoil = sqrt( PlasmaCof * g4mat_FoilMaterial->GetElectronDensity() );
196 m_WplasmaGas = sqrt( PlasmaCof * g4mat_Gas->GetElectronDensity() );
197
198 ATH_MSG_DEBUG("Foil material : " << g4mat_FoilMaterial->GetName()
199 << "; plasma energy : " << m_WplasmaFoil/CLHEP::eV << " eV");
200 ATH_MSG_DEBUG("Gas : " << g4mat_Gas->GetName()
201 << "; plasma energy : " << m_WplasmaGas/CLHEP::eV << " eV");
202
203 m_Omg = new G4double[m_NumBins];
204 m_om = new G4double[m_NumBins];
205 m_Ey = new G4double[m_NumBins];
206 m_Sr = new G4double[m_NumBins];
207 m_sigmaGas = new G4double[m_NumBins];
208 m_sigmaFoil = new G4double[m_NumBins];
209
210 for(G4int j(0); j<m_NumBins; ++j) {
211 m_Omg[j] = log( m_MinEnergyTR ) +
213 m_om [j] = exp( m_Omg[j] );
214 m_sigmaGas[j] = ComputePhotoAbsCof( g4mat_Gas, m_om[j] );
215 m_sigmaFoil[j] = ComputePhotoAbsCof( g4mat_FoilMaterial, m_om[j] );
216 }
217
218 const std::string file=PathResolver::find_file(m_xmlfilename,"DATAPATH");
219 if (!file.empty())
220 {
221 m_XMLhandler->Process(file);
222 }
223 else ATH_MSG_WARNING("File "<<m_xmlfilename<<" not found! Expect problems.");
224
225 ATH_MSG_DEBUG ( "Constructor TRTTransitionRadiation::Initialize() DONE!" );
226 return;
227
228}
229
231//
232// Returns condition for application of the model depending on particle type
233
234
235G4bool TRTTransitionRadiation::IsApplicable(const G4ParticleDefinition& particle) {
236 //return true if PDG Charge and PDG Mass are non-zero.
237 return ( std::abs(particle.GetPDGCharge())>std::numeric_limits<double>::epsilon() && std::abs(particle.GetPDGMass())>std::numeric_limits<double>::epsilon() ) ;
238}
239
241//
242// GetContinuousStepLimit
243//
244
245G4double TRTTransitionRadiation::GetMeanFreePath(const G4Track& aTrack,
246 G4double,
247 G4ForceCondition* condition) {
248 const G4DynamicParticle *aParticle(aTrack.GetDynamicParticle());
249 if(std::abs(aParticle->GetDefinition()->GetPDGCharge())< std::numeric_limits<double>::epsilon() ||
250 std::abs(aParticle->GetDefinition()->GetPDGMass())<std::numeric_limits<double>::epsilon() ) {
251 *condition = NotForced;
252 } else {
253 *condition = Forced;
254 }
255 return DBL_MAX;
256}
257
259//
260//
261G4VParticleChange* TRTTransitionRadiation::PostStepDoIt( const G4Track& aTrack,
262 const G4Step& aStep ) {
263 aParticleChange.Initialize(aTrack);
264
265 // Obtain kinetic energy
266 const G4DynamicParticle *aParticle(aTrack.GetDynamicParticle());
267 G4double KineticEnergy(aParticle->GetKineticEnergy());
268
269 // Kinetic energy too low for TR generation for any kind of particle
270 if ( KineticEnergy < m_EkinMin )
271 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
272
273 // Particle properties ...
274 G4double ParticleMass(aParticle->GetDefinition()->GetPDGMass());
275 G4double Gamma(1.0 + KineticEnergy/ParticleMass);
276
277 if( Gamma < m_GammaMin ) {
278 ATH_MSG_VERBOSE ( "Leaving PostStepDoIt(): Gamma < " << m_GammaMin );
279 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
280 }
281
282 G4double FoilThickness(0.);
283 G4double GasThickness(0.);
284 BEflag BEflag(BARREL); //enum?
285
286 // Check whether particle is inside radiator
287 std::vector<TRTRadiatorParameters>::const_iterator currentRadiator(m_radiators.begin());
288 const std::vector<TRTRadiatorParameters>::const_iterator endOfRadiators(m_radiators.end());
289 while(currentRadiator!=endOfRadiators) {
290 if ( currentRadiator->GetLogicalVolume() ==
291 aTrack.GetVolume()->GetLogicalVolume() ) {
292 FoilThickness = currentRadiator->GetFoilThickness();
293 GasThickness = currentRadiator->GetGasThickness();
294 BEflag = currentRadiator->GetBEflag();
295 break;
296 }
297 ++currentRadiator;
298 }
299 // Particle not inside radiator - return
300 if ( currentRadiator == endOfRadiators ) {
301 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
302 }
303
304 ATH_MSG_VERBOSE ( "PostStepDoIt(): Inside radiator "
305 << currentRadiator->GetLogicalVolume()->GetName() );
306
307 G4ThreeVector ParticleDirection(aParticle->GetMomentumDirection());
308
309 if ( BEflag == ENDCAP ) {
310 G4double costh(std::abs(ParticleDirection[2]));
311 FoilThickness /= costh;
312 GasThickness /= costh;
313 }
314
315 G4double StepLength(aStep.GetStepLength());
316 G4int FoilsTraversed(static_cast<G4int>(StepLength/(FoilThickness+GasThickness)+0.5));
317
318 // Steplength so short, that no foils are crossed.
319 // Forget it...
320 // Actually the treatment we do is only valid for "many" foils.
321 // Should we put this requirement higher?
322 if ( FoilsTraversed == 0 ) {
323 ATH_MSG_VERBOSE ( "Leaving PostStepDoIt(): No foils traversed." );
324 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
325 }
326
327 G4double Gy(0.), Gmany(0.), Neff(0.);
328
329 // Generate, in the language of Artru:
330 // 1) G_y
331 // 2) N_{eff}
332 // 3) G_{many}
333 // Divide by photon energy to get photon spectrum (m_Ey)
334
335 for( G4int i(0); i<m_NumBins; ++i ) {
336
337#ifdef ARTRU // Formulas from Artru
338 Gy = XEmitanArtru( m_om[i],
339 Gamma,
340 GasThickness,
341 FoilThickness );
342 Neff = NeffArtru( m_sigmaGas[i],
343 m_sigmaFoil[i],
344 GasThickness,
345 FoilThickness,
346 FoilsTraversed );
347#else
348 Gy = XEmitanNevski( m_om[i], // Formulas from Cherry via
349 Gamma, // Pavel Nevski in G3
350 GasThickness,
351 FoilThickness );
352 Neff = NeffNevski( m_sigmaGas[i],
353 m_sigmaFoil[i],
354 GasThickness,
355 FoilThickness,
356 FoilsTraversed );
357#endif
358
359 Gmany = Neff*Gy; // Energy spectrum
360 m_Ey[i] = Gmany / m_om[i]; // Photon spectrum
361
362 }
363
364 // Integrate spectrum
365 G4double Sph(XInteg( m_Ey, m_Sr ));
366 G4double Sph_org(Sph);
367
368 if ( ( BEflag == BARREL ) && (Gamma > 0.0) ) {
369 G4double tempLog10Gamma(log10(Gamma)); //Hardcoded numbers from TestBeam study (E.Klinkby). Better to move elsewhere.
370 G4double pHT_dEdxData = -0.001 + 0.0005*tempLog10Gamma;
371 G4double pHT_TRData = 0.15/(1.0 + exp(-(tempLog10Gamma - 3.3)/0.27));
372 G4double pHT_dEdxMC = -0.0031 + 0.00016*tempLog10Gamma;
373 G4double pHT_TRMC = 0.1289/(1.0 + exp(-(tempLog10Gamma - 3.01)/0.1253));
374
375 G4double fudge = std::min(2.0,(pHT_dEdxData+pHT_TRData)/( pHT_dEdxMC+pHT_TRMC));
376
377 Sph *= fudge;
378 }
379
380 // Now how many photons are generated?
381 G4long NumberOfPhotonsGenerated(CLHEP::RandPoisson::shoot(Sph));
382
383 //reset Sph to original value.
384 Sph = Sph_org;
385
386 ATH_MSG_VERBOSE ( ">>>> FoilsTraversed = " << FoilsTraversed
387 << ", Gamma = " << Gamma
388 << ", Sph = " << Sph
389 << ", Nph = " << NumberOfPhotonsGenerated );
390 // No photons, return
391 if ( NumberOfPhotonsGenerated <= 0 ) {
392 ATH_MSG_VERBOSE ( "Leaving TRTTransitionRadiation::PostStepDoIt: "
393 << "No photons generated." );
394
395 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
396 }
397
398 // Now, generate the photons
399 aParticleChange.SetNumberOfSecondaries( NumberOfPhotonsGenerated );
400 G4StepPoint *pPostStepPoint(aStep.GetPostStepPoint());
401 G4ThreeVector x0(pPostStepPoint->GetPosition());
402
403 G4double Esum(0.), XRannu(0.), Ephot(0.), Theta(0.), Phi(0.),
404 X(0.), Y(0.), Z(0.);
405
406 G4ThreeVector pos(0),TRPhotonDirection(0);
407
408 for( G4int i(0); i<NumberOfPhotonsGenerated; ++i ) {
409 XRannu = CLHEP::RandFlat::shoot();
410 Ephot = exp( XFinter( XRannu*Sph, m_Sr, m_Omg ) );
411
412 Esum += Ephot;
413
414 pos = pPostStepPoint->GetPosition();
415
416 ATH_MSG_VERBOSE ( " Photon generated with energy : " << Ephot/CLHEP::keV
417 << " keV; position: ( " << pos.x()/CLHEP::cm << " cm, "
418 << pos.y()/CLHEP::cm << " cm, " << pos.z()/CLHEP::cm << " cm )" );
419
420 // Angle w.r.t. electron direction (this is not correct but anyway...)
421 Theta = std::abs( CLHEP::RandGaussZiggurat::shoot( 0.0, M_PI/Gamma ) );
422 if( Theta >= 0.1 ) Theta = 0.1;
423
424 Phi = (2.*M_PI)*CLHEP::RandFlat::shoot();
425
426 X = sin(Theta)*cos(Phi);
427 Y = sin(Theta)*sin(Phi);
428 Z = cos(Theta);
429
430 TRPhotonDirection.set(X,Y,Z);
431 TRPhotonDirection.rotateUz( ParticleDirection );
432 TRPhotonDirection.unit();
433
434 G4DynamicParticle*
435 TRPhoton = new G4DynamicParticle( G4Gamma::Gamma(),
436 TRPhotonDirection,
437 Ephot );
438 aParticleChange.AddSecondary( TRPhoton );
439
440 }
441
442 aParticleChange.ProposeEnergy( KineticEnergy-Esum );
443
444 ATH_MSG_VERBOSE ( "Leaving TRTTransitionRadiation::PostStepDoIt; "
445 << "Number of TR photons generated: " << NumberOfPhotonsGenerated );
446
447 return &aParticleChange;
448}
449
450inline G4double TRTTransitionRadiation::NeffNevski( const G4double & sigGas,
451 const G4double & sigFoil,
452 const G4double & GasThickness,
453 const G4double & FoilThickness,
454 const G4int & FoilsTraversed ) const {
455 if ( FoilsTraversed <= 0 )
456 return 0.;
457
458 G4double C1 = sigGas * GasThickness*FoilsTraversed;
459 G4double C2 = sigFoil*FoilThickness*FoilsTraversed;
460
461 G4double W1 = exp(-C1);
462 G4double W2 = exp(-C2);
463
464 return FoilsTraversed * (1-W1*W2)/(C1+C2);
465}
466
467inline G4double TRTTransitionRadiation::NeffArtru( const G4double & sigGas,
468 const G4double & sigFoil,
469 const G4double & GasThickness,
470 const G4double & FoilThickness,
471 const G4int & FoilsTraversed ) const {
472 if ( FoilsTraversed <= 0 )
473 return 0.;
474
475 G4double sigma = FoilThickness*sigFoil + GasThickness*sigGas;
476 G4double xNF = FoilsTraversed;
477 return (1.-exp(-xNF*sigma))/(1.-exp(-sigma));
478}
479
481//
482// Nevskis XEmitan-function
483//
484//inline G4double max(G4double a, G4double b) {return a>b?a:b;}
485
486inline G4double TRTTransitionRadiation::XEmitanNevski( const G4double & PhotonEnergy,
487 const G4double & Gamma,
488 const G4double & Al1,
489 const G4double & Al2 ) const {
490 //FIXME hardcoded numbers are not good.
491 G4double eps = 0.05;
492
493 G4double Al = Al1 + Al2;
494
495 G4double CK = 1/(2*M_PI*CLHEP::hbarc);
496
497 // Plasma frequencies:
498 G4double om1 = m_WplasmaGas;
499 G4double om2 = m_WplasmaFoil;
500
501 G4double om12 = (om1*om1-om2*om2)/PhotonEnergy;
502 G4double oms = (Al1*om1*om1+Al2*om2*om2)/Al;
503 G4double Ao = 0.5*Al1*om12*CK;
504 G4double Bo = 0.5*Al2*om12*CK;
505 G4double ro = sqrt(oms)/Gamma;
506 G4double co = 0.5*(PhotonEnergy/(Gamma*Gamma)+oms/PhotonEnergy);
507 G4double cc = co*Al*CK;
508 G4double cb = M_PI*Al2/Al;
509 G4int i1 = static_cast<G4int>(1.0+(ro>co?ro:co)*Al*CK);//(G4int)(1.0+max(ro,co)*Al*CK);
510 G4int i2 = static_cast<G4int>(Gamma*ro*Al*CK);
511
512 G4double Sum = 0.0;
513 G4double Smx = 0.0;
514 G4double Sm1 = 0.0;
515 G4double S = 0.0;
516 G4double So = 0.0;
517
518 G4double Sx = 0.0;
519
520 for(G4int i(i1); i<=i2; ++i) {
521 Sx = So;
522 So = S;
523 S = (static_cast<G4double>(i)-cc)*(sin(cb*(static_cast<G4double>(i)-Ao))/((static_cast<G4double>(i)-Ao)*(static_cast<G4double>(i)+Bo)))*
524 (sin(cb*(static_cast<G4double>(i)-Ao))/((static_cast<G4double>(i)-Ao)*(static_cast<G4double>(i)+Bo)));
525 Smx=Smx>S?Smx:S;//max(Smx,S);
526 Sm1=Sm1>S?Sm1:S;//max(Sm1,S);
527 Sum += S;
528
529 if( (So<Sx)&&(So<=S) ) {
530 if( Sm1<eps*Smx ) break;
531 Sm1 = 0.0;
532 }
533 }
534
535 return ( CLHEP::fine_structure_const*om12*om12*Al*Al*CK*CK*Sum/M_PI );
536}
537
539//
540// Same as above but for the Artru way of doing things
541//
542inline G4double TRTTransitionRadiation::XEmitanArtru( const G4double & omega,
543 const G4double & Gamma,
544 const G4double & d2,
545 const G4double & d1 ) const {
546 G4double phimax = 4.*M_PI;
547 G4double ginv2 = 1./(Gamma*Gamma);
548 G4double tau = d2/d1;
549
550 G4double mom = omega/CLHEP::hbarc;
551
552 G4double xi1 = m_WplasmaFoil/omega;
553 G4double xi2 = m_WplasmaGas/omega;
554
555 G4double xi1sq = xi1*xi1;
556 G4double xi2sq = xi2*xi2;
557
558 G4double Z1 = 2./(mom*(ginv2+xi1sq));
559 G4double Z2 = 2./(mom*(ginv2+xi2sq));
560
561 G4double a = d1/Z2;
562 G4double V = d1*(1./Z1-1./Z2);
563
564 G4double p0 = mom*d1*(xi1sq-xi2sq)/(4.*M_PI);
565 G4double pmin = p0 + a*(1.+tau)/(2.*M_PI);
566 G4double ymax = tau*phimax;
567 G4double pmax = p0 + ymax*(1.+tau)/(2.*M_PI);
568
569 G4int ipmin = (int)pmin+1;
570 G4int ipmax = (int)pmax;
571
572 G4double sum=0.;
573 G4double y0 = -0.5*mom*d1*(xi1sq-xi2sq)/(1+tau);
574 G4double dy = (2.*M_PI)/(1+tau);
575 G4double y;
576
577 G4double term = 0.;
578 G4double t1 = 0.;
579 G4double t2 = 0.;
580 G4double tx = 0.;
581
582 for ( G4int ip(ipmin); ip<=ipmax; ++ip ) {
583 t2 = t1;
584 t1 = term;
585
586 y = y0 + ip*dy;
587 term = (1./y-1./(y+V))*sin(0.5*(y+V));
588 term *= term;
589 term *= y-a;
590 sum += term;
591
592 //G4double theta2 = 2.*y*tau/(mom*d2) - ginv2 - xi2sq;
593 //G4double theta = sqrt(theta2);
594
595 tx = (tx>term) ? tx : term;
596
597 // Locate local maxima (previous term)
598
599 if ( t2 > 0. && t1 > term && t1 > t2 ) {
600 if ( t1 < 0.0001*sum ) break;
601 }
602 }
603
604 return ( sum*4.*2.*CLHEP::fine_structure_const/(1.+tau));
605}
606
608//
609// Nevskis XInteg function
610// Spectrum integration
611//
612
613inline G4double TRTTransitionRadiation::XInteg( const G4double * yy, G4double * ss ) const {
614 G4double S(0.);
615 ss[0] = 0.;
616
617 for( G4int i(1); i<m_NumBins; ++i ) {
618 S += 0.5*( m_Omg[i]-m_Omg[i-1] )*( yy[i-1]*m_om[i-1] + yy[i]*m_om[i] );
619 ss[i] = S;
620 }
621
622 return S;
623}
624
626//
627// Nevskis XFINTER function
628
629G4double TRTTransitionRadiation::XFinter( const G4double & X, const G4double* A, const G4double* F ) const {
630 G4int K1(0);
631 G4int K2(0);
632
633 if( (K1>=m_NumBins) || ( X < A[K1] ) ) K1 = 0; //FIXME This line has no effect?
634 if( (K2>=m_NumBins) || ( X > A[K2] ) ) K2 = m_NumBins-1;
635
636 G4int K(0);
637 while( K2-K1>1 ) {
638 K = static_cast<G4int>((K1+K2)/2);
639 if( A[K]<X ) {
640 K1 = K;
641 } else {
642 K2 = K;
643 }
644 }
645
646 G4double X1(A[K1]);
647 G4double X2(A[K1+1]);
648
649 G4double xfinter((F[K1]*(X-X2)+F[K1+1]*(X1-X))/(X1-X2));
650 return xfinter;
651}
652
654//
655// Computes Sandia photo absorption cross section coefficients for Material
656
658 const G4double & PhotonEnergy ) const {
659
660 const G4double *SandiaCof(Material->GetSandiaTable()->GetSandiaCofForMaterial(PhotonEnergy));
661
662 const G4double Energy2(PhotonEnergy*PhotonEnergy);
663 const G4double Energy3(PhotonEnergy*Energy2);
664 const G4double Energy4(Energy2*Energy2);
665
666 return ( SandiaCof[0]/PhotonEnergy +
667 SandiaCof[1]/Energy2 +
668 SandiaCof[2]/Energy3 +
669 SandiaCof[3]/Energy4 );
670}
#define M_PI
#define ATH_MSG_FATAL(x)
#define ATH_MSG_INFO(x)
#define ATH_MSG_VERBOSE(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
Definition of the abstract IRDBAccessSvc interface.
std::shared_ptr< IRDBRecordset > IRDBRecordset_ptr
Definition of the abstract IRDBRecord interface.
Definition of the abstract IRDBRecordset interface.
std::unique_ptr< IRDBRecord > IRDBRecord_ptr
static Double_t a
static Double_t ss
#define F(x, y, z)
Definition MD5.cxx:112
@ Material
@ Phi
Definition RPCdef.h:8
#define y
AthMessaging(IMessageSvc *msgSvc, const std::string &name)
Constructor.
G4Material * Build(const GeoMaterial *)
static std::string find_file(const std::string &logical_file_name, const std::string &search_path)
This class holds one or more material managers and makes them storeable, under StoreGate.
virtual const GeoMaterial * getMaterial(std::string_view name)=0
G4bool IsApplicable(const G4ParticleDefinition &)
G4double NeffNevski(const G4double &sigGas, const G4double &sigFoil, const G4double &GasThickness, const G4double &FoilThickness, const G4int &FoilsTraversed) const
void AddRadiatorParameters(const TRTRadiatorParameters &p)
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep)
G4double XFinter(const G4double &X, const G4double *A, const G4double *F) const
G4double XInteg(const G4double *yy, G4double *ss) const
std::vector< TRTRadiatorParameters > m_radiators
G4double XEmitanNevski(const G4double &PhotonEnergy, const G4double &Gamma, const G4double &GasThickness, const G4double &FoilThickness) const
G4double NeffArtru(const G4double &sigGas, const G4double &sigFoil, const G4double &GasThickness, const G4double &FoilThickness, const G4int &FoilsTraversed) const
G4double XEmitanArtru(const G4double &PhotonEnergy, const G4double &Gamma, const G4double &GasThickness, const G4double &FoilThickness) const
TRTTransitionRadiation(const G4String &processName="TransitionRadiation", const std::string &xmlfilename="TRgeomodelgeometry.xml")
TRRegionXMLHandler * m_XMLhandler
G4double GetMeanFreePath(const G4Track &aTrack, G4double, G4ForceCondition *)
G4double ComputePhotoAbsCof(const G4Material *Material, const G4double &GammaEnergy) const
double ymax
Definition listroot.cxx:64
hold the test vectors and ease the comparison
TFile * file