Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ElNeutrinoNucleusProcess.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"
42#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 G4ElNeutrinoNucleusTotXsc();
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;
109 G4double totxsc =
111 aTrack.GetMaterial());
112
113 if ( rName == fEnvelopeName )
114 {
115 totxsc *= fNuNuclTotXscBias;
116 }
117 G4double res = (totxsc>0.0) ? 1.0/totxsc : DBL_MAX;
118 //G4cout << " xsection= " << totxsc << G4endl;
119 return res;
120}
121
122///////////////////////////////////////////////////
123
124void G4ElNeutrinoNucleusProcess::ProcessDescription(std::ostream& outFile) const
125{
126 outFile << "G4ElNeutrinoNucleusProcess handles the scattering of \n"
127 << "neutrino on electrons by invoking the following model(s) and \n"
128 << "cross section(s).\n";
129
130}
131
132///////////////////////////////////////////////////////////////////////
133
136{
137 // track.GetVolume()->GetLogicalVolume()->GetName()
138 // if( track.GetVolume()->GetLogicalVolume() != fEnvelope )
139
141
142 if( rName != fEnvelopeName )
143 {
144 if( verboseLevel > 0 )
145 {
146 G4cout<<"Go out from G4ElNeutrinoNucleusProcess::PostStepDoIt: wrong volume "<<G4endl;
147 }
148 return G4VDiscreteProcess::PostStepDoIt( track, step );
149 }
152 G4double weight = track.GetWeight();
154
155 if( track.GetTrackStatus() != fAlive )
156 {
157 return theTotalResult;
158 }
159 // Next check for illegal track status
160 //
161 if (track.GetTrackStatus() != fAlive &&
162 track.GetTrackStatus() != fSuspend)
163 {
164 if (track.GetTrackStatus() == fStopAndKill ||
167 {
169 ed << "G4HadronicProcess: track in unusable state - "
170 << track.GetTrackStatus() << G4endl;
171 ed << "G4HadronicProcess: returning unchanged track " << G4endl;
172 DumpState(track,"PostStepDoIt",ed);
173 G4Exception("G4HadronicProcess::PostStepDoIt", "had004", JustWarning, ed);
174 }
175 // No warning for fStopButAlive which is a legal status here
176 return theTotalResult;
177 }
178
179 // For elastic scattering, _any_ result is considered an interaction
181
182 G4double kineticEnergy = track.GetKineticEnergy();
183 const G4DynamicParticle* dynParticle = track.GetDynamicParticle();
184 const G4ParticleDefinition* part = dynParticle->GetDefinition();
185 const G4String pName = part->GetParticleName();
186
187 // NOTE: Very low energy scatters were causing numerical (FPE) errors
188 // in earlier releases; these limits have not been changed since.
189
190 if ( kineticEnergy <= lowestEnergy ) return theTotalResult;
191
192 const G4Material* material = track.GetMaterial();
193 G4Nucleus* targNucleus = GetTargetNucleusPointer();
194
195 //////////////// uniform random spread of the neutrino interaction point ////////////
196
197 const G4StepPoint* pPostStepPoint = step.GetPostStepPoint();
198 const G4DynamicParticle* aParticle = track.GetDynamicParticle();
199 G4ThreeVector position = pPostStepPoint->GetPosition(), newPosition=position;
200 G4ParticleMomentum direction = aParticle->GetMomentumDirection();
201
202 if( fNuNuclCcBias > 1.0 || fNuNuclNcBias > 1.0) // = true, if fBiasingfactor != 1., i.e. xsc is biased
203 {
204 const G4RotationMatrix* rotM = pPostStepPoint->GetTouchable()->GetRotation();
205 G4ThreeVector transl = pPostStepPoint->GetTouchable()->GetTranslation();
206 G4AffineTransform transform = G4AffineTransform(rotM,transl);
207 transform.Invert();
208
209 G4ThreeVector localP = transform.TransformPoint(position);
210 G4ThreeVector localV = transform.TransformAxis(direction);
211
212 G4double forward = track.GetVolume()->GetLogicalVolume()->GetSolid()->DistanceToOut(localP, localV);
213 G4double backward = track.GetVolume()->GetLogicalVolume()->GetSolid()->DistanceToOut(localP, -localV);
214
215 G4double distance = forward+backward;
216
217 // G4cout<<distance/cm<<", ";
218
219 // uniform sampling of nu-e interaction point
220 // along neutrino direction in current volume
221
222 G4double range = -backward+G4UniformRand()*distance;
223
224 newPosition = position + range*direction;
225
226 safetyHelper->ReLocateWithinVolume(newPosition);
227
228 theTotalResult->ProposePosition(newPosition); // G4Exception : GeomNav1002
229 }
230 G4HadProjectile theProj( track );
231 G4HadronicInteraction* hadi = nullptr;
232 G4HadFinalState* result = nullptr;
233 const G4Element* elm =
234 GetCrossSectionDataStore()->SampleZandA(dynParticle, material,
235 *targNucleus);
236 G4int ZZ = elm->GetZasInt();
237 fTotXsc->GetElementCrossSection(dynParticle, ZZ, material);
238 G4double ccTotRatio = fTotXsc->GetCcTotRatio();
239
240 if( G4UniformRand() < ccTotRatio ) // Cc-model
241 {
242 // Initialize the hadronic projectile from the track
243 thePro.Initialise(track);
244
245 if (pName == "nu_e" ) hadi = (GetHadronicInteractionList())[0];
246 else hadi = (GetHadronicInteractionList())[2];
247
248 result = hadi->ApplyYourself( thePro, *targNucleus);
249
251
253
254 FillResult(result, track);
255 }
256 else // Nc-model
257 {
258
259 if (pName == "nu_e" ) hadi = (GetHadronicInteractionList())[1];
260 else hadi = (GetHadronicInteractionList())[3];
261
262 std::size_t idx = track.GetMaterialCutsCouple()->GetIndex();
263
265
266 hadi->SetRecoilEnergyThreshold(tcut);
267
268 if( verboseLevel > 1 )
269 {
270 G4cout << "G4ElNeutrinoNucleusProcess::PostStepDoIt for "
271 << part->GetParticleName()
272 << " in " << material->GetName()
273 << " Target Z= " << targNucleus->GetZ_asInt()
274 << " A= " << targNucleus->GetA_asInt() << G4endl;
275 }
276 try
277 {
278 result = hadi->ApplyYourself( theProj, *targNucleus);
279 }
280 catch(G4HadronicException & aR)
281 {
283 aR.Report(ed);
284 ed << "Call for " << hadi->GetModelName() << G4endl;
285 ed << " Z= "
286 << targNucleus->GetZ_asInt()
287 << " A= " << targNucleus->GetA_asInt() << G4endl;
288 DumpState(track,"ApplyYourself",ed);
289 ed << " ApplyYourself failed" << G4endl;
290 G4Exception("G4ElNeutrinoNucleusProcess::PostStepDoIt", "had006",
291 FatalException, ed);
292 }
293 // directions
294
295 G4ThreeVector indir = track.GetMomentumDirection();
296 G4double phi = CLHEP::twopi*G4UniformRand();
297 G4ThreeVector it(0., 0., 1.);
298 G4ThreeVector outdir = result->GetMomentumChange();
299
300 if(verboseLevel>1)
301 {
302 G4cout << "Efin= " << result->GetEnergyChange()
303 << " de= " << result->GetLocalEnergyDeposit()
304 << " nsec= " << result->GetNumberOfSecondaries()
305 << " dir= " << outdir
306 << G4endl;
307 }
308 // energies
309
310 G4double edep = result->GetLocalEnergyDeposit();
311 G4double efinal = result->GetEnergyChange();
312
313 if(efinal < 0.0) { efinal = 0.0; }
314 if(edep < 0.0) { edep = 0.0; }
315
316 // NOTE: Very low energy scatters were causing numerical (FPE) errors
317 // in earlier releases; these limits have not been changed since.
318
319 if(efinal <= lowestEnergy)
320 {
321 edep += efinal;
322 efinal = 0.0;
323 }
324 // primary change
325
327
328 G4TrackStatus status = track.GetTrackStatus();
329
330 if(efinal > 0.0)
331 {
332 outdir.rotate(phi, it);
333 outdir.rotateUz(indir);
335 }
336 else
337 {
338 if( part->GetProcessManager()->GetAtRestProcessVector()->size() > 0)
339 {
340 status = fStopButAlive;
341 }
342 else
343 {
344 status = fStopAndKill;
345 }
347 }
348 //G4cout << "Efinal= " << efinal << " TrackStatus= " << status << G4endl;
349
351
352 // recoil
353
354 if( result->GetNumberOfSecondaries() > 0 )
355 {
356 G4DynamicParticle* p = result->GetSecondary(0)->GetParticle();
357
358 if(p->GetKineticEnergy() > tcut)
359 {
362
363 // G4cout << "recoil " << pdir << G4endl;
364 //!! is not needed for models inheriting G4ElNeutrinoNucleus
365
366 pdir.rotate(phi, it);
367 pdir.rotateUz(indir);
368
369 // G4cout << "recoil rotated " << pdir << G4endl;
370
371 p->SetMomentumDirection(pdir);
372
373 // in elastic scattering time and weight are not changed
374
375 G4Track* t = new G4Track(p, track.GetGlobalTime(),
376 track.GetPosition());
377 t->SetWeight(weight);
380 }
381 else
382 {
383 edep += p->GetKineticEnergy();
384 delete p;
385 }
386 }
389 result->Clear();
390 }
391 return theTotalResult;
392}
393
394void
396{
397 lowestEnergy = val;
398}
399
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
G4ElNeutrinoNucleusProcess(const G4String &anEnvelopeName, const G4String &procName="elNuNucleus")
G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep) override
void SetBiasingFactors(G4double bfCc, G4double bfNc)
void ProcessDescription(std::ostream &outFile) const override
G4double GetMeanFreePath(const G4Track &aTrack, G4double, G4ForceCondition *) override
G4double GetElementCrossSection(const G4DynamicParticle *dynPart, G4int Z, const G4Material *mat) override
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
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