Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/scoring/src/G4ParallelWorldScoringProcess.cc

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

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