Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/scoring/src/G4ParallelWorldProcess.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 ]

  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 "G4ios.hh"
 31 #include "G4ParallelWorldProcess.hh"
 32 #include "G4ParallelWorldProcessStore.hh"
 33 #include "G4Step.hh"
 34 #include "G4StepPoint.hh"
 35 #include "G4Navigator.hh"
 36 #include "G4VTouchable.hh"
 37 #include "G4VPhysicalVolume.hh"
 38 #include "G4ParticleChange.hh"
 39 #include "G4PathFinder.hh"
 40 #include "G4TransportationManager.hh"
 41 #include "G4ParticleChange.hh"
 42 #include "G4StepPoint.hh"
 43 #include "G4FieldTrackUpdator.hh"
 44 #include "G4Material.hh"
 45 #include "G4ProductionCuts.hh"
 46 #include "G4ProductionCutsTable.hh"
 47 
 48 #include "G4SDManager.hh"
 49 #include "G4VSensitiveDetector.hh"
 50 
 51 G4ThreadLocal G4Step* G4ParallelWorldProcess::fpHyperStep = 0;
 52 G4ThreadLocal G4int   G4ParallelWorldProcess::nParallelWorlds = 0;
 53 G4ThreadLocal G4int   G4ParallelWorldProcess::fNavIDHyp = 0;
 54 const G4Step* G4ParallelWorldProcess::GetHyperStep()
 55 { return fpHyperStep; }
 56 G4int G4ParallelWorldProcess::GetHypNavigatorID()
 57 { return fNavIDHyp; }
 58 
 59 G4ParallelWorldProcess::
 60 G4ParallelWorldProcess(const G4String& processName,G4ProcessType theType)
 61 :G4VProcess(processName,theType),fGhostWorld(nullptr),fGhostNavigator(nullptr),
 62  fNavigatorID(-1),fFieldTrack('0'),fGhostSafety(0.),fOnBoundary(false),
 63  layeredMaterialFlag(false)
 64 {
 65   SetProcessSubType(491);
 66   if(!fpHyperStep) fpHyperStep = new G4Step();
 67   iParallelWorld = ++nParallelWorlds;
 68 
 69   pParticleChange = &aDummyParticleChange;
 70 
 71   fGhostStep = new G4Step();
 72   fGhostPreStepPoint = fGhostStep->GetPreStepPoint();
 73   fGhostPostStepPoint = fGhostStep->GetPostStepPoint();
 74 
 75   fTransportationManager = G4TransportationManager::GetTransportationManager();
 76   fTransportationManager->GetNavigatorForTracking()->SetPushVerbosity(false);
 77   fPathFinder = G4PathFinder::GetInstance();
 78 
 79   fGhostWorldName = "** NotDefined **";
 80   G4ParallelWorldProcessStore::GetInstance()->SetParallelWorld(this,processName);
 81 
 82   if (verboseLevel>0)
 83   {
 84     G4cout << GetProcessName() << " is created " << G4endl;
 85   }
 86 }
 87 
 88 G4ParallelWorldProcess::~G4ParallelWorldProcess()
 89 {
 90   delete fGhostStep;
 91   nParallelWorlds--;
 92   if(nParallelWorlds==0)
 93   {
 94     delete fpHyperStep;
 95     fpHyperStep = 0;
 96   }
 97 }
 98 
 99 void G4ParallelWorldProcess::
