Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4NeutrinoElectronProcess.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.18 V.Grichine - PostStepDoIt implementation
34// 03.10.18 V. Grichine - G4Region name and optionally total cross section biased in the region only.
35#include <iostream>
36#include <typeinfo>
37
39#include "G4SystemOfUnits.hh"
40#include "G4Nucleus.hh"
41#include "G4ProcessManager.hh"
47#include "G4VDiscreteProcess.hh"
48
50//#include "G4NeutrinoElectronCcModel.hh"
51//#include "G4NeutrinoElectronNcModel.hh"
52
53#include "G4RotationMatrix.hh"
54#include "G4ThreeVector.hh"
55#include "G4AffineTransform.hh"
56#include "G4DynamicParticle.hh"
57#include "G4StepPoint.hh"
58#include "G4VSolid.hh"
59#include "G4LogicalVolume.hh"
60#include "G4SafetyHelper.hh"
62
63///////////////////////////////////////////////////////////////////////////////
64
65
67 : G4HadronicProcess( pName, fHadronElastic ), isInitialised(false), fBiased(true) // fHadronElastic???
68{
69 lowestEnergy = 1.*keV;
70 fEnvelope = nullptr;
71 fEnvelopeName = anEnvelopeName;
72 fTotXsc = nullptr; // new G4NeutrinoElectronTotXsc();
73 fNuEleCcBias=1.;
74 fNuEleNcBias=1.;
75 fNuEleTotXscBias=1.;
77 safetyHelper->InitialiseHelper();
78}
79
81{
82 // if( fTotXsc ) delete fTotXsc;
83}
84
85///////////////////////////////////////////////////////
86
88{
89 fNuEleTotXscBias = bf;
90
91 fTotXsc = new G4NeutrinoElectronTotXsc();
92 // fTotXsc->SetBiasingFactor(bf);
93}
94
95///////////////////////////////////////////////////////
96
98{
99 fNuEleCcBias=bfCc;
100 fNuEleNcBias=bfNc;
101
102 fTotXsc = new G4NeutrinoElectronTotXsc();
103 fTotXsc->SetBiasingFactors(bfCc, bfNc);
104}
105
106//////////////////////////////////////////////////
107
110{
111 //G4cout << "GetMeanFreePath " << aTrack.GetDefinition()->GetParticleName()
112 // << " Ekin= " << aTrack.GetKineticEnergy() << G4endl;
114 G4double totxsc(0.);
115 try
116 {
117 if( rName == fEnvelopeName && fNuEleTotXscBias > 1.)
118 {
119 totxsc = fNuEleTotXscBias*
121 aTrack.GetMaterial());
122 }
123 else
124 {
126 aTrack.GetMaterial());
127 }
128 }
129 catch(G4HadronicException & aR)
130 {
132 aR.Report(ed);
133 DumpState(aTrack,"GetMeanFreePath",ed);
134 ed << " Cross section is not available" << G4endl;
135 G4Exception("G4NeutrinoElectronProcess::GetMeanFreePath", "had002", FatalException,
136 ed);
137 }
138 G4double res = (totxsc>0.0) ? 1.0/totxsc : DBL_MAX;
139 //G4cout << " xsection= " << totxsc << G4endl;
140 return res;
141}
142
143///////////////////////////////////////////////////
144
145void G4NeutrinoElectronProcess::ProcessDescription(std::ostream& outFile) const
146{
147
148 outFile << "G4NeutrinoElectronProcess handles the scattering of \n"
149 << "neutrino on electrons by invoking the following model(s) and \n"
150 << "cross section(s).\n";
151
152}
153
154///////////////////////////////////////////////////////////////////////
155
158{
159 // track.GetVolume()->GetLogicalVolume()->GetName()
160 // if( track.GetVolume()->GetLogicalVolume() != fEnvelope )
161
163
164 if( rName != fEnvelopeName )
165 {
166 if( verboseLevel > 0 )
167 {
168 G4cout<<"Go out from G4NeutrinoElectronProcess::PostStepDoIt: wrong volume "<<G4endl;
169 }
170 return G4VDiscreteProcess::PostStepDoIt( track, step );
171 }
174 G4double weight = track.GetWeight();
176
177 if( track.GetTrackStatus() != fAlive )
178 {
179 return theTotalResult;
180 }
181 // Next check for illegal track status
182 //
183 if (track.GetTrackStatus() != fAlive &&
184 track.GetTrackStatus() != fSuspend)
185 {
186 if (track.GetTrackStatus() == fStopAndKill ||
189 {
191 ed << "G4HadronicProcess: track in unusable state - "
192 << track.GetTrackStatus() << G4endl;
193 ed << "G4HadronicProcess: returning unchanged track " << G4endl;
194 DumpState(track,"PostStepDoIt",ed);
195 G4Exception("G4HadronicProcess::PostStepDoIt", "had004", JustWarning, ed);
196 }
197 // No warning for fStopButAlive which is a legal status here
198 return theTotalResult;
199 }
200
201 // For elastic scattering, _any_ result is considered an interaction
203
204 G4double kineticEnergy = track.GetKineticEnergy();
205 const G4DynamicParticle* dynParticle = track.GetDynamicParticle();
206 const G4ParticleDefinition* part = dynParticle->GetDefinition();
207
208 // NOTE: Very low energy scatters were causing numerical (FPE) errors
209 // in earlier releases; these limits have not been changed since.
210
211 if ( kineticEnergy <= lowestEnergy ) return theTotalResult;
212
213 const G4Material* material = track.GetMaterial();
214 G4Nucleus* targNucleus = GetTargetNucleusPointer();
215
216 //////////////// uniform random spread of the neutrino interaction point ////////////
217
218 const G4StepPoint* pPostStepPoint = step.GetPostStepPoint();
219 const G4DynamicParticle* aParticle = track.GetDynamicParticle();
220 G4ThreeVector position = pPostStepPoint->GetPosition(), newPosition=position;
221 G4ParticleMomentum direction = aParticle->GetMomentumDirection();
222
223 if( fNuEleCcBias > 1.0 || fNuEleNcBias > 1.0) // = true, if fBiasingfactor != 1., i.e. xsc is biased
224 {
225 const G4RotationMatrix* rotM = pPostStepPoint->GetTouchable()->GetRotation();
226 G4ThreeVector transl = pPostStepPoint->GetTouchable()->GetTranslation();
227 G4AffineTransform transform = G4AffineTransform(rotM,transl);
228 transform.Invert();
229
230 G4ThreeVector localP = transform.TransformPoint(position);
231 G4ThreeVector localV = transform.TransformAxis(direction);
232
233 G4double forward = track.GetVolume()->GetLogicalVolume()->GetSolid()->DistanceToOut(localP, localV);
234 G4double backward = track.GetVolume()->GetLogicalVolume()->GetSolid()->DistanceToOut(localP, -localV);
235
236 G4double distance = forward+backward;
237
238 // G4cout<<distance/cm<<", ";
239
240 // uniform sampling of nu-e interaction point
241 // along neutrino direction in current volume
242
243 G4double range = -backward+G4UniformRand()*distance;
244
245 newPosition = position + range*direction;
246
247 safetyHelper->ReLocateWithinVolume(newPosition);
248
249 theTotalResult->ProposePosition(newPosition); // G4Exception : GeomNav1002
250 }
251 G4HadProjectile theProj( track );
252 G4HadronicInteraction* hadi = nullptr;
253 G4HadFinalState* result = nullptr;
254
255 // Select element
256 const G4Element* elm = nullptr;
257
258 try
259 {
260 elm = GetCrossSectionDataStore()->SampleZandA(dynParticle, material,
261 *targNucleus);
262 }
263 catch( G4HadronicException & aR )
264 {
266 aR.Report(ed);
267 DumpState(track,"SampleZandA",ed);
268 ed << " PostStepDoIt failed on element selection" << G4endl;
269 G4Exception("G4NeutrinoElectronProcess::PostStepDoIt", "had003",
270 FatalException, ed);
271 }
272
273 G4double ccTotRatio = fTotXsc->GetCcRatio();
274
275 if( G4UniformRand() < ccTotRatio ) // Cc-model
276 {
277 // Initialize the hadronic projectile from the track
278 thePro.Initialise(track);
279
280 hadi = (GetHadronicInteractionList())[0];
281
282 result = hadi->ApplyYourself( thePro, *targNucleus);
283
285
287
288 FillResult(result, track);
289 }
290 else // Nc-model, like 'elastic', 2->2 scattering
291 {
292
293 hadi = (GetHadronicInteractionList())[1];
294
295 size_t idx = track.GetMaterialCutsCouple()->GetIndex();
296
298
299 hadi->SetRecoilEnergyThreshold(tcut);
300
301 if( verboseLevel > 1 )
302 {
303 G4cout << "G4NeutrinoElectronProcess::PostStepDoIt for "
304 << part->GetParticleName()
305 << " in " << material->GetName()
306 << " Target Z= " << targNucleus->GetZ_asInt()
307 << " A= " << targNucleus->GetA_asInt() << G4endl;
308 }
309 try
310 {
311 result = hadi->ApplyYourself( theProj, *targNucleus);
312 }
313 catch(G4HadronicException & aR)
314 {
316 aR.Report(ed);
317 ed << "Call for " << hadi->GetModelName() << G4endl;
318 ed << "Target element "<< elm->GetName()<<" Z= "
319 << targNucleus->GetZ_asInt()
320 << " A= " << targNucleus->GetA_asInt() << G4endl;
321 DumpState(track,"ApplyYourself",ed);
322 ed << " ApplyYourself failed" << G4endl;
323 G4Exception("G4NeutrinoElectronProcess::PostStepDoIt", "had006",
324 FatalException, ed);
325 }
326 // directions
327
328 G4ThreeVector indir = track.GetMomentumDirection();
329 G4double phi = CLHEP::twopi*G4UniformRand();
330 G4ThreeVector it(0., 0., 1.);
331 G4ThreeVector outdir = result->GetMomentumChange();
332
333 if(verboseLevel>1)
334 {
335 G4cout << "Efin= " << result->GetEnergyChange()
336 << " de= " << result->GetLocalEnergyDeposit()
337 << " nsec= " << result->GetNumberOfSecondaries()
338 << " dir= " << outdir
339 << G4endl;
340 }
341 // energies
342
343 G4double edep = result->GetLocalEnergyDeposit();
344 G4double efinal = result->GetEnergyChange();
345
346 if(efinal < 0.0) { efinal = 0.0; }
347 if(edep < 0.0) { edep = 0.0; }
348
349 // NOTE: Very low energy scatters were causing numerical (FPE) errors
350 // in earlier releases; these limits have not been changed since.
351
352 if(efinal <= lowestEnergy)
353 {
354 edep += efinal;
355 efinal = 0.0;
356 }
357 // primary change
358
360
361 G4TrackStatus status = track.GetTrackStatus();
362
363 if(efinal > 0.0)
364 {
365 outdir.rotate(phi, it);
366 outdir.rotateUz(indir);
368 }
369 else
370 {
371 if( part->GetProcessManager()->GetAtRestProcessVector()->size() > 0)
372 {
373 status = fStopButAlive;
374 }
375 else
376 {
377 status = fStopAndKill;
378 }
380 }
381 //G4cout << "Efinal= " << efinal << " TrackStatus= " << status << G4endl;
382
384
385 // recoil
386
387 if( result->GetNumberOfSecondaries() > 0 )
388 {
389 G4DynamicParticle* p = result->GetSecondary(0)->GetParticle();
390
391 if(p->GetKineticEnergy() > tcut)
392 {
395
396 // G4cout << "recoil " << pdir << G4endl;
397 //!! is not needed for models inheriting G4NeutrinoElectron
398
399 pdir.rotate(phi, it);
400 pdir.rotateUz(indir);
401
402 // G4cout << "recoil rotated " << pdir << G4endl;
403
404 p->SetMomentumDirection(pdir);
405
406 // in elastic scattering time and weight are not changed
407
408 G4Track* t = new G4Track(p, track.GetGlobalTime(),
409 track.GetPosition());
410 t->SetWeight(weight);
413 }
414 else
415 {
416 edep += p->GetKineticEnergy();
417 delete p;
418 }
419 }
422 result->Clear();
423 }
424 return theTotalResult;
425}
426
427void
429{
430 if(!isInitialised) {
431 isInitialised = true;
432 if(G4Neutron::Neutron() == &part) { lowestEnergy = 1.e-6*eV; }
433 }
435}
436
437void
439{
440 lowestEnergy = val;
441}
442
@ 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
@ fHadronElastic
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 *)
const G4Element * SampleZandA(const G4DynamicParticle *, const G4Material *, G4Nucleus &target)
void SetMomentumDirection(const G4ThreeVector &aDirection)
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
const G4String & GetName() const
Definition: G4Element.hh:127
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
G4double GetMeanFreePath(const G4Track &aTrack, G4double, G4ForceCondition *)
virtual G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep)
G4NeutrinoElectronProcess(G4String anEnvelopeName, const G4String &procName="neutrino-electron")
void SetBiasingFactors(G4double bfCc, G4double bfNc)
virtual void ProcessDescription(std::ostream &outFile) const
virtual void SetLowestEnergy(G4double)
virtual void PreparePhysicsTable(const G4ParticleDefinition &)
void SetBiasingFactors(G4double bfCc, G4double bfNc)
static G4Neutron * Neutron()
Definition: G4Neutron.cc:103
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
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