Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4MuNeutrinoNucleusProcess.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 Elastic Scattering Process
28//
29// Created from G4HadronElasticProcess
30//
31// Modified:
32//
33// 2.2.19 V.Grichine - PostStepDoIt implementation
34// 24.04.19 V. Grichine - G4Region name and optionally total cross section biased in the region only.
35
36#include <iostream>
37#include <typeinfo>
38
40#include "G4SystemOfUnits.hh"
41#include "G4Nucleus.hh"
43#include "G4ProcessManager.hh"
48#include "G4VDiscreteProcess.hh"
49
51
52#include "G4RotationMatrix.hh"
53#include "G4ThreeVector.hh"
54#include "G4AffineTransform.hh"
55#include "G4DynamicParticle.hh"
56#include "G4StepPoint.hh"
57#include "G4VSolid.hh"
58#include "G4LogicalVolume.hh"
59#include "G4SafetyHelper.hh"
61
62///////////////////////////////////////////////////////////////////////////////
63
64
67{
68 lowestEnergy = 1.*keV;
69 fEnvelopeName = anEnvelopeName;
70 fTotXsc = new G4MuNeutrinoNucleusTotXsc();
72 safetyHelper->InitialiseHelper();
73}
74
75///////////////////////////////////////////////////////
76
78{
79 fNuNuclTotXscBias = bf;
80}
81
82///////////////////////////////////////////////////////
83
85{
86 fNuNuclCcBias = bfCc;
87 fNuNuclNcBias = bfNc;
88 fNuNuclTotXscBias = std::max(fNuNuclCcBias, fNuNuclNcBias);
89}
90
91///////////////////////////////////////////////
92
100
101//////////////////////////////////////////////////
102
105{
106 //G4cout << "GetMeanFreePath " << aTrack.GetDefinition()->GetParticleName()
107 // << " Ekin= " << aTrack.GetKineticEnergy() << G4endl;
110 aTrack.GetMaterial());
111
112 if( rName == fEnvelopeName && fNuNuclTotXscBias > 1.)
113 {
114 totxsc *= fNuNuclTotXscBias;
115 }
116 G4double res = (totxsc>0.0) ? 1.0/totxsc : DBL_MAX;
117 //G4cout << " xsection= " << totxsc << G4endl;
118 return res;
119}
120
121///////////////////////////////////////////////////
122
123void G4MuNeutrinoNucleusProcess::ProcessDescription(std::ostream& outFile) const
124{
125 outFile << "G4MuNeutrinoNucleusProcess handles the scattering of \n"
126 << "neutrino on electrons by invoking the following model(s) and \n"
127 << "cross section(s).\n";
128}
129
130///////////////////////////////////////////////////////////////////////
131
134{
136
137 if( rName != fEnvelopeName )
138 {
139 if( verboseLevel > 0 )
140 {
141 G4cout<<"Go out from G4MuNeutrinoNucleusProcess::PostStepDoIt: wrong volume "<<G4endl;
142 }
143 return G4VDiscreteProcess::PostStepDoIt( track, step );
144 }
147 G4double weight = track.GetWeight();
149
150 if( track.GetTrackStatus() != fAlive )
151 {
152 return theTotalResult;
153 }
154 // Next check for illegal track status
155 //
156 if (track.GetTrackStatus() != fAlive &&
157 track.GetTrackStatus() != fSuspend)
158 {
159 if (track.GetTrackStatus() == fStopAndKill ||
162 {
164 ed << "G4HadronicProcess: track in unusable state - "
165 << track.GetTrackStatus() << G4endl;
166 ed << "G4HadronicProcess: returning unchanged track " << G4endl;
167 DumpState(track,"PostStepDoIt",ed);
168 G4Exception("G4HadronicProcess::PostStepDoIt", "had004", JustWarning, ed);
169 }
170 // No warning for fStopButAlive which is a legal status here
171 return theTotalResult;
172 }
173
174 // For elastic scattering, _any_ result is considered an interaction
176
177 G4double kineticEnergy = track.GetKineticEnergy();
178 const G4DynamicParticle* dynParticle = track.GetDynamicParticle();
179 const G4ParticleDefinition* part = dynParticle->GetDefinition();
180 const G4String pName = part->GetParticleName();
181
182 // NOTE: Very low energy scatters were causing numerical (FPE) errors
183 // in earlier releases; these limits have not been changed since.
184
185 if ( kineticEnergy <= lowestEnergy ) return theTotalResult;
186
187 const G4Material* material = track.GetMaterial();
188 G4Nucleus* targNucleus = GetTargetNucleusPointer();
189
190 ///// uniform random spread of the neutrino interaction point ////////////
191
192 const G4StepPoint* pPostStepPoint = step.GetPostStepPoint();
193 const G4DynamicParticle* aParticle = track.GetDynamicParticle();
194 G4ThreeVector position = pPostStepPoint->GetPosition(), newPosition=position;
195 G4ParticleMomentum direction = aParticle->GetMomentumDirection();
196
197 if( fNuNuclCcBias > 1.0 || fNuNuclNcBias > 1.0) // = true, if fBiasingfactor != 1., i.e. xsc is biased
198 {
199 const G4RotationMatrix* rotM = pPostStepPoint->GetTouchable()->GetRotation();
200 G4ThreeVector transl = pPostStepPoint->GetTouchable()->GetTranslation();
201 G4AffineTransform transform = G4AffineTransform(rotM,transl);
202 transform.Invert();
203
204 G4ThreeVector localP = transform.TransformPoint(position);
205 G4ThreeVector localV = transform.TransformAxis(direction);
206
207 G4double forward = track.GetVolume()->GetLogicalVolume()->GetSolid()->DistanceToOut(localP, localV);
208 G4double backward = track.GetVolume()->GetLogicalVolume()->GetSolid()->DistanceToOut(localP, -localV);
209
210 G4double distance = forward+backward;
211
212 // G4cout<<distance/cm<<", ";
213
214 // uniform sampling of nu-e interaction point
215 // along neutrino direction in current volume
216
217 G4double range = -backward+G4UniformRand()*distance;
218
219 newPosition = position + range*direction;
220
221 safetyHelper->ReLocateWithinVolume(newPosition);
222
223 theTotalResult->ProposePosition(newPosition); // G4Exception : GeomNav1002
224 }
225 G4HadProjectile theProj( track );
226 G4HadronicInteraction* hadi = nullptr;
227 G4HadFinalState* result = nullptr;
228 const G4Element* elm =
229 GetCrossSectionDataStore()->SampleZandA(dynParticle, material, *targNucleus);
230 G4int ZZ = elm->GetZasInt();
231 fTotXsc->GetElementCrossSection(dynParticle, ZZ, material);
232 G4double ccTotRatio = fTotXsc->GetCcTotRatio();
233
234 if( G4UniformRand() < ccTotRatio ) // Cc-model
235 {
236 // Initialize the hadronic projectile from the track
237 thePro.Initialise(track);
238
239 if (pName == "nu_mu" ) hadi = (GetHadronicInteractionList())[0];
240 else hadi = (GetHadronicInteractionList())[2];
241
242 result = hadi->ApplyYourself( thePro, *targNucleus);
243
245
247
248 FillResult(result, track);
249 }
250 else // Nc-model
251 {
252 if (pName == "nu_mu" ) hadi = (GetHadronicInteractionList())[1];
253 else hadi = (GetHadronicInteractionList())[3];
254
255 std::size_t idx = track.GetMaterialCutsCouple()->GetIndex();
256
258
259 hadi->SetRecoilEnergyThreshold(tcut);
260
261 if( verboseLevel > 1 )
262 {
263 G4cout << "G4MuNeutrinoNucleusProcess::PostStepDoIt for "
264 << part->GetParticleName()
265 << " in " << material->GetName()
266 << " Target Z= " << targNucleus->GetZ_asInt()
267 << " A= " << targNucleus->GetA_asInt() << G4endl;
268 }
269 try
270 {
271 result = hadi->ApplyYourself( theProj, *targNucleus);
272 }
273 catch(G4HadronicException & aR)
274 {
276 aR.Report(ed);
277 ed << "Call for " << hadi->GetModelName() << G4endl;
278 ed << " Z= "
279 << targNucleus->GetZ_asInt()
280 << " A= " << targNucleus->GetA_asInt() << G4endl;
281 DumpState(track,"ApplyYourself",ed);
282 ed << " ApplyYourself failed" << G4endl;
283 G4Exception("G4MuNeutrinoNucleusProcess::PostStepDoIt", "had006",
284 FatalException, ed);
285 }
286 // directions
287
288 G4ThreeVector indir = track.GetMomentumDirection();
289 G4double phi = CLHEP::twopi*G4UniformRand();
290 G4ThreeVector it(0., 0., 1.);
291 G4ThreeVector outdir = result->GetMomentumChange();
292
293 if(verboseLevel>1)
294 {
295 G4cout << "Efin= " << result->GetEnergyChange()
296 << " de= " << result->GetLocalEnergyDeposit()
297 << " nsec= " << result->GetNumberOfSecondaries()
298 << " dir= " << outdir
299 << G4endl;
300 }
301 // energies
302
303 G4double edep = result->GetLocalEnergyDeposit();
304 G4double efinal = result->GetEnergyChange();
305
306 if(efinal < 0.0) { efinal = 0.0; }
307 if(edep < 0.0) { edep = 0.0; }
308
309 // NOTE: Very low energy scatters were causing numerical (FPE) errors
310 // in earlier releases; these limits have not been changed since.
311
312 if(efinal <= lowestEnergy)
313 {
314 edep += efinal;
315 efinal = 0.0;
316 }
317 // primary change
318
320
321 G4TrackStatus status = track.GetTrackStatus();
322
323 if(efinal > 0.0)
324 {
325 outdir.rotate(phi, it);
326 outdir.rotateUz(indir);
328 }
329 else
330 {
331 if( part->GetProcessManager()->GetAtRestProcessVector()->size() > 0)
332 {
333 status = fStopButAlive;
334 }
335 else
336 {
337 status = fStopAndKill;
338 }
340 }
341 //G4cout << "Efinal= " << efinal << " TrackStatus= " << status << G4endl;
342
344
345 // recoil
346
347 if( result->GetNumberOfSecondaries() > 0 )
348 {
349 G4DynamicParticle* p = result->GetSecondary(0)->GetParticle();
350
351 if(p->GetKineticEnergy() > tcut)
352 {
355
356 // G4cout << "recoil " << pdir << G4endl;
357 //!! is not needed for models inheriting G4MuNeutrinoNucleus
358
359 pdir.rotate(phi, it);
360 pdir.rotateUz(indir);
361
362 // G4cout << "recoil rotated " << pdir << G4endl;
363
364 p->SetMomentumDirection(pdir);
365
366 // in elastic scattering time and weight are not changed
367
368 G4Track* t = new G4Track(p, track.GetGlobalTime(),
369 track.GetPosition());
370 t->SetWeight(weight);
373 }
374 else
375 {
376 edep += p->GetKineticEnergy();
377 delete p;
378 }
379 }
382 result->Clear();
383 }
384 return theTotalResult;
385}
386
387void
389{
390 lowestEnergy = val;
391}
392
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
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep) override
G4double GetMeanFreePath(const G4Track &aTrack, G4double, G4ForceCondition *) override
void SetBiasingFactors(G4double bfCc, G4double bfNc)
void ProcessDescription(std::ostream &outFile) const override
G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
G4MuNeutrinoNucleusProcess(const G4String &anEnvelopeName, const G4String &procName="muNuNucleus")
virtual G4double GetElementCrossSection(const G4DynamicParticle *dynPart, G4int Z, const G4Material *mat)
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
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