100 SetParallelWorld(G4String parallelWorldName)
101 {
102   fGhostWorldName = parallelWorldName;
103   fGhostWorld = fTransportationManager->GetParallelWorld(fGhostWorldName);
104   fGhostNavigator = fTransportationManager->GetNavigator(fGhostWorld);
105   fGhostNavigator->SetPushVerbosity(false);
106 }
107 
108 void G4ParallelWorldProcess::
109 SetParallelWorld(G4VPhysicalVolume* parallelWorld)
110 {
111   fGhostWorldName = parallelWorld->GetName();
112   fGhostWorld = parallelWorld;
113   fGhostNavigator = fTransportationManager->GetNavigator(fGhostWorld);
114   fGhostNavigator->SetPushVerbosity(false);
115 }
116 
117 void G4ParallelWorldProcess::StartTracking(G4Track* trk)
118 {
119   if(fGhostNavigator)
120   { fNavigatorID = fTransportationManager->ActivateNavigator(fGhostNavigator); }
121   else
122   {
123     G4Exception("G4ParallelWorldProcess::StartTracking",
124        "ProcParaWorld000",FatalException,
125        "G4ParallelWorldProcess is used for tracking without having a parallel world assigned");
126   }
127   fPathFinder->PrepareNewTrack(trk->GetPosition(),trk->GetMomentumDirection());
128 
129   fOldGhostTouchable = fPathFinder->CreateTouchableHandle(fNavigatorID);
130   fGhostPreStepPoint->SetTouchableHandle(fOldGhostTouchable);
131   fNewGhostTouchable = fOldGhostTouchable;
132   fGhostPostStepPoint->SetTouchableHandle(fNewGhostTouchable);
133 
134   fGhostSafety = -1.;
135   fOnBoundary = false;
136   fGhostPreStepPoint->SetStepStatus(fUndefined);
137   fGhostPostStepPoint->SetStepStatus(fUndefined);
138 
139 //  G4VPhysicalVolume* thePhys = fNewGhostTouchable->GetVolume();
140 //  if(thePhys)
141 //  {
142 //    G4Material* ghostMaterial = thePhys->GetLogicalVolume()->GetMaterial();
143 //    if(ghostMaterial)
144 //    { G4cout << " --- Material : " << ghostMaterial->GetName() << G4endl; }
145 //  }
146 
147   *(fpHyperStep->GetPostStepPoint()) = *(trk->GetStep()->GetPostStepPoint());
148   if(layeredMaterialFlag)
149   {
150     G4StepPoint* realWorldPostStepPoint = trk->GetStep()->GetPostStepPoint();
151     SwitchMaterial(realWorldPostStepPoint);
152     G4StepPoint *realWorldPreStepPoint = trk->GetStep()->GetPreStepPoint();
153     SwitchMaterial(realWorldPreStepPoint);
154     G4double velocity = trk->CalculateVelocity();
155     realWorldPostStepPoint->SetVelocity(velocity);
156     realWorldPreStepPoint->SetVelocity(velocity);
157     trk->SetVelocity(velocity);
158   }
159   *(fpHyperStep->GetPreStepPoint()) = *(fpHyperStep->GetPostStepPoint());
160 }
161 
162 G4double 
163 G4ParallelWorldProcess::AtRestGetPhysicalInteractionLength(
164          const G4Track& /*track*/, 
165          G4ForceCondition* condition)
166 {
167 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
168 // At Rest must be registered ONLY for the particle which has other At Rest
169 // process(es).
170 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
171   *condition = Forced;
172   return DBL_MAX;
173 }
174 
175 G4VParticleChange* G4ParallelWorldProcess::AtRestDoIt(
176      const G4Track& track,
177      const G4Step& step)
178 { 
179 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
180 // At Rest must be registered ONLY for the particle which has other At Rest
181 // process(es).
182 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
183   fOldGhostTouchable = fGhostPostStepPoint->GetTouchableHandle();
184   G4VSensitiveDetector* aSD = 0;
185   if(fOldGhostTouchable->GetVolume())
186   { aSD = fOldGhostTouchable->GetVolume()->GetLogicalVolume()->GetSensitiveDetector(); }
187   fOnBoundary = false;
188   if(aSD)
189   {
190     CopyStep(step);
191     fGhostPreStepPoint->SetSensitiveDetector(aSD);
192 
193     fNewGhostTouchable = fOldGhostTouchable;
194   
195     fGhostPreStepPoint->SetTouchableHandle(fOldGhostTouchable);
196     fGhostPostStepPoint->SetTouchableHandle(fNewGhostTouchable);
197     if(fNewGhostTouchable->GetVolume())
198     {
199       fGhostPostStepPoint->SetSensitiveDetector(
200         fNewGhostTouchable->GetVolume()->GetLogicalVolume()->GetSensitiveDetector());
201     }
202     else
203     { fGhostPostStepPoint->SetSensitiveDetector(0); }
204 
205     aSD->Hit(fGhostStep);
206   }
207 
208   pParticleChange->Initialize(track);
209   return pParticleChange;
210 }
211 
212 G4double 
213 G4ParallelWorldProcess::PostStepGetPhysicalInteractionLength(
214          const G4Track& /*track*/, 
215          G4double   /*previousStepSize*/, 
216          G4ForceCondition* condition)
217 {
218   *condition = StronglyForced;
219   return DBL_MAX;
220 }
221 
222 G4VParticleChange* G4ParallelWorldProcess::PostStepDoIt(
223      const G4Track& track,
224      const G4Step& step)
225 { 
226   fOldGhostTouchable = fGhostPostStepPoint->GetTouchableHandle();
227   G4VSensitiveDetector* aSD = 0;
228   if(fOldGhostTouchable->GetVolume())
229   { aSD = fOldGhostTouchable->GetVolume()->GetLogicalVolume()->GetSensitiveDetector(); }
230   CopyStep(step);
231   fGhostPreStepPoint->SetSensitiveDetector(aSD);
232 
233   if(fOnBoundary)
234   {
235     fNewGhostTouchable = fPathFinder->CreateTouchableHandle(fNavigatorID);
236   }
237   else
238   {
239     fNewGhostTouchable = fOldGhostTouchable;
240   }
241     
242   fGhostPreStepPoint->SetTouchableHandle(fOldGhostTouchable);
243   fGhostPostStepPoint->SetTouchableHandle(fNewGhostTouchable);
244 
245   if(fNewGhostTouchable->GetVolume())
246   {
247     fGhostPostStepPoint->SetSensitiveDetector(
248       fNewGhostTouchable->GetVolume()->GetLogicalVolume()->GetSensitiveDetector());
249   }
250   else
251   { fGhostPostStepPoint->SetSensitiveDetector(0); }
252 
253   G4VSensitiveDetector* sd = fGhostPreStepPoint->GetSensitiveDetector();
254   if(sd)
255   {
256     sd->Hit(fGhostStep);
257   }
258 
259   pParticleChange->Initialize(track); 
260   if(layeredMaterialFlag)
261   {
262     G4StepPoint* realWorldPostStepPoint =
263      ((G4Step*)(track.GetStep()))->GetPostStepPoint();
264     SwitchMaterial(realWorldPostStepPoint);
265   }
266   return pParticleChange;
267 }
268 
269 G4double G4ParallelWorldProcess::AlongStepGetPhysicalInteractionLength(
270             const G4Track& track, G4double  previousStepSize, G4double  currentMinimumStep,
271             G4double& proposedSafety, G4GPILSelection* selection)
272 {
273   static G4ThreadLocal G4FieldTrack *endTrack_G4MT_TLS_ = 0 ; if (!endTrack_G4MT_TLS_) endTrack_G4MT_TLS_ = new  G4FieldTrack ('0') ;  G4FieldTrack &endTrack = *endTrack_G4MT_TLS_;
274   //static ELimited eLimited;
275   ELimited eLimited;
276   ELimited eLim = kUndefLimited;
277   
278   *selection = NotCandidateForSelection;
279   G4double returnedStep = DBL_MAX;
280 
281   if (previousStepSize > 0.)
282   { fGhostSafety -= previousStepSize; }
283   if (fGhostSafety < 0.) fGhostSafety = 0.0;
284       
285   if (currentMinimumStep <= fGhostSafety && currentMinimumStep > 0.)
286   {
287     // I have no chance to limit
288     returnedStep = currentMinimumStep;
289     fOnBoundary = false;
290     proposedSafety = fGhostSafety - currentMinimumStep;
291     eLim = kDoNot;
292   }
293   else 
294   {
295     G4FieldTrackUpdator::Update(&fFieldTrack,&track);
296 
297 #ifdef G4DEBUG_PARALLEL_WORLD_PROCESS    
298     if( verboseLevel > 0 ){
299        int localVerb = verboseLevel-1; 
300 
301        if( localVerb == 1 ) { 
302           G4cout << " Pll Wrl  proc::AlongStepGPIL " << this->GetProcessName() << G4endl;
303        }else if( localVerb > 1 ) { 
304           G4cout << "----------------------------------------------" << G4endl;
305           G4cout << " ParallelWorldProcess: field Track set to : " << G4endl;
306           G4cout << "----------------------------------------------" << G4endl;
307           G4cout << fFieldTrack << G4endl;
308           G4cout << "----------------------------------------------" << G4endl;
309        }
310     }
311 #endif
312     
313     returnedStep
314       = fPathFinder->ComputeStep(fFieldTrack,currentMinimumStep,fNavigatorID,
315                      track.GetCurrentStepNumber(),fGhostSafety,eLimited,
316                      endTrack,track.GetVolume());
317     if(eLimited == kDoNot)
318     {
319       fOnBoundary = false;
320       fGhostSafety = fGhostNavigator->ComputeSafety(endTrack.GetPosition());      
321     }
322     else
323     {
324       fOnBoundary = true;
325       // fGhostSafetyEnd = 0.0;    // At end-point of expected step only
326     }
327     proposedSafety = fGhostSafety;
328     if(eLimited == kUnique || eLimited == kSharedOther) {
329        *selection = CandidateForSelection;
330     }
331     else if (eLimited == kSharedTransport) { 
332        returnedStep *= (1.0 + 1.0e-9);  
333     }
334     eLim = eLimited;
335   }
336 
337   if(iParallelWorld==nParallelWorlds) fNavIDHyp = 0;
338   if(eLim == kUnique || eLim == kSharedOther) fNavIDHyp = fNavigatorID;
339   return returnedStep;
340 }
341 
342 G4VParticleChange* G4ParallelWorldProcess::AlongStepDoIt(
343     const G4Track& track, const G4Step& )
344 {
345   pParticleChange->Initialize(track);
346   return pParticleChange; 
347 }
348 
349 void G4ParallelWorldProcess::CopyStep(const G4Step & step)
350 {
351   G4StepStatus prevStat = fGhostPostStepPoint->GetStepStatus();
352 
353   fGhostStep->SetTrack(step.GetTrack());
354   fGhostStep->SetStepLength(step.GetStepLength());
355   fGhostStep->SetTotalEnergyDeposit(step.GetTotalEnergyDeposit());
356   fGhostStep->SetNonIonizingEnergyDeposit(step.GetNonIonizingEnergyDeposit());
357   fGhostStep->SetControlFlag(step.GetControlFlag());
358   fGhostStep->SetSecondary((const_cast<G4Step&>(step)).GetfSecondary());
359 
360   *fGhostPreStepPoint = *(step.GetPreStepPoint());
361   *fGhostPostStepPoint = *(step.GetPostStepPoint());
362 
363   fGhostPreStepPoint->SetStepStatus(prevStat);
364   if(fOnBoundary)
365   { fGhostPostStepPoint->SetStepStatus(fGeomBoundary); }
366   else if(fGhostPostStepPoint->GetStepStatus()==fGeomBoundary)
367   { fGhostPostStepPoint->SetStepStatus(fPostStepDoItProc); }
368 
369   if(iParallelWorld==1)
370   {
371     G4StepStatus prevStatHyp = fpHyperStep->GetPostStepPoint()->GetStepStatus();
372 
373     fpHyperStep->SetTrack(step.GetTrack());
374     fpHyperStep->SetStepLength(step.GetStepLength());
375     fpHyperStep->SetTotalEnergyDeposit(step.GetTotalEnergyDeposit());
376     fpHyperStep->SetNonIonizingEnergyDeposit(step.GetNonIonizingEnergyDeposit());
377     fpHyperStep->SetControlFlag(step.GetControlFlag());
378 
379     *(fpHyperStep->GetPreStepPoint()) = *(fpHyperStep->GetPostStepPoint());
380     *(fpHyperStep->GetPostStepPoint()) = *(step.GetPostStepPoint());
381   
382     fpHyperStep->GetPreStepPoint()->SetStepStatus(prevStatHyp);
383   }
384 
385   if(fOnBoundary)
386   { fpHyperStep->GetPostStepPoint()->SetStepStatus(fGeomBoundary); }
387 }
388 
389 void G4ParallelWorldProcess::SwitchMaterial(G4StepPoint* realWorldStepPoint)
390 {
391   if(realWorldStepPoint->GetStepStatus()==fWorldBoundary) return;
392   G4VPhysicalVolume* thePhys = fNewGhostTouchable->GetVolume();
393   if(thePhys)
394   {
395     G4Material* ghostMaterial = thePhys->GetLogicalVolume()->GetMaterial();
396     if(ghostMaterial)
397     {
398       G4Region* ghostRegion = thePhys->GetLogicalVolume()->GetRegion();
399       G4ProductionCuts* prodCuts =
400           realWorldStepPoint->GetMaterialCutsCouple()->GetProductionCuts();
401       if(ghostRegion)
402       {
403         G4ProductionCuts* ghostProdCuts = ghostRegion->GetProductionCuts();
404         if(ghostProdCuts) prodCuts = ghostProdCuts;
405       }
406       const G4MaterialCutsCouple* ghostMCCouple =
407           G4ProductionCutsTable::GetProductionCutsTable()
408           ->GetMaterialCutsCouple(ghostMaterial,prodCuts);
409       if(ghostMCCouple)
410       {
411         realWorldStepPoint->SetMaterial(ghostMaterial);
412         realWorldStepPoint->SetMaterialCutsCouple(ghostMCCouple);
413         *(fpHyperStep->GetPostStepPoint()) = *(fGhostPostStepPoint);
414         fpHyperStep->GetPostStepPoint()->SetMaterial(ghostMaterial);
415         fpHyperStep->GetPostStepPoint()->SetMaterialCutsCouple(ghostMCCouple);
416       }
417       else
418       {
419         G4cout << "!!! MaterialCutsCouple is not found for "
420                << ghostMaterial->GetName() << "." << G4endl
421                << "    Material in real world ("
422                << realWorldStepPoint->GetMaterial()->GetName()
423                << ") is used." << G4endl;
424       }
425     }
426   }
427 }
428 
429 G4bool G4ParallelWorldProcess::IsAtRestRequired(G4ParticleDefinition* partDef)
430 {
431   G4int pdgCode = partDef->GetPDGEncoding();
432   if(pdgCode==0)
433   {
434     G4String partName = partDef->GetParticleName();
435     if(partName=="geantino") return false;
436     if(partName=="chargedgeantino") return false;
437   }
438   else
439   {
440     if(pdgCode==11 || pdgCode==2212) return false; // electrons and proton
441     pdgCode = std::abs(pdgCode);
442     if(pdgCode==22) return false; // gamma and optical photons
443     if(pdgCode==12 || pdgCode==14 || pdgCode==16) return false; // all neutronos
444   }
445   return true;
446 }
447 
448