Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4TauNeutrinoNucleusProcess.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26//
27// Geant4 Hadron Inelastic Scattering Process
28//
29// Created from G4MuNeutrinoNucleusProcess
30//
31
32
33#include <iostream>
34#include <typeinfo>
35
37#include "G4SystemOfUnits.hh"
38#include "G4Nucleus.hh"
39#include "G4ProcessManager.hh"
45#include "G4VDiscreteProcess.hh"
46
48
49#include "G4RotationMatrix.hh"
50#include "G4ThreeVector.hh"
51#include "G4AffineTransform.hh"
52#include "G4DynamicParticle.hh"
53#include "G4StepPoint.hh"
54#include "G4VSolid.hh"
55#include "G4LogicalVolume.hh"
56#include "G4SafetyHelper.hh"
58
59///////////////////////////////////////////////////////////////////////////////
60
61
64{
65 lowestEnergy = 1.*keV;
66 fEnvelopeName = anEnvelopeName;
67 fTotXsc = new G4TauNeutrinoNucleusTotXsc();
69 safetyHelper->InitialiseHelper();
70}
71
72///////////////////////////////////////////////////////
73
75{
76 fNuNuclTotXscBias = bf;
77}
78
79///////////////////////////////////////////////////////
80
82{
83 fNuNuclCcBias = bfCc;
84 fNuNuclNcBias = bfNc;
85 fNuNuclTotXscBias = std::max(fNuNuclCcBias, fNuNuclNcBias);
86}
87
88///////////////////////////////////////////////
89
97
98//////////////////////////////////////////////////
99
102{
103 //G4cout << "GetMeanFreePath " << aTrack.GetDefinition()->GetParticleName()
104 // << " Ekin= " << aTrack.GetKineticEnergy() << G4endl;
106 G4double totxsc =
108 aTrack.GetMaterial());
109
110 if( rName == fEnvelopeName )
111 {
112 totxsc *= fNuNuclTotXscBias;
113 }
114 G4double res = (totxsc>0.0) ? 1.0/totxsc : DBL_MAX;
115 //G4cout << " xsection= " << totxsc << G4endl;
116 return res;
117}
118
119///////////////////////////////////////////////////
120
121void G4TauNeutrinoNucleusProcess::ProcessDescription(std::ostream& outFile) const
122{
123
124 outFile << "G4TauNeutrinoNucleusProcess handles the inelastic scattering of \n"
125 << "tau-neutrino on nucleus by invoking the following model(s) and \n"
126 << "cross section(s).\n";
127
128}
129
130///////////////////////////////////////////////////////////////////////
131
134{
135 // track.GetVolume()->GetLogicalVolume()->GetName()
136 // if( track.GetVolume()->GetLogicalVolume() != fEnvelope )
137
139
140 if( rName != fEnvelopeName )
141 {
142 if( verboseLevel > 0 )
143 {
144 G4cout<<"Go out from G4TauNeutrinoNucleusProcess::PostStepDoIt: wrong volume "<<G4endl;
145 }
146 return G4VDiscreteProcess::PostStepDoIt( track, step );
147 }
150 G4double weight = track.GetWeight();
152
153 if( track.GetTrackStatus() != fAlive )
154 {
155 return theTotalResult;
156 }
157 // Next check for illegal track status
158 //
159 if (track.GetTrackStatus() != fAlive &&
160 track.GetTrackStatus() != fSuspend)
161 {
162 if (track.GetTrackStatus() == fStopAndKill ||
165 {
167 ed << "G4TauNeutrinoNucleusProcess: track in unusable state - "
168 << track.GetTrackStatus() << G4endl;
169 ed << "G4TauNeutrinoNucleusProcess: returning unchanged track " << G4endl;
170 DumpState(track,"PostStepDoIt",ed);
171 G4Exception("G4TauNeutrinoNucleusProcess::PostStepDoIt", "had004", JustWarning, ed);
172 }
173 // No warning for fStopButAlive which is a legal status here
174 return theTotalResult;
175 }
176
177 // For elastic scattering, _any_ result is considered an interaction
179
180 G4double kineticEnergy = track.GetKineticEnergy();
181 const G4DynamicParticle* dynParticle = track.GetDynamicParticle();
182 const G4ParticleDefinition* part = dynParticle->GetDefinition();
183 const G4String pName = part->GetParticleName();
184
185 // NOTE: Very low energy scatters were causing numerical (FPE) errors
186 // in earlier releases; these limits have not been changed since.
187
188 if ( kineticEnergy <= lowestEnergy ) return theTotalResult;
189
190 const G4Material* material = track.GetMaterial();
191 G4Nucleus* targNucleus = GetTargetNucleusPointer();
192
193 //////////////// uniform random spread of the neutrino interaction point ////////////
194
195 const G4StepPoint* pPostStepPoint = step.GetPostStepPoint();
196 const G4DynamicParticle* aParticle = track.GetDynamicParticle();
197 G4ThreeVector position = pPostStepPoint->GetPosition(), newPosition=position;
198 G4ParticleMomentum direction = aParticle->GetMomentumDirection();
199
200 if( fNuNuclCcBias > 1.0 || fNuNuclNcBias > 1.0) // = true, if fBiasingfactor != 1., i.e. xsc is biased
201 {
202 const G4RotationMatrix* rotM = pPostStepPoint->GetTouchable()->GetRotation();
203 G4ThreeVector transl = pPostStepPoint->GetTouchable()->GetTranslation();
204 G4AffineTransform transform = G4AffineTransform(rotM,transl);
205 transform.Invert();
206
207 G4ThreeVector localP = transform.TransformPoint(position);
208 G4ThreeVector localV = transform.TransformAxis(direction);
209
210 G4double forward = track.GetVolume()->GetLogicalVolume()->GetSolid()->DistanceToOut(localP, localV);
211 G4double backward = track.GetVolume()->GetLogicalVolume()->GetSolid()->DistanceToOut(localP, -localV);
212
213 G4double distance = forward+backward;
214
215 // G4cout<<distance/cm<<", ";
216
217 // uniform sampling of nu-e interaction point
218 // along neutrino direction in current volume
219
220 G4double range = -backward+G4UniformRand()*distance;
221
222 newPosition = position + range*direction;
223
224 safetyHelper->ReLocateWithinVolume(newPosition);
225
226 theTotalResult->ProposePosition(newPosition); // G4Exception : GeomNav1002
227 }
228 G4HadProjectile theProj( track );
229 G4HadronicInteraction* hadi = nullptr;
230 G4HadFinalState* result = nullptr;
231 const G4Element* elm = GetCrossSectionDataStore()->SampleZandA(dynParticle, material,
232 *targNucleus);
233 G4int ZZ = elm->GetZasInt();
234 fTotXsc->GetElementCrossSection(dynParticle, ZZ, material);
235 G4double ccTotRatio = fTotXsc->GetCcTotRatio();
236
237 if( G4UniformRand() < ccTotRatio ) // Cc-model
238 {
239 // Initialize the hadronic projectile from the track
240 thePro.Initialise(track);
241
242 if (pName == "nu_tau" ) hadi = (GetHadronicInteractionList())[0];
243 else hadi = (GetHadronicInteractionList())[2];
244
245 result = hadi->ApplyYourself( thePro, *targNucleus);
246
248
250
251 FillResult(result, track);
252 }
253 else // Nc-model
254 {
255
256 if (pName == "nu_tau" ) hadi = (GetHadronicInteractionList())[1];
257 else hadi = (GetHadronicInteractionList())[3];
258
259 size_t idx = track.GetMaterialCutsCouple()->GetIndex();
260
262
263 hadi->SetRecoilEnergyThreshold(tcut);
264
265 if( verboseLevel > 1 )
266 {
267 G4cout << "G4TauNeutrinoNucleusProcess::PostStepDoIt for "
268 << part->GetParticleName()
269 << " in " << material->GetName()
270 << " Target Z= " << targNucleus->GetZ_asInt()
271 << " A= " << targNucleus->GetA_asInt() << G4endl;
272 }
273 try
274 {
275 result = hadi->ApplyYourself( theProj, *targNucleus);
276 }
277 catch(G4HadronicException & aR)
278 {
280 aR.Report(ed);
281 ed << "Call for " << hadi->GetModelName() << G4endl;
282 ed << " Z= "
283 << targNucleus->GetZ_asInt()
284 << " A= " << targNucleus->GetA_asInt() << G4endl;
285 DumpState(track,"ApplyYourself",ed);
286 ed << " ApplyYourself failed" << G4endl;
287 G4Exception("G4TauNeutrinoNucleusProcess::PostStepDoIt", "had006",
288 FatalException, ed);
289 }
290 // directions
291
292 G4ThreeVector indir = track.GetMomentumDirection();
293 G4double phi = CLHEP::twopi*G4UniformRand();
294 G4ThreeVector it(0., 0., 1.);
295 G4ThreeVector outdir = result->GetMomentumChange();
296
297 if(verboseLevel>1)
298 {
299 G4cout << "Efin= " << result->GetEnergyChange()
300 << " de= " << result->GetLocalEnergyDeposit()
301 << " nsec= " << result->GetNumberOfSecondaries()
302 << " dir= " << outdir
303 << G4endl;
304 }
305 // energies
306
307 G4double edep = result->GetLocalEnergyDeposit();
308 G4double efinal = result->GetEnergyChange();
309
310 if(efinal < 0.0) { efinal = 0.0; }
311 if(edep < 0.0) { edep = 0.0; }
312
313 // NOTE: Very low energy scatters were causing numerical (FPE) errors
314 // in earlier releases; these limits have not been changed since.
315
316 if(efinal <= lowestEnergy)
317 {
318 edep += efinal;
319 efinal = 0.0;
320 }
321 // primary change
322
324
325 G4TrackStatus status = track.GetTrackStatus();
326
327 if(efinal > 0.0)
328 {
329 outdir.rotate(phi, it);
330 outdir.rotateUz(indir);
332 }
333 else
334 {
335 if( part->GetProcessManager()->GetAtRestProcessVector()->size() > 0)
336 {
337 status = fStopButAlive;
338 }
339 else
340 {
341 status = fStopAndKill;
342 }
344 }
345 //G4cout << "Efinal= " << efinal << " TrackStatus= " << status << G4endl;
346
348
349 // recoil
350
351 if( result->GetNumberOfSecondaries() > 0 )
352 {
353 G4DynamicParticle* p = result->GetSecondary(0)->GetParticle();
354
355 if(p->GetKineticEnergy() > tcut)
356 {
359
360 // G4cout << "recoil " << pdir << G4endl;
361 //!! is not needed for models inheriting G4TauNeutrinoNucleus
362
363 pdir.rotate(phi, it);
364 pdir.rotateUz(indir);
365
366 // G4cout << "recoil rotated " << pdir << G4endl;
367
368 p->SetMomentumDirection(pdir);
369
370 // in elastic scattering time and weight are not changed
371
372 G4Track* t = new G4Track(p, track.GetGlobalTime(),
373 track.GetPosition());
374 t->SetWeight(weight);
377 }
378 else
379 {
380 edep += p->GetKineticEnergy();
381 delete p;
382 }
383 }
386 result->Clear();
387 }
388 return theTotalResult;
389}
390
391void
393{
394 lowestEnergy = val;
395}
396
G4double condition(const G4ErrorSymMatrix &m)
@ JustWarning
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
G4ForceCondition
G4TrackStatus
@ fKillTrackAndSecondaries
@ fSuspend
@ fAlive
@ fStopAndKill
@ fStopButAlive
@ fPostponeToNextEvent
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition Randomize.hh:52
Hep3Vector & rotateUz(const Hep3Vector &)
Hep3Vector & rotate(double, const Hep3Vector &)
G4AffineTransform & Invert()
G4ThreeVector TransformPoint(const G4ThreeVector &vec) const
G4ThreeVector TransformAxis(const G4ThreeVector &axis) const
G4double ComputeCrossSection(const G4DynamicParticle *, const G4Material *)
const G4Element * SampleZandA(const G4DynamicParticle *, const G4Material *, G4Nucleus &target)
void SetMomentumDirection(const G4ThreeVector &aDirection)
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4int GetZasInt() const
Definition G4Element.hh:120
G4double GetEnergyChange() const
void SetTrafoToLab(const G4LorentzRotation &aT)
G4double GetLocalEnergyDeposit() const
const G4ThreeVector & GetMomentumChange() const
std::size_t GetNumberOfSecondaries() const
G4HadSecondary * GetSecondary(size_t i)
void Initialise(const G4Track &aT)
G4LorentzRotation & GetTrafoToLab()
G4DynamicParticle * GetParticle()
void Report(std::ostream &aS) const
virtual G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
const G4String & GetModelName() const
void SetRecoilEnergyThreshold(G4double val)
void FillResult(G4HadFinalState *aR, const G4Track &aT)
G4HadProjectile thePro
G4Nucleus * GetTargetNucleusPointer()
G4ParticleChange * theTotalResult
std::vector< G4HadronicInteraction * > & GetHadronicInteractionList()
G4CrossSectionDataStore * GetCrossSectionDataStore()
void DumpState(const G4Track &, const G4String &, G4ExceptionDescription &)
G4VSolid * GetSolid() const
G4Region * GetRegion() const
const G4String & GetName() const
G4int GetA_asInt() const
Definition G4Nucleus.hh:99
G4int GetZ_asInt() const
Definition G4Nucleus.hh:105
void AddSecondary(G4Track *aSecondary)
void ProposePosition(G4double x, G4double y, G4double z)
void Initialize(const G4Track &) override
void ProposeEnergy(G4double finalEnergy)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
G4ProcessManager * GetProcessManager() const
const G4String & GetParticleName() const
G4ProcessVector * GetAtRestProcessVector(G4ProcessVectorTypeIndex typ=typeGPIL) const
std::size_t size() const
const std::vector< G4double > * GetEnergyCutsVector(std::size_t pcIdx) const
static G4ProductionCutsTable * GetProductionCutsTable()
const G4String & GetName() const
void ReLocateWithinVolume(const G4ThreeVector &pGlobalPoint)
const G4VTouchable * GetTouchable() const
const G4ThreeVector & GetPosition() const
G4VPhysicalVolume * GetPhysicalVolume() const
G4StepPoint * GetPreStepPoint() const
G4StepPoint * GetPostStepPoint() const
G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
void SetBiasingFactors(G4double bfCc, G4double bfNc)
void ProcessDescription(std::ostream &outFile) const override
G4TauNeutrinoNucleusProcess(const G4String &anEnvelopeName, const G4String &procName="tau-neutrino-nucleus")
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep) override
G4double GetMeanFreePath(const G4Track &aTrack, G4double, G4ForceCondition *) override
virtual G4double GetElementCrossSection(const G4DynamicParticle *dynPart, G4int Z, const G4Material *mat)
virtual const G4RotationMatrix * GetRotation(G4int depth=0) const
virtual const G4ThreeVector & GetTranslation(G4int depth=0) const
G4TrackStatus GetTrackStatus() const
G4VPhysicalVolume * GetVolume() const
G4double GetWeight() const
void SetWeight(G4double aValue)
const G4ThreeVector & GetPosition() const
void SetTouchableHandle(const G4TouchableHandle &apValue)
G4double GetGlobalTime() const
G4Material * GetMaterial() const
const G4DynamicParticle * GetDynamicParticle() const
const G4TouchableHandle & GetTouchableHandle() const
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
const G4Step * GetStep() const
static G4TransportationManager * GetTransportationManager()
G4SafetyHelper * GetSafetyHelper() const
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
void ProposeTrackStatus(G4TrackStatus status)
void ProposeNonIonizingEnergyDeposit(G4double anEnergyPart)
void ProposeWeight(G4double finalWeight)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void SetNumberOfSecondaries(G4int totSecondaries)
G4LogicalVolume * GetLogicalVolume() const
void ClearNumberOfInteractionLengthLeft()
G4int verboseLevel
virtual G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=nullptr, G4ThreeVector *n=nullptr) const =0
#define DBL_MAX
Definition templates.hh:62