Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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//#include "G4NuMuNucleusCcModel.hh"
49//#include "G4NuMuNucleusNcModel.hh"
50
51#include "G4RotationMatrix.hh"
52#include "G4ThreeVector.hh"
53#include "G4AffineTransform.hh"
54#include "G4DynamicParticle.hh"
55#include "G4StepPoint.hh"
56#include "G4VSolid.hh"
57#include "G4LogicalVolume.hh"
58#include "G4SafetyHelper.hh"
60
61///////////////////////////////////////////////////////////////////////////////
62
63
65 : G4HadronicProcess( pName, fHadronInelastic ), isInitialised(false), fBiased(true) // fHadronElastic???
66{
67 lowestEnergy = 1.*keV;
68 fEnvelope = nullptr;
69 fEnvelopeName = anEnvelopeName;
70 fTotXsc = nullptr; // new G4TauNeutrinoNucleusTotXsc();
71 fNuNuclCcBias=1.;
72 fNuNuclNcBias=1.;
73 fNuNuclTotXscBias=1.;
75 safetyHelper->InitialiseHelper();
76}
77
79{
80 if( fTotXsc ) delete fTotXsc;
81}
82
83///////////////////////////////////////////////////////
84
86{
87 fNuNuclTotXscBias = bf;
88
89 fTotXsc = new G4TauNeutrinoNucleusTotXsc();
90 fTotXsc->SetBiasingFactor(bf);
91}
92
93///////////////////////////////////////////////////////
94
96{
97 fNuNuclCcBias = bfCc;
98 fNuNuclNcBias = bfNc;
99
100 fTotXsc = new G4TauNeutrinoNucleusTotXsc();
101 // fTotXsc->SetBiasingFactors(bfCc, bfNc);
102}
103
104//////////////////////////////////////////////////
105
108{
109 //G4cout << "GetMeanFreePath " << aTrack.GetDefinition()->GetParticleName()
110 // << " Ekin= " << aTrack.GetKineticEnergy() << G4endl;
112 G4double totxsc(0.);
113
114 if( rName == fEnvelopeName && fNuNuclTotXscBias > 1.)
115 {
116 totxsc = fNuNuclTotXscBias*
118 aTrack.GetMaterial());
119 }
120 else
121 {
123 aTrack.GetMaterial());
124 }
125 G4double res = (totxsc>0.0) ? 1.0/totxsc : DBL_MAX;
126 //G4cout << " xsection= " << totxsc << G4endl;
127 return res;
128}
129
130///////////////////////////////////////////////////
131
132void G4TauNeutrinoNucleusProcess::ProcessDescription(std::ostream& outFile) const
133{
134
135 outFile << "G4TauNeutrinoNucleusProcess handles the inelastic scattering of \n"
136 << "tau-neutrino on nucleus by invoking the following model(s) and \n"
137 << "cross section(s).\n";
138
139}
140
141///////////////////////////////////////////////////////////////////////
142
145{
146 // track.GetVolume()->GetLogicalVolume()->GetName()
147 // if( track.GetVolume()->GetLogicalVolume() != fEnvelope )
148
150
151 if( rName != fEnvelopeName )
152 {
153 if( verboseLevel > 0 )
154 {
155 G4cout<<"Go out from G4TauNeutrinoNucleusProcess::PostStepDoIt: wrong volume "<<G4endl;
156 }
157 return G4VDiscreteProcess::PostStepDoIt( track, step );
158 }
161 G4double weight = track.GetWeight();
163
164 if( track.GetTrackStatus() != fAlive )
165 {
166 return theTotalResult;
167 }
168 // Next check for illegal track status
169 //
170 if (track.GetTrackStatus() != fAlive &&
171 track.GetTrackStatus() != fSuspend)
172 {
173 if (track.GetTrackStatus() == fStopAndKill ||
176 {
178 ed << "G4TauNeutrinoNucleusProcess: track in unusable state - "
179 << track.GetTrackStatus() << G4endl;
180 ed << "G4TauNeutrinoNucleusProcess: returning unchanged track " << G4endl;
181 DumpState(track,"PostStepDoIt",ed);
182 G4Exception("G4TauNeutrinoNucleusProcess::PostStepDoIt", "had004", JustWarning, ed);
183 }
184 // No warning for fStopButAlive which is a legal status here
185 return theTotalResult;
186 }
187
188 // For elastic scattering, _any_ result is considered an interaction
190
191 G4double kineticEnergy = track.GetKineticEnergy();
192 const G4DynamicParticle* dynParticle = track.GetDynamicParticle();
193 const G4ParticleDefinition* part = dynParticle->GetDefinition();
194 const G4String pName = part->GetParticleName();
195
196 // NOTE: Very low energy scatters were causing numerical (FPE) errors
197 // in earlier releases; these limits have not been changed since.
198
199 if ( kineticEnergy <= lowestEnergy ) return theTotalResult;
200
201 const G4Material* material = track.GetMaterial();
202 G4Nucleus* targNucleus = GetTargetNucleusPointer();
203
204 //////////////// uniform random spread of the neutrino interaction point ////////////
205
206 const G4StepPoint* pPostStepPoint = step.GetPostStepPoint();
207 const G4DynamicParticle* aParticle = track.GetDynamicParticle();
208 G4ThreeVector position = pPostStepPoint->GetPosition(), newPosition=position;
209 G4ParticleMomentum direction = aParticle->GetMomentumDirection();
210
211 if( fNuNuclCcBias > 1.0 || fNuNuclNcBias > 1.0) // = true, if fBiasingfactor != 1., i.e. xsc is biased
212 {
213 const G4RotationMatrix* rotM = pPostStepPoint->GetTouchable()->GetRotation();
214 G4ThreeVector transl = pPostStepPoint->GetTouchable()->GetTranslation();
215 G4AffineTransform transform = G4AffineTransform(rotM,transl);
216 transform.Invert();
217
218 G4ThreeVector localP = transform.TransformPoint(position);
219 G4ThreeVector localV = transform.TransformAxis(direction);
220
221 G4double forward = track.GetVolume()->GetLogicalVolume()->GetSolid()->DistanceToOut(localP, localV);
222 G4double backward = track.GetVolume()->GetLogicalVolume()->GetSolid()->DistanceToOut(localP, -localV);
223
224 G4double distance = forward+backward;
225
226 // G4cout<<distance/cm<<", ";
227
228 // uniform sampling of nu-e interaction point
229 // along neutrino direction in current volume
230
231 G4double range = -backward+G4UniformRand()*distance;
232
233 newPosition = position + range*direction;
234
235 safetyHelper->ReLocateWithinVolume(newPosition);
236
237 theTotalResult->ProposePosition(newPosition); // G4Exception : GeomNav1002
238 }
239 G4HadProjectile theProj( track );
240 G4HadronicInteraction* hadi = nullptr;
241 G4HadFinalState* result = nullptr;
242
243 G4double ccTotRatio = fTotXsc->GetCcTotRatio();
244
245 if( G4UniformRand() < ccTotRatio ) // Cc-model
246 {
247 // Initialize the hadronic projectile from the track
248 thePro.Initialise(track);
249
250 if (pName == "nu_tau" ) hadi = (GetHadronicInteractionList())[0];
251 else hadi = (GetHadronicInteractionList())[2];
252
253 result = hadi->ApplyYourself( thePro, *targNucleus);
254
256
258
259 FillResult(result, track);
260 }
261 else // Nc-model
262 {
263
264 if (pName == "nu_tau" ) hadi = (GetHadronicInteractionList())[1];
265 else hadi = (GetHadronicInteractionList())[3];
266
267 size_t idx = track.GetMaterialCutsCouple()->GetIndex();
268
270
271 hadi->SetRecoilEnergyThreshold(tcut);
272
273 if( verboseLevel > 1 )
274 {
275 G4cout << "G4TauNeutrinoNucleusProcess::PostStepDoIt for "
276 << part->GetParticleName()
277 << " in " << material->GetName()
278 << " Target Z= " << targNucleus->GetZ_asInt()
279 << " A= " << targNucleus->GetA_asInt() << G4endl;
280 }
281 try
282 {
283 result = hadi->ApplyYourself( theProj, *targNucleus);
284 }
285 catch(G4HadronicException & aR)
286 {
288 aR.Report(ed);
289 ed << "Call for " << hadi->GetModelName() << G4endl;
290 ed << " Z= "
291 << targNucleus->GetZ_asInt()
292 << " A= " << targNucleus->GetA_asInt() << G4endl;
293 DumpState(track,"ApplyYourself",ed);
294 ed << " ApplyYourself failed" << G4endl;
295 G4Exception("G4TauNeutrinoNucleusProcess::PostStepDoIt", "had006",
296 FatalException, ed);
297 }
298 // directions
299
300 G4ThreeVector indir = track.GetMomentumDirection();
301 G4double phi = CLHEP::twopi*G4UniformRand();
302 G4ThreeVector it(0., 0., 1.);
303 G4ThreeVector outdir = result->GetMomentumChange();
304
305 if(verboseLevel>1)
306 {
307 G4cout << "Efin= " << result->GetEnergyChange()
308 << " de= " << result->GetLocalEnergyDeposit()
309 << " nsec= " << result->GetNumberOfSecondaries()
310 << " dir= " << outdir
311 << G4endl;
312 }
313 // energies
314
315 G4double edep = result->GetLocalEnergyDeposit();
316 G4double efinal = result->GetEnergyChange();
317
318 if(efinal < 0.0) { efinal = 0.0; }
319 if(edep < 0.0) { edep = 0.0; }
320
321 // NOTE: Very low energy scatters were causing numerical (FPE) errors
322 // in earlier releases; these limits have not been changed since.
323
324 if(efinal <= lowestEnergy)
325 {
326 edep += efinal;
327 efinal = 0.0;
328 }
329 // primary change
330
332
333 G4TrackStatus status = track.GetTrackStatus();
334
335 if(efinal > 0.0)
336 {
337 outdir.rotate(phi, it);
338 outdir.rotateUz(indir);
340 }
341 else
342 {
343 if( part->GetProcessManager()->GetAtRestProcessVector()->size() > 0)
344 {
345 status = fStopButAlive;
346 }
347 else
348 {
349 status = fStopAndKill;
350 }
352 }
353 //G4cout << "Efinal= " << efinal << " TrackStatus= " << status << G4endl;
354
356
357 // recoil
358
359 if( result->GetNumberOfSecondaries() > 0 )
360 {
361 G4DynamicParticle* p = result->GetSecondary(0)->GetParticle();
362
363 if(p->GetKineticEnergy() > tcut)
364 {
367
368 // G4cout << "recoil " << pdir << G4endl;
369 //!! is not needed for models inheriting G4TauNeutrinoNucleus
370
371 pdir.rotate(phi, it);
372 pdir.rotateUz(indir);
373
374 // G4cout << "recoil rotated " << pdir << G4endl;
375
376 p->SetMomentumDirection(pdir);
377
378 // in elastic scattering time and weight are not changed
379
380 G4Track* t = new G4Track(p, track.GetGlobalTime(),
381 track.GetPosition());
382 t->SetWeight(weight);
385 }
386 else
387 {
388 edep += p->GetKineticEnergy();
389 delete p;
390 }
391 }
394 result->Clear();
395 }
396 return theTotalResult;
397}
398
399void
401{
402 if(!isInitialised) {
403 isInitialised = true;
404 // if(G4Neutron::Neutron() == &part) { lowestEnergy = 1.e-6*eV; }
405 }
407}
408
409void
411{
412 lowestEnergy = val;
413}
414
@ JustWarning
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
G4ForceCondition
@ fHadronInelastic
G4TrackStatus
@ fKillTrackAndSecondaries
@ fSuspend
@ fAlive
@ fStopAndKill
@ fStopButAlive
@ fPostponeToNextEvent
double G4double
Definition: G4Types.hh:83
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:33
Hep3Vector & rotate(double, const Hep3Vector &)
Definition: ThreeVectorR.cc:24
G4AffineTransform & Invert()
G4ThreeVector TransformPoint(const G4ThreeVector &vec) const
G4ThreeVector TransformAxis(const G4ThreeVector &axis) const
G4double ComputeCrossSection(const G4DynamicParticle *, const G4Material *)
void SetMomentumDirection(const G4ThreeVector &aDirection)
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
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()
void PreparePhysicsTable(const G4ParticleDefinition &) override
G4CrossSectionDataStore * GetCrossSectionDataStore()
void DumpState(const G4Track &, const G4String &, G4ExceptionDescription &)
G4VSolid * GetSolid() const
G4Region * GetRegion() const
const G4String & GetName() const
Definition: G4Material.hh:172
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)
void InitialiseHelper()
const G4VTouchable * GetTouchable() const
const G4ThreeVector & GetPosition() const
G4VPhysicalVolume * GetPhysicalVolume() const
Definition: G4Step.hh:62
G4StepPoint * GetPreStepPoint() const
G4StepPoint * GetPostStepPoint() const
void PreparePhysicsTable(const G4ParticleDefinition &) override
G4TauNeutrinoNucleusProcess(G4String anEnvelopeName, const G4String &procName="tau-neutrino-nucleus")
void SetBiasingFactors(G4double bfCc, G4double bfNc)
void ProcessDescription(std::ostream &outFile) const override
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep) override
G4double GetMeanFreePath(const G4Track &aTrack, G4double, G4ForceCondition *) override
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 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()
Definition: G4VProcess.hh:428
G4int verboseLevel
Definition: G4VProcess.hh:360
virtual G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=nullptr, G4ThreeVector *n=nullptr) const =0
virtual const G4ThreeVector & GetTranslation(G4int depth=0) const =0
virtual const G4RotationMatrix * GetRotation(G4int depth=0) const =0
#define DBL_MAX
Definition: templates.hh:62