CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
GflashHadronWrapperProcess.cc
Go to the documentation of this file.
1 //
2 // S.Y. Jun, August 2007
3 //
5 
6 #include "G4Track.hh"
7 #include "G4VParticleChange.hh"
8 #include "G4ProcessManager.hh"
9 #include "G4ProcessVector.hh"
10 #include "G4VProcess.hh"
11 #include "G4GPILSelection.hh"
12 
13 using namespace CLHEP;
14 
16  particleChange(0),
17  pmanager(0),
18  fProcessVector(0),
19  fProcess(0)
20 {
21  theProcessName = processName;
22 }
23 
25 }
26 
27 G4VParticleChange* GflashHadronWrapperProcess::PostStepDoIt(const G4Track& track, const G4Step& step)
28 {
29  // process PostStepDoIt for the original process
30 
31  particleChange = pRegProcess->PostStepDoIt(track, step);
32 
33  // specific actions of the wrapper process
34 
35  // update step/track information after PostStep of the original process is done
36  // these steps will be repeated again without additional conflicts even if the
37  // parameterized physics process doesn't take over the original process
38 
39  particleChange->UpdateStepForPostStep(const_cast<G4Step *> (&step));
40 
41  // we may not want to update G4Track during this wrapper process if
42  // G4Track accessors are not used in the G4VFastSimulationModel model
43  // (const_cast<G4Step *> (&step))->UpdateTrack();
44 
45  // update safety after each invocation of PostStepDoIts
46  // G4double safety = std::max(endpointSafety - (endpointSafOrigin - fPostStepPoint->GetPosition()).mag(),0.);
47  // step.GetPostStepPoint()->SetSafety(safety);
48 
49  // the secondaries from the original process are also created at the wrapper
50  // process, but they must be deleted whether the parameterized physics is
51  // an 'ExclusivelyForced' process or not. Normally the secondaries will be
52  // created after this wrapper process in G4SteppingManager::InvokePSDIP
53 
54  // store the secondaries from ParticleChange to SecondaryList
55 
56  G4TrackVector* fSecondary = (const_cast<G4Step *> (&step))->GetfSecondary();
57  G4int nSecondarySave = fSecondary->size();
58 
59  G4int num2ndaries = particleChange->GetNumberOfSecondaries();
60  G4Track* tempSecondaryTrack;
61 
62  for(G4int DSecLoop=0 ; DSecLoop< num2ndaries; DSecLoop++){
63  tempSecondaryTrack = particleChange->GetSecondary(DSecLoop);
64 
65  // Set the process pointer which created this track
66  tempSecondaryTrack->SetCreatorProcess( pRegProcess );
67 
68  //add secondaries from this wrapper process to the existing list
69  fSecondary->push_back( tempSecondaryTrack );
70 
71  } //end of loop on secondary
72 
73  // Now we can still impose conditions on the secondaries from ModelTrigger,
74  // such as the number of secondaries produced by the original process as well as
75  // those by other continuous processes - for the hadronic inelastic interaction,
76  // it is always true that particleChange->GetNumberOfSecondaries() > 0
77 
78  // at this stage, updated step information after all processes involved
79  // should be available
80 
81  // Print(step);
82 
83  // call ModelTrigger for the second time - the primary loop of PostStepGPIL
84  // is inside G4SteppingManager::DefinePhysicalStepLength().
85 
86  pmanager = track.GetDefinition()->GetProcessManager();
87 
88  G4double testGPIL = DBL_MAX;
89  G4double fStepLength = 0.0;
90  G4ForceCondition fForceCondition = InActivated;
91 
92  fStepLength = step.GetStepLength();
93 
94  fProcessVector = pmanager->GetPostStepProcessVector(typeDoIt);
95 
96  // keep the current status of track, use fPostponeToNextEvent word for
97  // this particular PostStep GPIL and then restore G4TrackStatus if the
98  // paramterized physics doesn't meet trigger conditions in ModelTrigger
99 
100  const G4TrackStatus keepStatus = track.GetTrackStatus();
101 
102  (const_cast<G4Track *> (&track))->SetTrackStatus(fPostponeToNextEvent);
103 
104  for(G4int ipm = 0 ; ipm < fProcessVector->entries() ; ipm++) {
105  fProcess = (*fProcessVector)(ipm);
106 
107  if ( fProcess->GetProcessType() == fParameterisation ) {
108  // test ModelTrigger via PostStepGPIL
109 
110  testGPIL = fProcess->PostStepGPIL(track,fStepLength,&fForceCondition );
111 
112  // if G4FastSimulationModel:: ModelTrigger is true, then the parameterized
113  // physics process takes over the current process
114 
115  if( fForceCondition == ExclusivelyForced) {
116 
117  // clean up memory for changing the process - counter clean up for
118  // the secondaries created by new G4Track in G4HadronicProcess::FillTotalResult
119  G4int nsec = particleChange->GetNumberOfSecondaries();
120  for(G4int DSecLoop=0 ; DSecLoop< nsec ; DSecLoop++){
121  G4Track* tempSecondaryTrack = particleChange->GetSecondary(DSecLoop);
122  delete tempSecondaryTrack;
123  }
124  particleChange->Clear();
125 
126  // updating G4Step between PostStepGPIL and PostStepDoIt for the parameterized
127  // process may not be necessary, but do it anyway
128 
129  (const_cast<G4Step *> (&step))->SetStepLength(testGPIL);
130  (const_cast<G4Track *> (&track))->SetStepLength(testGPIL);
131 
132  step.GetPostStepPoint()->SetStepStatus(fExclusivelyForcedProc);;
133  step.GetPostStepPoint()->SetProcessDefinedStep(fProcess);
134  step.GetPostStepPoint()->SetSafety(0.0);
135 
136  // invoke PostStepDoIt: equivalent steps for G4SteppingManager::InvokePSDIP
137  particleChange = fProcess->PostStepDoIt(track,step);
138 
139  // update PostStepPoint of Step according to ParticleChange
140  particleChange->UpdateStepForPostStep(const_cast<G4Step *> (&step));
141 
142  // update G4Track according to ParticleChange after each PostStepDoIt
143  (const_cast<G4Step *> (&step))->UpdateTrack();
144 
145  // update safety after each invocation of PostStepDoIts - acutally this
146  // is not necessary for the parameterized physics process, but do it anyway
147  step.GetPostStepPoint()->SetSafety(0.0);
148 
149  // additional nullification
150  (const_cast<G4Track *> (&track))->SetTrackStatus( particleChange->GetTrackStatus() );
151  }
152  else {
153  //restore TrackStatus if fForceCondition != ExclusivelyForced
154  (const_cast<G4Track *> (&track))->SetTrackStatus(keepStatus);
155  }
156  // assume that there is one and only one parameterized physics
157  break;
158  }
159  }
160 
161  // remove secondaries of this wrapper process that were added to the secondary list
162  // since they will be added in the normal stepping procedure after this->PostStepDoIt
163  // in G4SteppingManager::InvokePSDIP
164 
165  //move the iterator to the (nSecondarySave+1)th element in the secondary list
166  G4TrackVector::iterator itv = fSecondary->begin();
167  itv += nSecondarySave;
168 
169  //delete next num2ndaries tracks from the secondary list
170  fSecondary->erase(itv,itv+num2ndaries);
171 
172  //end of specific actions of this wrapper process
173 
174  return particleChange;
175 }
176 
178 
179  std::cout << " GflashHadronWrapperProcess ProcessName, PreStepPosition, preStepPoint KE, PostStepPoint KE, DeltaEnergy Nsec \n "
180  << step.GetPostStepPoint()->GetProcessDefinedStep()->GetProcessName() << " "
181  << step.GetPostStepPoint()->GetPosition() << " "
182  << step.GetPreStepPoint()->GetKineticEnergy()/GeV << " "
183  << step.GetPostStepPoint()->GetKineticEnergy()/GeV << " "
184  << step.GetDeltaEnergy()/GeV << " "
185  << particleChange->GetNumberOfSecondaries() << std::endl;
186 }
GflashHadronWrapperProcess(G4String processName)
G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &step)
tuple cout
Definition: gather_cfg.py:121