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