1 #define AVO 6.022e23 //Avogadro's number (#/mol) 2 #define EMASS 9.109e-31*CLHEP::kg 3 #define MillerDriftSpeed true 5 #define GASGAP 0.25*CLHEP::cm //S2 generation region 6 #define BORDER 0*CLHEP::cm //liquid-gas border z-coordinate 8 #define QE_EFF 1 //a base or maximum quantum efficiency 9 #define phe_per_e 1 //S2 gain for quick studies 12 #define WIN 0*CLHEP::mm //top Cu block (also, quartz window) 13 #define TOP 0 //top grid wires 14 #define ANE 0 //anode mesh 15 #define SRF 0 //liquid-gas interface 16 #define GAT 0 //gate grid 17 #define CTH 0 //cathode grid 18 #define BOT 0 //bottom PMT grid 19 #define PMT 0 //bottom Cu block and PMTs 20 #define MIN_ENE -1*CLHEP::eV //lets you turn NEST off BELOW a certain energy 21 #define MAX_ENE 1.*CLHEP::TeV //lets you turn NEST off ABOVE a certain energy 22 #define HIENLIM 5*CLHEP::MeV //energy at which Doke model used exclusively 23 #define R_TOL 0.2*CLHEP::mm //tolerance (for edge events) 24 #define R_MAX 1*CLHEP::km //for corraling diffusing electrons 25 #define Density_LXe 2.9 //reference density for density-dep. effects 26 #define Density_LAr 1.393 27 #define Density_LNe 1.207 28 #define Density_LKr 2.413 30 #include "Geant4/Randomize.hh" 31 #include "Geant4/G4Ions.hh" 32 #include "Geant4/G4OpticalPhoton.hh" 33 #include "Geant4/G4VProcess.hh" 38 #include "CLHEP/Random/RandGauss.h" 39 #include "CLHEP/Random/RandFlat.h" 85 CLHEP::RandGauss GaussGen(
fEngine);
86 CLHEP::RandFlat UniformGen(
fEngine);
100 bool fTrackSecondariesFirst =
false;
101 bool fExcitedNucleus =
false;
102 bool fVeryHighEnergy =
false;
104 bool fMultipleScattering =
false;
107 if( aTrack.GetParentID() == 0 && aTrack.GetCurrentStepNumber() == 1 ) {
108 fExcitedNucleus =
false;
109 fVeryHighEnergy =
false;
111 fMultipleScattering =
false;
114 const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
115 G4ParticleDefinition *pDef = aParticle->GetDefinition();
116 G4String particleName = pDef->GetParticleName();
117 const G4Material* aMaterial = aStep.GetPreStepPoint()->GetMaterial();
118 const G4Material* bMaterial = aStep.GetPostStepPoint()->GetMaterial();
120 if((particleName ==
"neutron" || particleName ==
"antineutron") &&
121 aStep.GetTotalEnergyDeposit() <= 0)
129 G4Element *ElementA = NULL, *ElementB = NULL;
131 const G4ElementVector* theElementVector1 =
132 aMaterial->GetElementVector();
133 ElementA = (*theElementVector1)[0];
136 const G4ElementVector* theElementVector2 =
137 bMaterial->GetElementVector();
138 ElementB = (*theElementVector2)[0];
140 G4int z1,z2,j=1; G4bool NobleNow=
false,NobleLater=
false;
141 if (ElementA) z1 = (G4int)(ElementA->GetZ());
else z1 = -1;
142 if (ElementB) z2 = (G4int)(ElementB->GetZ());
else z2 = -1;
143 if ( z1==2 || z1==10 || z1==18 || z1==36 || z1==54 ) {
153 if ( z2==2 || z2==10 || z2==18 || z2==36 || z2==54 ) {
164 if ( !NobleNow && !NobleLater )
169 G4StepPoint* pPreStepPoint = aStep.GetPreStepPoint();
170 G4StepPoint* pPostStepPoint = aStep.GetPostStepPoint();
171 G4ThreeVector
x1 = pPostStepPoint->GetPosition();
172 G4ThreeVector x0 = pPreStepPoint->GetPosition();
173 G4double evtStrt = pPreStepPoint->GetGlobalTime();
174 G4double
t0 = pPreStepPoint->GetLocalTime();
175 G4double
t1 = pPostStepPoint->GetLocalTime();
181 G4bool outside =
false, inside =
false, InsAndOuts =
false;
182 G4MaterialPropertiesTable* aMaterialPropertiesTable =
183 aMaterial->GetMaterialPropertiesTable();
184 if ( NobleNow && !NobleLater ) outside =
true;
185 if ( !NobleNow && NobleLater ) {
186 aMaterial = bMaterial; inside =
true; z1 = z2;
187 aMaterialPropertiesTable = bMaterial->GetMaterialPropertiesTable();
189 if ( NobleNow && NobleLater &&
190 aMaterial->GetDensity() != bMaterial->GetDensity() )
194 G4double Density = aMaterial->GetDensity()/(CLHEP::g/CLHEP::cm3);
195 G4double nDensity = Density*
AVO;
196 G4int Phase = aMaterial->GetState();
197 G4double ElectricField(0.), FieldSign(0.);
198 G4bool GlobalFields =
false;
201 ElectricField = aMaterialPropertiesTable->GetConstProperty(
"ELECTRICFIELD");
205 if ( x1[2] <
WIN && x1[2] >
TOP && Phase == kStateGas )
206 ElectricField = aMaterialPropertiesTable->GetConstProperty(
"ELECTRICFIELDWINDOW");
207 else if ( x1[2] <
TOP && x1[2] >
ANE && Phase == kStateGas )
208 ElectricField = aMaterialPropertiesTable->GetConstProperty(
"ELECTRICFIELDTOP");
209 else if ( x1[2] <
ANE && x1[2] >
SRF && Phase == kStateGas )
210 ElectricField = aMaterialPropertiesTable->
211 GetConstProperty(
"ELECTRICFIELDANODE");
212 else if ( x1[2] <
SRF && x1[2] >
GAT && Phase == kStateLiquid )
213 ElectricField = aMaterialPropertiesTable->
214 GetConstProperty(
"ELECTRICFIELDSURFACE");
215 else if ( x1[2] <
GAT && x1[2] >
CTH && Phase == kStateLiquid )
216 ElectricField = aMaterialPropertiesTable->
217 GetConstProperty(
"ELECTRICFIELDGATE");
218 else if ( x1[2] <
CTH && x1[2] >
BOT && Phase == kStateLiquid )
219 ElectricField = aMaterialPropertiesTable->
220 GetConstProperty(
"ELECTRICFIELDCATHODE");
221 else if ( x1[2] <
BOT && x1[2] >
PMT && Phase == kStateLiquid )
222 ElectricField = aMaterialPropertiesTable->
223 GetConstProperty(
"ELECTRICFIELDBOTTOM");
225 ElectricField = aMaterialPropertiesTable->
226 GetConstProperty(
"ELECTRICFIELD");
228 if ( ElectricField >= 0 ) FieldSign = 1;
else FieldSign = -1;
229 ElectricField = std::abs((1e3*ElectricField)/(CLHEP::kilovolt/CLHEP::cm));
230 G4double Temperature = aMaterial->GetTemperature();
231 G4double ScintillationYield, ResolutionScale, R0 = 1.0*CLHEP::um,
232 DokeBirks[3], ThomasImel = 0.00, delta = 1*CLHEP::mm;
233 DokeBirks[0] = 0.00; DokeBirks[2] = 1.00;
234 G4double PhotMean = 7*CLHEP::eV, PhotWidth = 1.0*CLHEP::eV;
235 G4double SingTripRatioR, SingTripRatioX, tau1, tau3, tauR = 0*CLHEP::ns;
238 ScintillationYield = 1 / (41.3*CLHEP::eV);
240 ResolutionScale = 0.2;
241 PhotMean = 15.9*CLHEP::eV;
242 tau1 = GaussGen.fire(10.0*CLHEP::ns,0.0*CLHEP::ns);
243 tau3 = 1.6e3*CLHEP::ns;
247 ScintillationYield = 1 / (29.2*CLHEP::eV);
249 ResolutionScale = 0.13;
250 PhotMean = 15.5*CLHEP::eV; PhotWidth = 0.26*CLHEP::eV;
251 tau1 = GaussGen.fire(10.0*CLHEP::ns,10.*CLHEP::ns);
252 tau3 = GaussGen.fire(15.4e3*CLHEP::ns,200*CLHEP::ns);
255 ScintillationYield = 1 / (19.5*CLHEP::eV);
257 ResolutionScale = 0.107;
258 R0 = 1.568*CLHEP::um;
260 ThomasImel = 0.156977*pow(ElectricField,-0.1);
261 DokeBirks[0] = 0.07*pow((ElectricField/1.0e3),-0.85);
266 DokeBirks[0] = 0.0003;
268 } nDensity /= 39.948;
269 PhotMean = 9.69*CLHEP::eV; PhotWidth = 0.22*CLHEP::eV;
270 tau1 = GaussGen.fire(6.5*CLHEP::ns,0.8*CLHEP::ns);
271 tau3 = GaussGen.fire(1300*CLHEP::ns,50*CLHEP::ns);
272 tauR = GaussGen.fire(0.8*CLHEP::ns,0.2*CLHEP::ns);
275 if ( Phase == kStateGas ) ScintillationYield = 1 / (30.0*CLHEP::eV);
276 else ScintillationYield = 1 / (15.0*CLHEP::eV);
278 ResolutionScale = 0.05;
279 PhotMean = 8.43*CLHEP::eV;
280 tau1 = GaussGen.fire(2.8*CLHEP::ns,.04*CLHEP::ns);
281 tau3 = GaussGen.fire(93.*CLHEP::ns,1.1*CLHEP::ns);
282 tauR = GaussGen.fire(12.*CLHEP::ns,.76*CLHEP::ns);
285 default: nDensity /= 131.293;
286 ScintillationYield = 48.814+0.80877*Density+2.6721*pow(Density,2.);
287 ScintillationYield /= CLHEP::keV;
291 ResolutionScale = 1.00 *
292 (0.12724-0.032152*Density-0.0013492*pow(Density,2.));
294 if ( Phase == kStateLiquid ) {
295 ResolutionScale *= 1.5;
299 R0 = 69.492*pow(ElectricField,-0.50422)*CLHEP::um;
301 DokeBirks[0]= 19.171*pow(ElectricField+25.552,-0.83057)+0.026772;
311 delta = 0.4*CLHEP::mm;
312 PhotMean = 6.97*CLHEP::eV; PhotWidth = 0.23*CLHEP::eV;
317 tau1 = GaussGen.fire(3.1*CLHEP::ns,.7*CLHEP::ns);
318 tau3 = GaussGen.fire(24.*CLHEP::ns,1.*CLHEP::ns);
320 else if ( Phase == kStateGas ) {
323 ScintillationYield = 1 / (12.98*CLHEP::eV); }
325 G4double Townsend = (ElectricField/nDensity)*1e17;
326 DokeBirks[0] = 0.0000;
327 DokeBirks[2] = 0.1933*pow(Density,2.6199)+0.29754 -
328 (0.045439*pow(Density,2.4689)+0.066034)*log10(ElectricField);
329 if ( ElectricField>6990 ) DokeBirks[2]=0.0;
330 if ( ElectricField<1000 ) DokeBirks[2]=0.2;
331 if ( ElectricField<100. ) { DokeBirks[0]=0.18; DokeBirks[2]=0.58; }
332 if( Density < 0.061 ) ThomasImel = 0.041973*pow(Density,1.8105);
333 else if( Density >= 0.061 && Density <= 0.167 )
334 ThomasImel=5.9583e-5+0.0048523*Density-0.023109*pow(Density,2.);
335 else ThomasImel = 6.2552e-6*pow(Density,-1.9963);
336 if(ElectricField)ThomasImel=1.2733e-5*pow(Townsend/Density,-0.68426);
338 PhotMean = 7.1*CLHEP::eV; PhotWidth = 0.2*CLHEP::eV;
339 tau1 = GaussGen.fire(5.18*CLHEP::ns,1.55*CLHEP::ns);
340 tau3 = GaussGen.fire(100.1*CLHEP::ns,7.9*CLHEP::ns);
343 tau1 = 3.5*CLHEP::ns; tau3 = 20.*CLHEP::ns; tauR = 40.*CLHEP::ns;
348 G4double anExcitationEnergy = ((
const G4Ions*)(pDef))->
349 GetExcitationEnergy();
350 G4double TotalEnergyDeposit =
351 aMaterialPropertiesTable->GetConstProperty(
"ENERGY_DEPOSIT_TOT" );
352 G4bool convert =
false, annihil =
false;
354 if(pPreStepPoint->GetKineticEnergy()>=(2*CLHEP::electron_mass_c2) &&
355 !pPostStepPoint->GetKineticEnergy() &&
356 !aStep.GetTotalEnergyDeposit() && aParticle->GetPDGcode()==22) {
357 convert =
true; TotalEnergyDeposit = CLHEP::electron_mass_c2;
359 if(pPreStepPoint->GetKineticEnergy() &&
360 !pPostStepPoint->GetKineticEnergy() &&
361 aParticle->GetPDGcode()==-11) {
362 annihil =
true; TotalEnergyDeposit += aStep.GetTotalEnergyDeposit();
364 G4bool either =
false;
365 if(inside || outside || convert || annihil || InsAndOuts) either=
true;
367 if( anExcitationEnergy<100*CLHEP::eV && aStep.GetTotalEnergyDeposit()<1*CLHEP::eV &&
368 !either && !fExcitedNucleus )
371 if ( !annihil ) TotalEnergyDeposit += aStep.GetTotalEnergyDeposit();
372 if ( !convert ) aMaterialPropertiesTable->
373 AddConstProperty(
"ENERGY_DEPOSIT_TOT", TotalEnergyDeposit );
375 TotalEnergyDeposit = aStep.GetTotalEnergyDeposit();
380 G4double InitialKinetEnergy = aMaterialPropertiesTable->
381 GetConstProperty(
"ENERGY_DEPOSIT_GOL" );
383 if ( InitialKinetEnergy == 0 ) {
384 G4double tE = pPreStepPoint->GetKineticEnergy()+anExcitationEnergy;
385 if ( (fabs(tE-1.8*CLHEP::keV) < CLHEP::eV || fabs(tE-9.4*CLHEP::keV) < CLHEP::eV) &&
386 Phase == kStateLiquid && z1 == 54 ) tE = 9.4*CLHEP::keV;
387 if ( fKr83m && ElectricField != 0 )
389 aMaterialPropertiesTable->
390 AddConstProperty (
"ENERGY_DEPOSIT_GOL", tE );
394 if ( anExcitationEnergy ) fExcitedNucleus =
true;
398 if(outside){ aMaterialPropertiesTable->
399 AddConstProperty(
"ENERGY_DEPOSIT_GOL",
400 InitialKinetEnergy-pPostStepPoint->GetKineticEnergy());
401 if(aMaterialPropertiesTable->
402 GetConstProperty(
"ENERGY_DEPOSIT_GOL")<0)
403 aMaterialPropertiesTable->AddConstProperty(
"ENERGY_DEPOSIT_GOL",0);
407 if(inside) { aMaterialPropertiesTable->
408 AddConstProperty(
"ENERGY_DEPOSIT_GOL",
409 InitialKinetEnergy+pPreStepPoint->GetKineticEnergy());
410 if ( TotalEnergyDeposit > 0 && InitialKinetEnergy == 0 ) {
411 aMaterialPropertiesTable->AddConstProperty(
"ENERGY_DEPOSIT_GOL",0);
412 TotalEnergyDeposit = .000000;
418 aMaterialPropertiesTable->
419 AddConstProperty(
"ENERGY_DEPOSIT_GOL",(-0.1*CLHEP::keV)+
420 InitialKinetEnergy-pPostStepPoint->GetKineticEnergy());
421 InitialKinetEnergy = bMaterial->GetMaterialPropertiesTable()->
422 GetConstProperty(
"ENERGY_DEPOSIT_GOL");
423 bMaterial->GetMaterialPropertiesTable()->
424 AddConstProperty(
"ENERGY_DEPOSIT_GOL",(-0.1*CLHEP::keV)+
425 InitialKinetEnergy+pPreStepPoint->GetKineticEnergy());
426 if(aMaterialPropertiesTable->
427 GetConstProperty(
"ENERGY_DEPOSIT_GOL")<0)
428 aMaterialPropertiesTable->AddConstProperty(
"ENERGY_DEPOSIT_GOL",0);
429 if ( bMaterial->GetMaterialPropertiesTable()->
430 GetConstProperty(
"ENERGY_DEPOSIT_GOL") < 0 )
431 bMaterial->GetMaterialPropertiesTable()->
432 AddConstProperty (
"ENERGY_DEPOSIT_GOL", 0 );
434 InitialKinetEnergy = aMaterialPropertiesTable->
435 GetConstProperty(
"ENERGY_DEPOSIT_GOL");
437 InitialKinetEnergy += 2*CLHEP::electron_mass_c2;
441 InitialKinetEnergy -= 2*CLHEP::electron_mass_c2;
443 aMaterialPropertiesTable->
444 AddConstProperty(
"ENERGY_DEPOSIT_GOL",InitialKinetEnergy);
445 if (anExcitationEnergy < 1
e-100 && aStep.GetTotalEnergyDeposit()==0 &&
446 aMaterialPropertiesTable->GetConstProperty(
"ENERGY_DEPOSIT_GOL")==0 &&
447 aMaterialPropertiesTable->GetConstProperty(
"ENERGY_DEPOSIT_TOT")==0)
451 if ( aTrack.GetCreatorProcess() )
452 procName = aTrack.GetCreatorProcess()->GetProcessName();
456 fMultipleScattering =
true;
459 if ( fAlpha ) delta = 1000.*CLHEP::km;
460 G4int i, k, counter = 0; G4double pos[3];
463 fMultipleScattering =
true;
467 char xCoord[80];
char yCoord[80];
char zCoord[80];
468 G4bool exists =
false;
471 sprintf(xCoord,
"POS_X_%d",i); sprintf(yCoord,
"POS_Y_%d",i);
472 sprintf(zCoord,
"POS_Z_%d",i);
473 pos[0] = aMaterialPropertiesTable->GetConstProperty(xCoord);
474 pos[1] = aMaterialPropertiesTable->GetConstProperty(yCoord);
475 pos[2] = aMaterialPropertiesTable->GetConstProperty(zCoord);
476 if ( sqrt(pow(x1[0]-pos[0],2.)+pow(x1[1]-pos[1],2.)+
477 pow(x1[2]-pos[2],2.)) < delta ) {
478 exists =
true;
break;
481 if(!exists && TotalEnergyDeposit) {
483 sprintf(xCoord,
"POS_X_%i",j); sprintf(yCoord,
"POS_Y_%i",j);
484 sprintf(zCoord,
"POS_Z_%i",j);
486 aMaterialPropertiesTable->AddConstProperty( xCoord, x1[0] );
487 aMaterialPropertiesTable->AddConstProperty( yCoord, x1[1] );
488 aMaterialPropertiesTable->AddConstProperty( zCoord, x1[2] );
490 aMaterialPropertiesTable->
491 AddConstProperty(
"TOTALNUM_INT_SITES", j );
499 G4double
a1 = ElementA->GetA();
500 z2 = pDef->GetAtomicNumber();
501 G4double
a2 = (G4double)(pDef->GetAtomicMass());
502 if ( particleName ==
"alpha" || (z2 == 2 && a2 == 4) )
504 if ( fAlpha || abs(aParticle->GetPDGcode()) == 2112 )
506 G4double epsilon = 11.5*(TotalEnergyDeposit/CLHEP::keV)*pow(z1,(-7./3.));
507 G4double gamma = 3.*pow(epsilon,0.15)+0.7*pow(epsilon,0.6)+epsilon;
508 G4double kappa = 0.133*pow(z1,(2./3.))*pow(a2,(-1./2.))*(2./3.);
510 if ( (z1 == z2 && pDef->GetParticleType() ==
"nucleus") ||
511 particleName ==
"neutron" || particleName ==
"antineutron" ) {
513 if ( z1 == 18 && Phase == kStateLiquid )
518 if ( ElectricField == 0 && Phase == kStateLiquid ) {
519 if ( z1 == 54 ) ThomasImel = 0.19;
520 if ( z1 == 18 ) ThomasImel = 0.25;
523 0.3065*exp(-0.008806*pow(ElectricField,0.76313));
528 G4double MeanNumberOfQuanta =
529 ScintillationYield*TotalEnergyDeposit;
532 G4double sigma = sqrt(ResolutionScale*MeanNumberOfQuanta);
534 G4int(floor(GaussGen.fire(MeanNumberOfQuanta,sigma)+0.5));
536 if (LeffVar > 1) { LeffVar = 1.00000; }
537 if (LeffVar < 0) { LeffVar = 0; }
541 if(TotalEnergyDeposit < 1/ScintillationYield || NumQuanta < 0)
547 G4int NumIons = NumQuanta - NumExcitons;
553 G4double dE, dx=0, LET=0, recombProb;
554 dE = TotalEnergyDeposit/CLHEP::MeV;
555 if ( particleName !=
"e-" && particleName !=
"e+" && z1 != z2 &&
556 particleName !=
"mu-" && particleName !=
"mu+" ) {
562 if(LET) dx = dE/(Density*LET);
563 if(abs(aParticle->GetPDGcode())==2112) dx=0;
566 dx = aStep.GetStepLength()/CLHEP::cm;
567 if(dx) LET = (dE/dx)*(1/Density);
568 if ( LET > 0 && dE > 0 && dx > 0 ) {
571 particleName ==
"e-" ) {
572 dx /= ratio; LET *= ratio; }}
574 DokeBirks[1] = DokeBirks[0]/(1-DokeBirks[2]);
576 recombProb = (DokeBirks[0]*LET)/(1+DokeBirks[1]*LET)+DokeBirks[2];
577 if ( Phase == kStateLiquid ) {
578 if ( z1 == 54 ) recombProb *= (Density/
Density_LXe);
579 if ( z1 == 18 ) recombProb *= (Density/
Density_LAr);
582 if(recombProb<0) recombProb=0;
583 if(recombProb>1) recombProb=1;
587 G4int NumPhotons = NumExcitons +
BinomFluct(NumIons,recombProb);
588 G4int NumElectrons = NumQuanta - NumPhotons;
599 char numExc[80];
char numIon[80];
char numPho[80];
char numEle[80];
600 sprintf(numExc,
"N_EXC_%i",counter); sprintf(numIon,
"N_ION_%i",counter);
601 aMaterialPropertiesTable->AddConstProperty( numExc, NumExcitons );
602 aMaterialPropertiesTable->AddConstProperty( numIon, NumIons );
603 sprintf(numPho,
"N_PHO_%i",counter); sprintf(numEle,
"N_ELE_%i",counter);
604 aMaterialPropertiesTable->AddConstProperty( numPho, NumPhotons );
605 aMaterialPropertiesTable->AddConstProperty( numEle, NumElectrons );
609 char trackL[80];
char time00[80];
char time01[80];
char energy[80];
610 sprintf(trackL,
"TRACK_%i",counter); sprintf(energy,
"ENRGY_%i",counter);
611 sprintf(time00,
"TIME0_%i",counter); sprintf(time01,
"TIME1_%i",counter);
612 delta = aMaterialPropertiesTable->GetConstProperty( trackL );
613 G4double energ = aMaterialPropertiesTable->GetConstProperty( energy );
614 delta += dx*CLHEP::cm; energ += dE*CLHEP::MeV;
615 aMaterialPropertiesTable->AddConstProperty( trackL, delta );
616 aMaterialPropertiesTable->AddConstProperty( energy, energ );
617 if ( TotalEnergyDeposit > 0 ) {
618 G4double deltaTime = aMaterialPropertiesTable->
619 GetConstProperty( time00 );
622 if (aParticle->GetCharge() != 0) {
624 aMaterialPropertiesTable->AddConstProperty( time00, t0 );
628 aMaterialPropertiesTable->AddConstProperty( time00, t1 );
630 deltaTime = aMaterialPropertiesTable->GetConstProperty( time01 );
633 aMaterialPropertiesTable->AddConstProperty( time01, t1 );
637 TotalEnergyDeposit=aMaterialPropertiesTable->
638 GetConstProperty(
"ENERGY_DEPOSIT_TOT");
639 InitialKinetEnergy=aMaterialPropertiesTable->
640 GetConstProperty(
"ENERGY_DEPOSIT_GOL");
641 if(InitialKinetEnergy >
HIENLIM &&
642 abs(aParticle->GetPDGcode()) != 2112) fVeryHighEnergy=
true;
644 if (fVeryHighEnergy && !fExcitedNucleus) safety = 0.2*CLHEP::keV;
645 else safety = 2.*CLHEP::eV;
648 if( !anExcitationEnergy && pDef->GetParticleType() ==
"nucleus" &&
649 aTrack.GetTrackStatus() != fAlive && !fAlpha )
650 InitialKinetEnergy = TotalEnergyDeposit;
651 if ( particleName ==
"neutron" || particleName ==
"antineutron" )
652 InitialKinetEnergy = TotalEnergyDeposit;
657 if( std::abs(TotalEnergyDeposit-InitialKinetEnergy)<safety ||
658 TotalEnergyDeposit>=InitialKinetEnergy ){
662 NumPhotons = 0; NumElectrons = 0;
664 sprintf(numPho,
"N_PHO_%d",i); sprintf(numEle,
"N_ELE_%d",i);
665 sprintf(trackL,
"TRACK_%d",i); sprintf(energy,
"ENRGY_%d",i);
667 dx += aMaterialPropertiesTable->GetConstProperty(trackL);
668 dE += aMaterialPropertiesTable->GetConstProperty(energy);
671 G4int buffer = 100;
if ( fVeryHighEnergy ) buffer = 1;
672 fParticleChange.SetNumberOfSecondaries(buffer*(NumPhotons+NumElectrons));
674 if (fTrackSecondariesFirst) {
675 if (aTrack.GetTrackStatus() == fAlive )
684 sprintf(xCoord,
"POS_X_%d",i); sprintf(yCoord,
"POS_Y_%d",i);
685 sprintf(zCoord,
"POS_Z_%d",i);
686 sprintf(numExc,
"N_EXC_%d",i); sprintf(numIon,
"N_ION_%d",i);
687 sprintf(numPho,
"N_PHO_%d",i); sprintf(numEle,
"N_ELE_%d",i);
688 NumExcitons = (G4int)aMaterialPropertiesTable->
689 GetConstProperty( numExc );
690 NumIons = (G4int)aMaterialPropertiesTable->
691 GetConstProperty( numIon );
692 sprintf(trackL,
"TRACK_%d",i); sprintf(energy,
"ENRGY_%d",i);
693 sprintf(time00,
"TIME0_%d",i); sprintf(time01,
"TIME1_%d",i);
694 delta = aMaterialPropertiesTable->GetConstProperty( trackL );
695 energ = aMaterialPropertiesTable->GetConstProperty( energy );
696 t0 = aMaterialPropertiesTable->GetConstProperty( time00 );
697 t1 = aMaterialPropertiesTable->GetConstProperty( time01 );
703 if ( (delta < R0 && !fVeryHighEnergy) || z2 == z1 || fAlpha ) {
704 if( z1 == 54 && ElectricField &&
705 Phase == kStateLiquid ) {
706 if ( abs ( z1 - z2 ) &&
707 abs ( aParticle->GetPDGcode() ) != 2112 ) {
708 ThomasImel = 0.056776*pow(ElectricField,-0.11844);
710 ThomasImel=0.057675*pow(ElectricField,-0.49362);
716 -0.15169*pow((ElectricField+215.25),0.01811)+0.20952;
719 if (ThomasImel > 0.19) ThomasImel = 0.19;
720 if (ThomasImel < 0.00) ThomasImel = 0.00;
722 if ( Phase == kStateLiquid ) {
723 if ( z1 == 54 ) ThomasImel *= pow((Density/2.84),0.3);
724 if ( z1 == 18 ) ThomasImel *= pow((Density/
Density_LAr),0.3);
731 xi = (G4double(NumIons)/4.)*ThomasImel;
732 if ( InitialKinetEnergy == 9.4*CLHEP::keV ) {
733 G4double NumIonsEff = 1.066e7*pow(t0/CLHEP::ns,-1.303)*
734 (0.17163+162.32/(ElectricField+191.39));
735 if ( NumIonsEff > 1e6 ) NumIonsEff = 1e6;
736 xi = (G4double(NumIonsEff)/4.)*ThomasImel;
738 if ( fKr83m && ElectricField==0 )
739 xi = (G4double(1.3*NumIons)/4.)*ThomasImel;
740 recombProb = 1-log(1+xi)/xi;
741 if(recombProb<0) recombProb=0;
742 if(recombProb>1) recombProb=1;
745 NumPhotons = NumExcitons +
BinomFluct(NumIons,recombProb);
746 NumElectrons = (NumExcitons + NumIons) - NumPhotons;
748 aMaterialPropertiesTable->
749 AddConstProperty( numPho, NumPhotons );
750 aMaterialPropertiesTable->
751 AddConstProperty( numEle, NumElectrons );
756 NumPhotons = (G4int)aMaterialPropertiesTable->
757 GetConstProperty( numPho );
758 NumElectrons =(G4int)aMaterialPropertiesTable->
759 GetConstProperty( numEle );
761 G4double FanoFactor =0;
764 2575.9*pow((ElectricField+15.154),-0.64064)-1.4707;
765 if ( fKr83m ) TotalEnergyDeposit = 4*CLHEP::keV;
766 if ( (dE/CLHEP::keV) <= 100 && ElectricField >= 0 ) {
767 G4double keVee = (TotalEnergyDeposit/(100.*CLHEP::keV));
769 FanoFactor *= -0.00075+0.4625*keVee+34.375*pow(keVee,2.);
771 FanoFactor *= 0.069554+1.7322*keVee-.80215*pow(keVee,2.);
774 if ( Phase == kStateGas && Density>0.5 ) FanoFactor =
775 0.42857-4.7857*Density+7.8571*pow(Density,2.);
776 if( FanoFactor <= 0 || fVeryHighEnergy ) FanoFactor = 0;
777 NumQuanta = NumPhotons + NumElectrons;
778 if(z1==54 && FanoFactor) NumElectrons = G4int(
779 floor(GaussGen.fire(NumElectrons,
780 sqrt(FanoFactor*NumElectrons))+0.5));
781 NumPhotons = NumQuanta - NumElectrons;
782 if ( NumElectrons <= 0 ) NumElectrons = 0;
783 if ( NumPhotons <= 0 ) NumPhotons = 0;
788 } NumElectrons = G4int(floor(NumElectrons*
phe_per_e+0.5));
794 if(fKr83m || InitialKinetEnergy==9.4*CLHEP::keV) fKr83m += dE/CLHEP::keV;
795 if(fKr83m > 41) fKr83m = 0;
801 aMaterialPropertiesTable->AddConstProperty( numExc, 0 );
802 aMaterialPropertiesTable->AddConstProperty( numIon, 0 );
803 aMaterialPropertiesTable->AddConstProperty( numPho, 0 );
804 aMaterialPropertiesTable->AddConstProperty( numEle, 0 );
807 if( InitialKinetEnergy < MAX_ENE && InitialKinetEnergy >
MIN_ENE &&
808 !fMultipleScattering )
809 NumQuanta = NumPhotons + NumElectrons;
811 for(k = 0; k < NumQuanta; k++) {
812 G4double sampledEnergy;
813 std::unique_ptr<G4DynamicParticle> aQuantum;
816 G4double cost = 1. - 2.*UniformGen.fire();
817 G4double sint = std::sqrt((1.-cost)*(1.+cost));
818 G4double phi = CLHEP::twopi*UniformGen.fire();
819 G4double sinp = std::sin(phi); G4double cosp = std::cos(phi);
820 G4double px = sint*cosp; G4double py = sint*sinp;
824 G4ParticleMomentum photonMomentum(px, py, pz);
827 if (k < NumPhotons) {
829 G4double sx = cost*cosp;
830 G4double sy = cost*sinp;
832 G4ThreeVector photonPolarization(sx, sy, sz);
833 G4ThreeVector perp = photonMomentum.cross(photonPolarization);
834 phi = CLHEP::twopi*UniformGen.fire();
835 sinp = std::sin(phi);
836 cosp = std::cos(phi);
837 photonPolarization = cosp * photonPolarization + sinp * perp;
838 photonPolarization = photonPolarization.unit();
841 sampledEnergy = GaussGen.fire(PhotMean,PhotWidth);
842 aQuantum = std::make_unique<G4DynamicParticle>(G4OpticalPhoton::OpticalPhoton(),
844 aQuantum->SetPolarization(photonPolarization.x(),
845 photonPolarization.y(),
846 photonPolarization.z());
852 G4ParticleMomentum electronMomentum(0, 0, -FieldSign);
855 if ( Phase == kStateGas ) {
866 sampledEnergy = 1.38e-23*(CLHEP::joule/CLHEP::kelvin)*Temperature;
871 aQuantum->SetKineticEnergy(sampledEnergy);
880 G4double aSecondaryTime = t0+UniformGen.fire()*(t1-t0)+evtStrt;
881 if (tau1<0) { tau1=0; }
882 if (tau3<0) { tau3=0; }
883 if (tauR<0) { tauR=0; }
884 if ( aQuantum->GetDefinition()->
885 GetParticleName()==
"opticalphoton" ) {
886 if ( abs(z2-z1) && !fAlpha &&
887 abs(aParticle->GetPDGcode()) != 2112 ) {
888 LET = (energ/CLHEP::MeV)/(delta/CLHEP::cm)*(1/Density);
894 if ( Phase == kStateLiquid && z1 == 54 )
895 tauR = 3.5*((1+0.41*LET)/(0.18*LET))*CLHEP::ns
896 *exp(-0.00900*ElectricField);
903 SingTripRatioX = GaussGen.fire(0.17,0.05);
904 SingTripRatioR = GaussGen.fire(0.8,0.2);
906 SingTripRatioR = 0.2701+0.003379*LET-4.7338e-5*pow(LET,2.)
907 +8.1449e-6*pow(LET,3.); SingTripRatioX = SingTripRatioR;
909 SingTripRatioX = GaussGen.fire(0.36,0.06);
910 SingTripRatioR = GaussGen.fire(0.5,0.2); }
914 SingTripRatioR = GaussGen.fire(2.3,0.51);
917 if (z1==18) SingTripRatioR = (-0.065492+1.9996
918 *exp(-dE/CLHEP::MeV))/(1+0.082154/pow(dE/CLHEP::MeV,2.)) + 2.1811;
919 SingTripRatioX = SingTripRatioR;
924 SingTripRatioR = GaussGen.fire(7.8,1.5);
925 if (z1==18) SingTripRatioR = 0.22218*pow(energ/CLHEP::keV,0.48211);
926 SingTripRatioX = SingTripRatioR;
931 if ( k > NumExcitons ) {
934 aSecondaryTime += tauR*(1./UniformGen.fire()-1);
935 if(UniformGen.fire()<SingTripRatioR/(1+SingTripRatioR))
936 aSecondaryTime -= tau1*log(UniformGen.fire());
937 else aSecondaryTime -= tau3*log(UniformGen.fire());
940 if(UniformGen.fire()<SingTripRatioX/(1+SingTripRatioX))
941 aSecondaryTime -= tau1*log(UniformGen.fire());
942 else aSecondaryTime -= tau3*log(UniformGen.fire());
946 G4double gainField = 12;
947 G4double tauTrap = 884.83-62.069*gainField;
948 if ( Phase == kStateLiquid )
949 aSecondaryTime -= tauTrap*CLHEP::ns*log(UniformGen.fire());
958 std::string sx(xCoord);
959 std::string sy(yCoord);
960 std::string sz(zCoord);
961 x0[0] = aMaterialPropertiesTable->GetConstProperty( xCoord );
962 x0[1] = aMaterialPropertiesTable->GetConstProperty( yCoord );
963 x0[2] = aMaterialPropertiesTable->GetConstProperty( zCoord );
964 G4double
radius = sqrt(pow(x0[0],2.)+pow(x0[1],2.));
967 if ( radius >=
R_TOL ) {
968 if (x0[0] == 0) { x0[0] = 1*CLHEP::nm; }
969 if (x0[1] == 0) { x0[1] = 1*CLHEP::nm; }
970 radius -=
R_TOL; phi = atan ( x0[1] / x0[0] );
971 x0[0] = fabs(radius*cos(phi))*((fabs(x0[0]))/(x0[0]));
972 x0[1] = fabs(radius*sin(phi))*((fabs(x0[1]))/(x0[1]));
975 G4ThreeVector aSecondaryPosition = x0;
976 if ( k >= NumPhotons &&
diffusion && ElectricField > 0 ) {
977 G4double D_T = 64*pow(1
e-3*ElectricField,-.17);
979 G4double D_L = 13.859*pow(1
e-3*ElectricField,-0.58559);
981 if ( Phase == kStateLiquid && z1 == 18 ) {
982 D_T = 93.342*pow(ElectricField/nDensity,0.041322);
984 if ( Phase == kStateGas && z1 == 54 ) {
985 D_L=4.265+19097/ElectricField-1.7397e6/pow(ElectricField,2.)+
986 1.2477e8/pow(ElectricField,3.); D_T *= 0.01;
988 if (ElectricField<100 && Phase == kStateLiquid) D_L = 8*CLHEP::cm2/
CLHEP::s;
989 G4double vDrift = sqrt((2*sampledEnergy)/(
EMASS));
990 if (
BORDER == 0 ) x0[2] = 0;
991 G4double sigmaDT = sqrt(2*D_T*fabs(
BORDER-x0[2])/vDrift);
992 G4double sigmaDL = sqrt(2*D_L*fabs(
BORDER-x0[2])/vDrift);
993 G4double dr = std::abs(GaussGen.fire(0.,sigmaDT));
994 phi = CLHEP::twopi * UniformGen.fire();
995 aSecondaryPosition[0] += cos(phi) * dr;
996 aSecondaryPosition[1] += sin(phi) * dr;
997 aSecondaryPosition[2] += GaussGen.fire(0.,sigmaDL);
998 radius = std::sqrt(std::pow(aSecondaryPosition[0],2.)+
999 std::pow(aSecondaryPosition[1],2.));
1000 if(aSecondaryPosition[2] >=
BORDER && Phase == kStateLiquid) {
1002 if(aSecondaryPosition[2] <=
PMT && !GlobalFields)
1003 aSecondaryPosition[2] =
PMT +
R_TOL;
1007 if ( aSecondaryTime < 0 ) aSecondaryTime = 0;
1018 aMaterialPropertiesTable->AddConstProperty( xCoord, 999*CLHEP::km );
1019 aMaterialPropertiesTable->AddConstProperty( yCoord, 999*CLHEP::km );
1020 aMaterialPropertiesTable->AddConstProperty( zCoord, 999*CLHEP::km );
1021 aMaterialPropertiesTable->AddConstProperty( trackL, 0*CLHEP::um );
1022 aMaterialPropertiesTable->AddConstProperty( energy, 0*CLHEP::eV );
1023 aMaterialPropertiesTable->AddConstProperty( time00, DBL_MAX );
1024 aMaterialPropertiesTable->AddConstProperty( time01, -1*CLHEP::ns );
1029 aMaterialPropertiesTable->
1030 AddConstProperty(
"TOTALNUM_INT_SITES", 0 );
1031 aMaterialPropertiesTable->
1032 AddConstProperty(
"ENERGY_DEPOSIT_TOT", 0*CLHEP::keV );
1033 aMaterialPropertiesTable->
1034 AddConstProperty(
"ENERGY_DEPOSIT_GOL", 0*CLHEP::MeV );
1035 fExcitedNucleus =
false;
1052 std::cout <<
"WARNING: NestAlg::GetGasElectronDriftSpeed(G4double, G4double) " 1053 <<
"is not defined, returning bogus value of -999." << std::endl;
1061 G4double efieldinput,
1064 if(efieldinput<0) efieldinput *= (-1);
1066 G4double onea=144623.235704015,
1067 oneb=850.812714257629,
1068 onec=1192.87056676815,
1069 oned=-395969.575204061,
1070 onef=-355.484170008875,
1071 oneg=-227.266219627672,
1072 oneh=223831.601257495,
1073 onei=6.1778950907965,
1074 onej=18.7831533426398,
1075 onek=-76132.6018884368;
1077 G4double twoa=17486639.7118995,
1078 twob=-113.174284723134,
1079 twoc=28.005913193763,
1080 twod=167994210.094027,
1081 twof=-6766.42962575088,
1082 twog=901.474643115395,
1083 twoh=-185240292.471665,
1084 twoi=-633.297790813084,
1085 twoj=87.1756135457949;
1087 G4double thra=10626463726.9833,
1088 thrb=224025158.134792,
1089 thrc=123254826.300172,
1090 thrd=-4563.5678061122,
1091 thrf=-1715.269592063,
1092 thrg=-694181.921834368,
1093 thrh=-50.9753281079838,
1094 thri=58.3785811395493,
1095 thrj=201512.080026704;
1097 t1=0,
t2=0,slope=0,intercept=0;
1100 f1=onea/(1+exp(-(efieldinput-oneb)/onec))+oned/
1101 (1+exp(-(efieldinput-onef)/oneg))+
1102 oneh/(1+exp(-(efieldinput-onei)/onej))+onek;
1103 f2=twoa/(1+exp(-(efieldinput-twob)/twoc))+twod/
1104 (1+exp(-(efieldinput-twof)/twog))+
1105 twoh/(1+exp(-(efieldinput-twoi)/twoj));
1106 f3=thra*exp(-thrb*efieldinput)+thrc*exp(-(pow(efieldinput-thrd,2))/
1108 thrg*exp(-(pow(efieldinput-thrh,2)/(thri*thri)))+thrj;
1110 if(efieldinput<20 && efieldinput>=0) {
1111 f1=2951*efieldinput;
1112 f2=5312*efieldinput;
1113 f3=7101*efieldinput;
1116 if(tempinput<200.0 && tempinput>165.0) {
1122 if(tempinput<230.0 && tempinput>200.0) {
1128 if((tempinput>230.0 || tempinput<165.0) && !Miller) {
1129 G4cout <<
"\nWARNING: TEMPERATURE OUT OF RANGE (165-230 K)\n";
1132 if (tempinput == 165.0) edrift =
f1;
1133 else if (tempinput == 200.0) edrift =
f2;
1134 else if (tempinput == 230.0) edrift =
f3;
1138 intercept=y1-slope*
t1;
1139 edrift=slope*tempinput+intercept;
1143 if ( efieldinput <= 40. )
1144 edrift = -0.13274+0.041082*efieldinput-0.0006886*pow(efieldinput,2.)+
1145 5.5503e-6*pow(efieldinput,3.);
1147 edrift = 0.060774*efieldinput/pow(1+0.11336*pow(efieldinput,0.5218),2.);
1148 if ( efieldinput >= 1e5 ) edrift = 2.7;
1149 if ( efieldinput >= 100 )
1150 edrift -= 0.017 * ( tempinput - 163 );
1152 edrift += 0.017 * ( tempinput - 163 );
1155 if ( Z == 18 ) edrift = 1e5 * (.097384*pow(log10(efieldinput),3.0622)-.018614*sqrt(efieldinput) );
1156 if ( edrift < 0 ) edrift = 0.;
1167 if ( E >= 1 ) LET = 58.482-61.183*log10(E)+19.749*pow(log10(E),2)+
1168 2.3101*pow(log10(E),3)-3.3469*pow(log10(E),4)+
1169 0.96788*pow(log10(E),5)-0.12619*pow(log10(E),6)+0.0065108*pow(log10(E),7);
1173 else if ( E>0 && E<1 ) LET = 6.9463+815.98*E-4828*pow(E,2)+17079*pow(E,3)-
1174 36394*pow(E,4)+44553*pow(E,5)-28659*pow(E,6)+7483.8*pow(E,7);
1179 if ( E >= 1 ) LET = 116.70-162.97*log10(E)+99.361*pow(log10(E),2)-
1180 33.405*pow(log10(E),3)+6.5069*pow(log10(E),4)-
1181 0.69334*pow(log10(E),5)+.031563*pow(log10(E),6);
1182 else if ( E>0 && E<1 ) LET = 100;
1191 CLHEP::RandGauss GaussGen(
fEngine);
1192 CLHEP::RandFlat UniformGen(
fEngine);
1194 G4double
mean = N0*prob;
1195 G4double sigma = sqrt(N0*prob*(1-prob));
1197 if ( prob == 0.00 )
return N1;
1198 if ( prob == 1.00 )
return N0;
1201 for(G4int i = 0; i < N0; i++) {
1202 if(UniformGen.fire() < prob) N1++;
1206 N1 = G4int(floor(GaussGen.fire(mean,sigma)+0.5));
1208 if ( N1 > N0 ) N1 = N0;
1209 if ( N1 < 0 ) N1 = 0;
1216 char xCoord[80];
char yCoord[80];
char zCoord[80];
1217 char numExc[80];
char numIon[80];
char numPho[80];
char numEle[80];
1218 char trackL[80];
char time00[80];
char time01[80];
char energy[80];
1221 for( G4int i=0; i<10000; i++ ) {
1222 sprintf(xCoord,
"POS_X_%d",i); sprintf(yCoord,
"POS_Y_%d",i);
1223 sprintf(zCoord,
"POS_Z_%d",i);
1224 nobleElementMat->AddConstProperty( xCoord, 999*CLHEP::km );
1225 nobleElementMat->AddConstProperty( yCoord, 999*CLHEP::km );
1226 nobleElementMat->AddConstProperty( zCoord, 999*CLHEP::km );
1227 sprintf(numExc,
"N_EXC_%d",i); sprintf(numIon,
"N_ION_%d",i);
1228 sprintf(numPho,
"N_PHO_%d",i); sprintf(numEle,
"N_ELE_%d",i);
1229 nobleElementMat->AddConstProperty( numExc, 0 );
1230 nobleElementMat->AddConstProperty( numIon, 0 );
1231 nobleElementMat->AddConstProperty( numPho, 0 );
1232 nobleElementMat->AddConstProperty( numEle, 0 );
1233 sprintf(trackL,
"TRACK_%d",i); sprintf(energy,
"ENRGY_%d",i);
1234 sprintf(time00,
"TIME0_%d",i); sprintf(time01,
"TIME1_%d",i);
1235 nobleElementMat->AddConstProperty( trackL, 0*CLHEP::um );
1236 nobleElementMat->AddConstProperty( energy, 0*CLHEP::eV );
1237 nobleElementMat->AddConstProperty( time00, DBL_MAX );
1238 nobleElementMat->AddConstProperty( time01,-1*CLHEP::ns );
1244 nobleElementMat->AddConstProperty(
"TOTALNUM_INT_SITES", 0 );
1245 nobleElementMat->AddConstProperty(
"ENERGY_DEPOSIT_TOT", 0*CLHEP::keV );
1246 nobleElementMat->AddConstProperty(
"ENERGY_DEPOSIT_GOL", 0*CLHEP::MeV );
1255 G4double a_0 = 5.29e-11*CLHEP::m; G4double a = 0.626*a_0*pow(Z,(-1./3.));
1256 G4double epsilon_0 = 8.854e-12*(CLHEP::farad/CLHEP::m);
1257 G4double epsilon = (a*E*2.*CLHEP::twopi*epsilon_0)/(2*pow(CLHEP::eplus,2.)*pow(Z,2.));
1258 G4double zeta_0 = pow(Z,(1./6.)); G4double m_N = A*1.66e-27*CLHEP::kg;
1259 G4double hbar = 6.582e-16*CLHEP::eV*
CLHEP::s;
1263 G4double s_n = log(1+1.1383*epsilon)/(2.*(epsilon +
1264 0.01321*pow(epsilon,0.21226) +
1265 0.19593*sqrt(epsilon)));
1266 G4double s_e = (a_0*zeta_0/a)*hbar*sqrt(8*epsilon*2.*CLHEP::twopi*epsilon_0/
1267 (a*m_N*pow(CLHEP::eplus,2.)));
1268 return 1.38e5*0.5*(1+tanh(50*epsilon-0.25))*epsilon*(s_e/s_n);
code to link reconstructed objects back to the MC truth information
Float_t y1[n_points_granero]
G4double GetLiquidElectronDriftSpeed(double T, double F, G4bool M, G4int Z)
Float_t x1[n_points_granero]
NestAlg(CLHEP::HepRandomEngine &engine)
double fYieldFactor
turns scint. on/off
G4VParticleChange fParticleChange
pointer to G4VParticleChange
int fNumIonElectrons
number of ionization electrons produced by step
G4int BinomFluct(G4int N0, G4double prob)
Float_t y2[n_points_geant4]
int fNumScintPhotons
number of photons produced by the step
static G4ThermalElectron * ThermalElectron()
const G4VParticleChange & CalculateIonizationAndScintillation(G4Track const &aTrack, G4Step const &aStep)
G4double UnivScreenFunc(G4double E, G4double Z, G4double A)
G4double GetGasElectronDriftSpeed(G4double efieldinput, G4double density)
CLHEP::HepRandomEngine & fEngine
random engine
G4double CalculateElectronLET(G4double E, G4int Z)
double mean(const std::vector< short > &wf, size_t start, size_t nsample)
void InitMatPropValues(G4MaterialPropertiesTable *nobleElementMat, int z)
std::map< int, bool > fElementPropInit
double fEnergyDep
energy deposited by the step