Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ParallelWorldProcess.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
30#include "G4ios.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"
41#include "G4ParticleChange.hh"
42#include "G4StepPoint.hh"
44#include "G4Material.hh"
45#include "G4ProductionCuts.hh"
47
48#include "G4SDManager.hh"
50
51G4ThreadLocal G4Step* G4ParallelWorldProcess::fpHyperStep = 0;
52G4ThreadLocal G4int G4ParallelWorldProcess::nParallelWorlds = 0;
53G4ThreadLocal G4int G4ParallelWorldProcess::fNavIDHyp = 0;
55{ return fpHyperStep; }
57{ return fNavIDHyp; }
58
60G4ParallelWorldProcess(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{
66 if(!fpHyperStep) fpHyperStep = new G4Step();
67 iParallelWorld = ++nParallelWorlds;
68
70
71 fGhostStep = new G4Step();
74
78
79 fGhostWorldName = "** NotDefined **";
81
82 if (verboseLevel>0)
83 {
84 G4cout << GetProcessName() << " is created " << G4endl;
85 }
86}
87
89{
90 delete fGhostStep;
91 nParallelWorlds--;
92 if(nParallelWorlds==0)
93 {
94 delete fpHyperStep;
95 fpHyperStep = 0;
96 }
97}
98
100SetParallelWorld(G4String parallelWorldName)
101{
102 fGhostWorldName = parallelWorldName;
106}
107
110{
111 fGhostWorldName = parallelWorld->GetName();
112 fGhostWorld = parallelWorld;
115}
116
118{
121 else
122 {
123 G4Exception("G4ParallelWorldProcess::StartTracking",
124 "ProcParaWorld000",FatalException,
125 "G4ParallelWorldProcess is used for tracking without having a parallel world assigned");
126 }
128
133
134 fGhostSafety = -1.;
135 fOnBoundary = false;
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());
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
164 const G4Track& /*track*/,
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
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//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
184 G4VSensitiveDetector* aSD = 0;
187 fOnBoundary = false;
188 if(aSD)
189 {
190 CopyStep(step);
192
194
198 {
201 }
202 else
204
205 aSD->Hit(fGhostStep);
206 }
207
209 return pParticleChange;
210}
211
214 const G4Track& /*track*/,
215 G4double /*previousStepSize*/,
217{
219 return DBL_MAX;
220}
221
223 const G4Track& track,
224 const G4Step& step)
225{
227 G4VSensitiveDetector* aSD = 0;
230 CopyStep(step);
232
233 if(fOnBoundary)
234 {
236 }
237 else
238 {
240 }
241
244
246 {
249 }
250 else
252
254 if(sd)
255 {
256 sd->Hit(fGhostStep);
257 }
258
261 {
262 G4StepPoint* realWorldPostStepPoint =
263 ((G4Step*)(track.GetStep()))->GetPostStepPoint();
264 SwitchMaterial(realWorldPostStepPoint);
265 }
266 return pParticleChange;
267}
268
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 {
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;
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;
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
343 const G4Track& track, const G4Step& )
344{
346 return pParticleChange;
347}
348
350{
352
358 fGhostStep->SetSecondary((const_cast<G4Step&>(step)).GetfSecondary());
359
362
364 if(fOnBoundary)
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());
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
390{
391 if(realWorldStepPoint->GetStepStatus()==fWorldBoundary) return;
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 =
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
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
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
@ kUndefLimited
@ kUnique
@ kSharedOther
@ kSharedTransport
G4ProcessType
G4StepStatus
Definition: G4StepStatus.hh:40
@ fGeomBoundary
Definition: G4StepStatus.hh:43
@ fWorldBoundary
Definition: G4StepStatus.hh:41
@ 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
G4Region * GetRegion() const
G4Material * GetMaterial() const
G4ProductionCuts * GetProductionCuts() const
const G4String & GetName() const
Definition: G4Material.hh:175
void SetPushVerbosity(G4bool mode)
virtual G4double ComputeSafety(const G4ThreeVector &globalpoint, const G4double pProposedMaxLength=DBL_MAX, const G4bool keepState=true)
static G4ParallelWorldProcessStore * GetInstance()
void SetParallelWorld(G4ParallelWorldProcess *proc, G4String parallelWorldName)
G4TouchableHandle fNewGhostTouchable
G4VParticleChange * AlongStepDoIt(const G4Track &, const G4Step &)
G4double AtRestGetPhysicalInteractionLength(const G4Track &, G4ForceCondition *)
G4VParticleChange aDummyParticleChange
static const G4Step * GetHyperStep()
G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
G4TransportationManager * fTransportationManager
void CopyStep(const G4Step &step)
G4VPhysicalVolume * fGhostWorld
G4TouchableHandle fOldGhostTouchable
void SwitchMaterial(G4StepPoint *)
G4double AlongStepGetPhysicalInteractionLength(const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
G4ParallelWorldProcess(const G4String &processName="ParaWorld", G4ProcessType theType=fParallel)
void SetParallelWorld(G4String parallelWorldName)
G4double PostStepGetPhysicalInteractionLength(const G4Track &, G4double, G4ForceCondition *)
G4bool IsAtRestRequired(G4ParticleDefinition *)
G4VParticleChange * AtRestDoIt(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
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
static G4ProductionCutsTable * GetProductionCutsTable()
G4ProductionCuts * GetProductionCuts() const
void SetSensitiveDetector(G4VSensitiveDetector *)
void SetMaterial(G4Material *)
G4StepStatus GetStepStatus() const
void SetStepStatus(const G4StepStatus aValue)
void SetVelocity(G4double v)
void SetTouchableHandle(const G4TouchableHandle &apValue)
G4Material * GetMaterial() const
const G4TouchableHandle & GetTouchableHandle() const
void SetMaterialCutsCouple(const G4MaterialCutsCouple *)
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
G4VSensitiveDetector * GetSensitiveDetector() 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 SetSecondary(G4TrackVector *value)
void SetControlFlag(G4SteppingControl StepControlFlag)
void SetTotalEnergyDeposit(G4double value)
G4StepPoint * GetPostStepPoint() const
void SetTrack(G4Track *value)
void SetVelocity(G4double val)
G4VPhysicalVolume * GetVolume() const
const G4ThreeVector & GetPosition() const
G4int GetCurrentStepNumber() const
const G4ThreeVector & GetMomentumDirection() const
G4double CalculateVelocity() const
const G4Step * GetStep() const
G4VPhysicalVolume * GetParallelWorld(const G4String &worldName)
static G4TransportationManager * GetTransportationManager()
G4Navigator * GetNavigatorForTracking() const
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
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:406
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
#define DBL_MAX
Definition: templates.hh:62
#define G4ThreadLocal
Definition: tls.hh:77