Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ParallelWorldScoringProcess.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//
28//
29
31
32#include "G4ios.hh"
33#include "G4SystemOfUnits.hh"
34#include "G4Step.hh"
35#include "G4Navigator.hh"
36#include "G4VTouchable.hh"
37#include "G4VPhysicalVolume.hh"
38#include "G4ParticleChange.hh"
39#include "G4PathFinder.hh"
41#include "G4ParticleChange.hh"
42#include "G4StepPoint.hh"
45
46#include "G4SDManager.hh"
48
49//--------------------------------
50// Constructor with name and type:
51//--------------------------------
54:G4VProcess(processName,theType), fGhostNavigator(0), fNavigatorID(-1), fFieldTrack('0')
55{
56 pParticleChange = &aDummyParticleChange;
57
58 fGhostStep = new G4Step();
59 fGhostPreStepPoint = fGhostStep->GetPreStepPoint();
60 fGhostPostStepPoint = fGhostStep->GetPostStepPoint();
61
63 fPathFinder = G4PathFinder::GetInstance();
64
65 fGhostWorld = 0;
66 fGhostSafety = 0.;
67 fOnBoundary = false;
68
69 if (verboseLevel>0)
70 {
71 G4cout << GetProcessName() << " is created " << G4endl;
72 }
73}
74
75// -----------
76// Destructor:
77// -----------
79{
80 delete fGhostStep;
81}
82
83//------------------------------------------------------
84//
85// SetParallelWorld
86//
87//------------------------------------------------------
89SetParallelWorld(G4String parallelWorldName)
90{
91//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
92// Get pointers of the parallel world and its navigator
93//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
94 fGhostWorldName = parallelWorldName;
95 fGhostWorld = fTransportationManager->GetParallelWorld(fGhostWorldName);
96 fGhostNavigator = fTransportationManager->GetNavigator(fGhostWorld);
97//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
98}
99
102{
103//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
104// Get pointer of navigator
105//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
106 fGhostWorldName = parallelWorld->GetName();
107 fGhostWorld = parallelWorld;
108 fGhostNavigator = fTransportationManager->GetNavigator(fGhostWorld);
109//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
110}
111
114{
115 G4int pdgCode = partDef->GetPDGEncoding();
116 if(pdgCode==0)
117 {
118 G4String partName = partDef->GetParticleName();
119 if(partName=="geantino") return false;
120 if(partName=="chargedgeantino") return false;
121 }
122 else
123 {
124 if(pdgCode==11 || pdgCode==2212) return false; // electrons and proton
125 pdgCode = std::abs(pdgCode);
126 if(pdgCode==22) return false; // gamma and optical photons
127 if(pdgCode==12 || pdgCode==14 || pdgCode==16) return false; // all neutronos
128 }
129 return true;
130}
131
132
133//------------------------------------------------------
134//
135// StartTracking
136//
137//------------------------------------------------------
139{
140//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
141// Activate navigator and get the navigator ID
142//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
143// G4cout << " G4ParallelWorldScoringProcess::StartTracking" << G4endl;
144 if(fGhostNavigator)
145 { fNavigatorID = fTransportationManager->ActivateNavigator(fGhostNavigator); }
146 else
147 {
148 G4Exception("G4ParallelWorldScoringProcess::StartTracking",
149 "ProcParaWorld000",FatalException,
150 "G4ParallelWorldScoringProcess is used for tracking without having a parallel world assigned");
151 }
152//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
153
154// G4cout << "G4ParallelWorldScoringProcess::StartTracking <<<<<<<<<<<<<<<<<< " << G4endl;
155//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
156// Let PathFinder initialize
157//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
158 fPathFinder->PrepareNewTrack(trk->GetPosition(),trk->GetMomentumDirection());
159//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
160
161//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
162// Setup initial touchables for the first step
163//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
164 fOldGhostTouchable = fPathFinder->CreateTouchableHandle(fNavigatorID);
165 fGhostPreStepPoint->SetTouchableHandle(fOldGhostTouchable);
166 fNewGhostTouchable = fOldGhostTouchable;
167 fGhostPostStepPoint->SetTouchableHandle(fNewGhostTouchable);
168
169 // Initialize
170 fGhostSafety = -1.;
171 fOnBoundary = false;
172 fGhostPreStepPoint->SetStepStatus(fUndefined);
173 fGhostPostStepPoint->SetStepStatus(fUndefined);
174}
175
176//----------------------------------------------------------
177//
178// AtRestGetPhysicalInteractionLength()
179//
180//----------------------------------------------------------
183 const G4Track& /*track*/,
185{
186 *condition = Forced;
187 return DBL_MAX;
188}
189
190//------------------------------------
191//
192// AtRestDoIt()
193//
194//------------------------------------
196 const G4Track& track,
197 const G4Step& step)
198{
199 fOldGhostTouchable = fGhostPostStepPoint->GetTouchableHandle();
200 G4VSensitiveDetector* aSD = 0;
201 if(fOldGhostTouchable->GetVolume())
202 { aSD = fOldGhostTouchable->GetVolume()->GetLogicalVolume()->GetSensitiveDetector(); }
203 fOnBoundary = false;
204 CopyStep(step);
205 fGhostPreStepPoint->SetSensitiveDetector(aSD);
206
207 fNewGhostTouchable = fOldGhostTouchable;
208
209 fGhostPreStepPoint->SetTouchableHandle(fOldGhostTouchable);
210 fGhostPostStepPoint->SetTouchableHandle(fNewGhostTouchable);
211 if(fNewGhostTouchable->GetVolume())
212 {
213 fGhostPostStepPoint->SetSensitiveDetector(
214 fNewGhostTouchable->GetVolume()->GetLogicalVolume()->GetSensitiveDetector());
215 }
216 else
217 { fGhostPostStepPoint->SetSensitiveDetector(0); }
218
219 if (verboseLevel>1) Verbose(step);
220
221 G4VSensitiveDetector* sd = fGhostPreStepPoint->GetSensitiveDetector();
222 if(sd)
223 {
224 sd->Hit(fGhostStep);
225 }
226
228 return pParticleChange;
229}
230
231//----------------------------------------------------------
232//
233// PostStepGetPhysicalInteractionLength()
234//
235//----------------------------------------------------------
238 const G4Track& /*track*/,
239 G4double /*previousStepSize*/,
241{
242 // I must be invoked anyway to score the hit.
244 return DBL_MAX;
245}
246
247//------------------------------------
248//
249// PostStepDoIt()
250//
251//------------------------------------
253 const G4Track& track,
254 const G4Step& step)
255{
256 fOldGhostTouchable = fGhostPostStepPoint->GetTouchableHandle();
257 G4VSensitiveDetector* aSD = 0;
258 if(fOldGhostTouchable->GetVolume())
259 { aSD = fOldGhostTouchable->GetVolume()->GetLogicalVolume()->GetSensitiveDetector(); }
260 CopyStep(step);
261 fGhostPreStepPoint->SetSensitiveDetector(aSD);
262
263 // fPathFinder->Locate( track.GetPosition(),
264 // track.GetMomentumDirection(),
265 // true);
266
267 // fPathFinder->Locate(step.GetPostStepPoint()->GetPosition(),
268 // step.GetPostStepPoint()->GetMomentumDirection());
269
270 if(fOnBoundary)
271 {
272//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
273// Locate the point and get new touchable
274//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
275 //?? fPathFinder->Locate(step.GetPostStepPoint()->GetPosition(),
276 //?? step.GetPostStepPoint()->GetMomentumDirection());
277 fNewGhostTouchable = fPathFinder->CreateTouchableHandle(fNavigatorID);
278//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
279 }
280 else
281 {
282// Do I need this ??????????????????????????????????????????????????????????
283// fGhostNavigator->LocateGlobalPointWithinVolume(track.GetPosition());
284// ?????????????????????????????????????????????????????????????????????????
285
286 // fPathFinder->ReLocate(track.GetPosition());
287
288 // reuse the touchable
289 fNewGhostTouchable = fOldGhostTouchable;
290 }
291
292 fGhostPreStepPoint->SetTouchableHandle(fOldGhostTouchable);
293 fGhostPostStepPoint->SetTouchableHandle(fNewGhostTouchable);
294
295 if(fNewGhostTouchable->GetVolume())
296 {
297 fGhostPostStepPoint->SetSensitiveDetector(
298 fNewGhostTouchable->GetVolume()->GetLogicalVolume()->GetSensitiveDetector());
299 }
300 else
301 { fGhostPostStepPoint->SetSensitiveDetector(0); }
302
303 if (verboseLevel>1) Verbose(step);
304
305 G4VSensitiveDetector* sd = fGhostPreStepPoint->GetSensitiveDetector();
306 if(sd)
307 {
308 sd->Hit(fGhostStep);
309 }
310
311 pParticleChange->Initialize(track); // Does not change the track properties
312 return pParticleChange;
313}
314
315
316//---------------------------------------
317//
318// AlongStepGetPhysicalInteractionLength
319//
320//---------------------------------------
322 const G4Track& track, G4double previousStepSize, G4double currentMinimumStep,
323 G4double& proposedSafety, G4GPILSelection* selection)
324{
325 static G4ThreadLocal G4FieldTrack *endTrack_G4MT_TLS_ = 0 ; if (!endTrack_G4MT_TLS_) endTrack_G4MT_TLS_ = new G4FieldTrack ('0') ; G4FieldTrack &endTrack = *endTrack_G4MT_TLS_;
326 static G4ThreadLocal ELimited *eLimited_G4MT_TLS_ = 0 ; if (!eLimited_G4MT_TLS_) eLimited_G4MT_TLS_ = new ELimited ; ELimited &eLimited = *eLimited_G4MT_TLS_;
327
328 *selection = NotCandidateForSelection;
329 G4double returnedStep = DBL_MAX;
330
331 if (previousStepSize > 0.)
332 { fGhostSafety -= previousStepSize; }
333// else
334// { fGhostSafety = -1.; }
335 if (fGhostSafety < 0.) fGhostSafety = 0.0;
336
337 // ------------------------------------------
338 // Determination of the proposed STEP LENGTH:
339 // ------------------------------------------
340 if (currentMinimumStep <= fGhostSafety && currentMinimumStep > 0.)
341 {
342 // I have no chance to limit
343 returnedStep = currentMinimumStep;
344 fOnBoundary = false;
345 proposedSafety = fGhostSafety - currentMinimumStep;
346 }
347 else // (currentMinimumStep > fGhostSafety: I may limit the Step)
348 {
349 G4FieldTrackUpdator::Update(&fFieldTrack,&track);
350//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
351// ComputeStep
352//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
353 returnedStep
354 = fPathFinder->ComputeStep(fFieldTrack,currentMinimumStep,fNavigatorID,
355 track.GetCurrentStepNumber(),fGhostSafety,eLimited,
356 endTrack,track.GetVolume());
357//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
358 if(eLimited == kDoNot)
359 {
360 // Track is not on the boundary
361 fOnBoundary = false;
362 fGhostSafety = fGhostNavigator->ComputeSafety(endTrack.GetPosition());
363 }
364 else
365 {
366 // Track is on the boundary
367 fOnBoundary = true;
368 // fEndGhostSafety = 0.0; // Will apply at the end of the step ...
369 }
370 proposedSafety = fGhostSafety;
371 if(eLimited == kUnique || eLimited == kSharedOther) {
372 *selection = CandidateForSelection;
373 }else if (eLimited == kSharedTransport) {
374 returnedStep *= (1.0 + 1.0e-9);
375 // Expand to disable its selection in Step Manager comparison
376 }
377 }
378
379 // ----------------------------------------------
380 // Returns the fGhostSafety as the proposedSafety
381 // The SteppingManager will take care of keeping
382 // the smallest one.
383 // ----------------------------------------------
384 return returnedStep;
385}
386
388 const G4Track& track, const G4Step& )
389{
390 // Dummy ParticleChange ie: does nothing
391 // Expecting G4Transportation to move the track
393 return pParticleChange;
394}
395
396
397void G4ParallelWorldScoringProcess::CopyStep(const G4Step & step)
398{
399 G4StepStatus prevStat = fGhostPostStepPoint->GetStepStatus();
400
401 fGhostStep->SetTrack(step.GetTrack());
402 fGhostStep->SetStepLength(step.GetStepLength());
405 fGhostStep->SetControlFlag(step.GetControlFlag());
406
407 *fGhostPreStepPoint = *(step.GetPreStepPoint());
408 *fGhostPostStepPoint = *(step.GetPostStepPoint());
409
410//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
411// Set StepStatus for ghost world
412//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
413 fGhostPreStepPoint->SetStepStatus(prevStat);
414 if(fOnBoundary)
415 { fGhostPostStepPoint->SetStepStatus(fGeomBoundary); }
416 else if(fGhostPostStepPoint->GetStepStatus()==fGeomBoundary)
417 { fGhostPostStepPoint->SetStepStatus(fPostStepDoItProc); }
418//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
419}
420
422{
423 G4cout << "In mass geometry ------------------------------------------------" << G4endl;
424 G4cout << " StepLength : " << step.GetStepLength()/mm << " TotalEnergyDeposit : "
425 << step.GetTotalEnergyDeposit()/MeV << G4endl;
426 G4cout << " PreStepPoint : "
427 << step.GetPreStepPoint()->GetPhysicalVolume()->GetName() << " - ";
430 else
431 { G4cout << "NoProcessAssigned"; }
432 G4cout << G4endl;
433 G4cout << " " << step.GetPreStepPoint()->GetPosition() << G4endl;
434 G4cout << " PostStepPoint : ";
437 else
438 { G4cout << "OutOfWorld"; }
439 G4cout << " - ";
442 else
443 { G4cout << "NoProcessAssigned"; }
444 G4cout << G4endl;
445 G4cout << " " << step.GetPostStepPoint()->GetPosition() << G4endl;
446
447 G4cout << "In ghost geometry ------------------------------------------------" << G4endl;
448 G4cout << " StepLength : " << fGhostStep->GetStepLength()/mm
449 << " TotalEnergyDeposit : "
450 << fGhostStep->GetTotalEnergyDeposit()/MeV << G4endl;
451 G4cout << " PreStepPoint : "
452 << fGhostStep->GetPreStepPoint()->GetPhysicalVolume()->GetName() << " ["
453 << fGhostStep->GetPreStepPoint()->GetTouchable()->GetReplicaNumber()
454 << " ]" << " - ";
455 if(fGhostStep->GetPreStepPoint()->GetProcessDefinedStep())
457 else
458 { G4cout << "NoProcessAssigned"; }
459 G4cout << G4endl;
460 G4cout << " " << fGhostStep->GetPreStepPoint()->GetPosition() << G4endl;
461 G4cout << " PostStepPoint : ";
462 if(fGhostStep->GetPostStepPoint()->GetPhysicalVolume())
463 {
464 G4cout << fGhostStep->GetPostStepPoint()->GetPhysicalVolume()->GetName() << " ["
465 << fGhostStep->GetPostStepPoint()->GetTouchable()->GetReplicaNumber()
466 << " ]";
467 }
468 else
469 { G4cout << "OutOfWorld"; }
470 G4cout << " - ";
471 if(fGhostStep->GetPostStepPoint()->GetProcessDefinedStep())
473 else
474 { G4cout << "NoProcessAssigned"; }
475 G4cout << G4endl;
476 G4cout << " " << fGhostStep->GetPostStepPoint()->GetPosition() << " == "
477 << fGhostStep->GetTrack()->GetMomentumDirection()
478 << G4endl;
479
480}
481
G4double condition(const G4ErrorSymMatrix &m)
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
G4ForceCondition
@ StronglyForced
@ Forced
G4GPILSelection
@ CandidateForSelection
@ NotCandidateForSelection
ELimited
@ kDoNot
@ kUnique
@ kSharedOther
@ kSharedTransport
G4ProcessType
G4StepStatus
Definition: G4StepStatus.hh:40
@ fGeomBoundary
Definition: G4StepStatus.hh:43
@ fUndefined
Definition: G4StepStatus.hh:55
@ fPostStepDoItProc
Definition: G4StepStatus.hh:49
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
static void Update(G4FieldTrack *, const G4Track *)
G4ThreeVector GetPosition() const
G4VSensitiveDetector * GetSensitiveDetector() const
virtual G4double ComputeSafety(const G4ThreeVector &globalpoint, const G4double pProposedMaxLength=DBL_MAX, const G4bool keepState=true)
G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
G4double AtRestGetPhysicalInteractionLength(const G4Track &, G4ForceCondition *)
G4bool IsAtRestRequired(G4ParticleDefinition *partDef)
void SetParallelWorld(G4String parallelWorldName)
G4double AlongStepGetPhysicalInteractionLength(const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
G4ParallelWorldScoringProcess(const G4String &processName="ParaWorldScore", G4ProcessType theType=fParameterisation)
G4VParticleChange * AtRestDoIt(const G4Track &, const G4Step &)
G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
G4VParticleChange * AlongStepDoIt(const G4Track &, const G4Step &)
const G4String & GetParticleName() const
G4double ComputeStep(const G4FieldTrack &pFieldTrack, G4double pCurrentProposedStepLength, G4int navigatorId, G4int stepNo, G4double &pNewSafety, ELimited &limitedStep, G4FieldTrack &EndState, G4VPhysicalVolume *currentVolume)
void PrepareNewTrack(const G4ThreeVector &position, const G4ThreeVector &direction, G4VPhysicalVolume *massStartVol=nullptr)
static G4PathFinder * GetInstance()
Definition: G4PathFinder.cc:52
G4TouchableHandle CreateTouchableHandle(G4int navId) const
void SetSensitiveDetector(G4VSensitiveDetector *)
G4StepStatus GetStepStatus() const
const G4VTouchable * GetTouchable() const
void SetStepStatus(const G4StepStatus aValue)
const G4VProcess * GetProcessDefinedStep() const
void SetTouchableHandle(const G4TouchableHandle &apValue)
const G4ThreeVector & GetPosition() const
const G4TouchableHandle & GetTouchableHandle() const
G4VSensitiveDetector * GetSensitiveDetector() const
G4VPhysicalVolume * GetPhysicalVolume() const
Definition: G4Step.hh:62
G4SteppingControl GetControlFlag() const
G4Track * GetTrack() const
void SetStepLength(G4double value)
void SetNonIonizingEnergyDeposit(G4double value)
G4StepPoint * GetPreStepPoint() const
G4double GetNonIonizingEnergyDeposit() const
G4double GetStepLength() const
G4double GetTotalEnergyDeposit() const
void SetControlFlag(G4SteppingControl StepControlFlag)
void SetTotalEnergyDeposit(G4double value)
G4StepPoint * GetPostStepPoint() const
void SetTrack(G4Track *value)
G4VPhysicalVolume * GetVolume() const
const G4ThreeVector & GetPosition() const
G4int GetCurrentStepNumber() const
const G4ThreeVector & GetMomentumDirection() const
G4VPhysicalVolume * GetParallelWorld(const G4String &worldName)
static G4TransportationManager * GetTransportationManager()
G4int ActivateNavigator(G4Navigator *aNavigator)
G4Navigator * GetNavigator(const G4String &worldName)
virtual void Initialize(const G4Track &)
G4LogicalVolume * GetLogicalVolume() const
const G4String & GetName() const
G4int verboseLevel
Definition: G4VProcess.hh:356
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:321
const G4String & GetProcessName() const
Definition: G4VProcess.hh:382
G4bool Hit(G4Step *aStep)
virtual G4VPhysicalVolume * GetVolume(G4int depth=0) const
Definition: G4VTouchable.cc:41
virtual G4int GetReplicaNumber(G4int depth=0) const
Definition: G4VTouchable.cc:55
#define DBL_MAX
Definition: templates.hh:62
#define G4ThreadLocal
Definition: tls.hh:77