Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4RichTrajectoryPoint.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: G4RichTrajectoryPoint.cc 69003 2013-04-15 09:25:23Z gcosmo $
28 //
29 //
30 // ---------------------------------------------------------------
31 //
32 // G4RichTrajectoryPoint.cc
33 //
34 // Contact:
35 // Questions and comments on G4TrajectoryPoint, on which this is based,
36 // should be sent to
37 // Katsuya Amako (e-mail: Katsuya.Amako@kek.jp)
38 // Makoto Asai (e-mail: asai@kekvax.kek.jp)
39 // Takashi Sasaki (e-mail: Takashi.Sasaki@kek.jp)
40 // and on the extended code to:
41 // John Allison (e-mail: John.Allison@manchester.ac.uk)
42 // Joseph Perl (e-mail: perl@slac.stanford.edu)
43 //
44 // ---------------------------------------------------------------
45 
46 #include "G4RichTrajectoryPoint.hh"
47 
48 #include "G4Track.hh"
49 #include "G4Step.hh"
50 #include "G4VProcess.hh"
51 
52 #include "G4AttDefStore.hh"
53 #include "G4AttDef.hh"
54 #include "G4AttValue.hh"
55 #include "G4UnitsTable.hh"
56 
57 //#define G4ATTDEBUG
58 #ifdef G4ATTDEBUG
59 #include "G4AttCheck.hh"
60 #endif
61 
62 #include <sstream>
63 
65 
67  fpAuxiliaryPointVector(0),
68  fTotEDep(0.),
69  fRemainingEnergy(0.),
70  fpProcess(0),
71  fPreStepPointStatus(fUndefined),
72  fPostStepPointStatus(fUndefined),
73  fPreStepPointGlobalTime(0),
74  fPostStepPointGlobalTime(0),
75  fPreStepPointWeight(1.),
76  fPostStepPointWeight(1.)
77 {}
78 
80  G4TrajectoryPoint(aTrack->GetPosition()),
81  fpAuxiliaryPointVector(0),
82  fTotEDep(0.),
83  fRemainingEnergy(aTrack->GetKineticEnergy()),
84  fpProcess(0),
85  fPreStepPointStatus(fUndefined),
86  fPostStepPointStatus(fUndefined),
87  fPreStepPointGlobalTime(aTrack->GetGlobalTime()),
88  fPostStepPointGlobalTime(aTrack->GetGlobalTime()),
89  fpPreStepPointVolume(aTrack->GetTouchableHandle()),
90  fpPostStepPointVolume(aTrack->GetNextTouchableHandle()),
91  fPreStepPointWeight(aTrack->GetWeight()),
92  fPostStepPointWeight(aTrack->GetWeight())
93 {}
94 
96  G4TrajectoryPoint(aStep->GetPostStepPoint()->GetPosition()),
97  fpAuxiliaryPointVector(aStep->GetPointerToVectorOfAuxiliaryPoints()),
98  fTotEDep(aStep->GetTotalEnergyDeposit())
99 {
100  G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
101  G4StepPoint* postStepPoint = aStep->GetPostStepPoint();
102  if (aStep->GetTrack()->GetCurrentStepNumber() <= 0) { // First step
103  fRemainingEnergy = aStep->GetTrack()->GetKineticEnergy();
104  } else {
105  fRemainingEnergy = preStepPoint->GetKineticEnergy() - fTotEDep;
106  }
107  fpProcess = postStepPoint->GetProcessDefinedStep();
108  fPreStepPointStatus = preStepPoint->GetStepStatus();
109  fPostStepPointStatus = postStepPoint->GetStepStatus();
110  fPreStepPointGlobalTime = preStepPoint->GetGlobalTime();
111  fPostStepPointGlobalTime = postStepPoint->GetGlobalTime();
112  fpPreStepPointVolume = preStepPoint->GetTouchableHandle();
113  fpPostStepPointVolume = postStepPoint->GetTouchableHandle();
114  fPreStepPointWeight = preStepPoint->GetWeight();
115  fPostStepPointWeight = postStepPoint->GetWeight();
116 }
117 
120  G4TrajectoryPoint(right),
121  fpAuxiliaryPointVector(right.fpAuxiliaryPointVector),
122  fTotEDep(right.fTotEDep),
123  fRemainingEnergy(right.fRemainingEnergy),
124  fpProcess(right.fpProcess),
125  fPreStepPointStatus(right.fPreStepPointStatus),
126  fPostStepPointStatus(right.fPostStepPointStatus),
127  fPreStepPointGlobalTime(right.fPreStepPointGlobalTime),
128  fPostStepPointGlobalTime(right.fPostStepPointGlobalTime),
129  fpPreStepPointVolume(right.fpPreStepPointVolume),
130  fpPostStepPointVolume(right.fpPostStepPointVolume),
131  fPreStepPointWeight(right.fPreStepPointWeight),
132  fPostStepPointWeight(right.fPostStepPointWeight)
133 {}
134 
136 {
137  if(fpAuxiliaryPointVector) {
138  /*
139  G4cout << "Deleting fpAuxiliaryPointVector at "
140  << (void*) fpAuxiliaryPointVector
141  << G4endl;
142  */
143  delete fpAuxiliaryPointVector;
144  }
145 }
146 
147 const std::map<G4String,G4AttDef>*
149 {
150  G4bool isNew;
151  std::map<G4String,G4AttDef>* store
152  = G4AttDefStore::GetInstance("G4RichTrajectoryPoint",isNew);
153  if (isNew) {
154 
155  // Copy base class att defs...
156  *store = *(G4TrajectoryPoint::GetAttDefs());
157 
158  G4String ID;
159 
160  ID = "Aux";
161  (*store)[ID] = G4AttDef(ID, "Auxiliary Point Position",
162  "Physics","G4BestUnit","G4ThreeVector");
163  ID = "TED";
164  (*store)[ID] = G4AttDef(ID,"Total Energy Deposit",
165  "Physics","G4BestUnit","G4double");
166  ID = "RE";
167  (*store)[ID] = G4AttDef(ID,"Remaining Energy",
168  "Physics","G4BestUnit","G4double");
169  ID = "PDS";
170  (*store)[ID] = G4AttDef(ID,"Process Defined Step",
171  "Physics","","G4String");
172  ID = "PTDS";
173  (*store)[ID] = G4AttDef(ID,"Process Type Defined Step",
174  "Physics","","G4String");
175  ID = "PreStatus";
176  (*store)[ID] = G4AttDef(ID,"Pre-step-point status",
177  "Physics","","G4String");
178  ID = "PostStatus";
179  (*store)[ID] = G4AttDef(ID,"Post-step-point status",
180  "Physics","","G4String");
181  ID = "PreT";
182  (*store)[ID] = G4AttDef(ID,"Pre-step-point global time",
183  "Physics","G4BestUnit","G4double");
184  ID = "PostT";
185  (*store)[ID] = G4AttDef(ID,"Post-step-point global time",
186  "Physics","G4BestUnit","G4double");
187  ID = "PreVPath";
188  (*store)[ID] = G4AttDef(ID,"Pre-step Volume Path",
189  "Physics","","G4String");
190  ID = "PostVPath";
191  (*store)[ID] = G4AttDef(ID,"Post-step Volume Path",
192  "Physics","","G4String");
193  ID = "PreW";
194  (*store)[ID] = G4AttDef(ID,"Pre-step-point weight",
195  "Physics","","G4double");
196  ID = "PostW";
197  (*store)[ID] = G4AttDef(ID,"Post-step-point weight",
198  "Physics","","G4double");
199  }
200  return store;
201 }
202 
203 static G4String Status(G4StepStatus stps)
204 {
206  switch (stps) {
207  case fWorldBoundary: status = "fWorldBoundary"; break;
208  case fGeomBoundary: status = "fGeomBoundary"; break;
209  case fAtRestDoItProc: status = "fAtRestDoItProc"; break;
210  case fAlongStepDoItProc: status = "fAlongStepDoItProc"; break;
211  case fPostStepDoItProc: status = "fPostStepDoItProc"; break;
212  case fUserDefinedLimit: status = "fUserDefinedLimit"; break;
213  case fExclusivelyForcedProc: status = "fExclusivelyForcedProc"; break;
214  case fUndefined: status = "fUndefined"; break;
215  default: status = "Not recognised"; break;
216  }
217  return status;
218 }
219 
220 static G4String Path(const G4TouchableHandle& th)
221 {
222  std::ostringstream oss;
223  G4int depth = th->GetHistoryDepth();
224  for (G4int i = depth; i >= 0; --i) {
225  oss << th->GetVolume(i)->GetName()
226  << ':' << th->GetCopyNumber(i);
227  if (i != 0) oss << '/';
228  }
229  return oss.str();
230 }
231 
232 std::vector<G4AttValue>* G4RichTrajectoryPoint::CreateAttValues() const
233 {
234  // Create base class att values...
235  std::vector<G4AttValue>* values = G4TrajectoryPoint::CreateAttValues();
236 
237  if (fpAuxiliaryPointVector) {
238  std::vector<G4ThreeVector>::iterator iAux;
239  for (iAux = fpAuxiliaryPointVector->begin();
240  iAux != fpAuxiliaryPointVector->end(); ++iAux) {
241  values->push_back(G4AttValue("Aux",G4BestUnit(*iAux,"Length"),""));
242  }
243  }
244 
245  values->push_back(G4AttValue("TED",G4BestUnit(fTotEDep,"Energy"),""));
246 
247  values->push_back(G4AttValue("RE",G4BestUnit(fRemainingEnergy,"Energy"),""));
248 
249  if (fpProcess) {
250  values->push_back
251  (G4AttValue("PDS",fpProcess->GetProcessName(),""));
252  values->push_back
253  (G4AttValue
254  ("PTDS",G4VProcess::GetProcessTypeName(fpProcess->GetProcessType()),
255  ""));
256  } else {
257  values->push_back(G4AttValue("PDS","None",""));
258  values->push_back(G4AttValue("PTDS","None",""));
259  }
260 
261  values->push_back
262  (G4AttValue("PreStatus",Status(fPreStepPointStatus),""));
263 
264  values->push_back
265  (G4AttValue("PostStatus",Status(fPostStepPointStatus),""));
266 
267  values->push_back
268  (G4AttValue("PreT",G4BestUnit(fPreStepPointGlobalTime,"Time"),""));
269 
270  values->push_back
271  (G4AttValue("PostT",G4BestUnit(fPostStepPointGlobalTime,"Time"),""));
272 
273  if (fpPreStepPointVolume && fpPreStepPointVolume->GetVolume()) {
274  values->push_back(G4AttValue("PreVPath",Path(fpPreStepPointVolume),""));
275  } else {
276  values->push_back(G4AttValue("PreVPath","None",""));
277  }
278 
279  if (fpPostStepPointVolume && fpPostStepPointVolume->GetVolume()) {
280  values->push_back(G4AttValue("PostVPath",Path(fpPostStepPointVolume),""));
281  } else {
282  values->push_back(G4AttValue("PostVPath","None",""));
283  }
284 
285  {
286  std::ostringstream oss;
287  oss << fPreStepPointWeight;
288  values->push_back
289  (G4AttValue("PreW",oss.str(),""));
290  }
291 
292  {
293  std::ostringstream oss;
294  oss << fPostStepPointWeight;
295  values->push_back
296  (G4AttValue("PostW",oss.str(),""));
297  }
298 
299 #ifdef G4ATTDEBUG
300  G4cout << G4AttCheck(values,GetAttDefs());
301 #endif
302 
303  return values;
304 }
static const G4String & GetProcessTypeName(G4ProcessType)
Definition: G4VProcess.cc:141
G4double GetWeight() const
G4StepStatus GetStepStatus() const
virtual std::vector< G4AttValue > * CreateAttValues() const
G4ThreadLocal G4Allocator< G4RichTrajectoryPoint > * aRichTrajectoryPointAllocator
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
G4ProcessType GetProcessType() const
Definition: G4VProcess.hh:414
#define G4ThreadLocal
Definition: tls.hh:52
G4int GetCopyNumber(G4int depth=0) const
int G4int
Definition: G4Types.hh:78
G4StepStatus
Definition: G4StepStatus.hh:49
G4StepPoint * GetPreStepPoint() const
G4double GetKineticEnergy() const
G4GLOB_DLL std::ostream G4cout
G4int GetCurrentStepNumber() const
const G4String & GetName() const
bool G4bool
Definition: G4Types.hh:79
virtual std::vector< G4AttValue > * CreateAttValues() const
Definition: G4Step.hh:76
virtual G4int GetHistoryDepth() const
Definition: G4VTouchable.cc:79
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
int status
Definition: tracer.cxx:24
const G4VProcess * GetProcessDefinedStep() const
virtual G4VPhysicalVolume * GetVolume(G4int depth=0) const
Definition: G4VTouchable.cc:44
G4StepPoint * GetPostStepPoint() const
G4double GetGlobalTime() const
std::map< G4String, G4AttDef > * GetInstance(const G4String &storeKey, G4bool &isNew)
G4double GetKineticEnergy() const
G4Track * GetTrack() const
const G4TouchableHandle & GetTouchableHandle() const
virtual const std::map< G4String, G4AttDef > * GetAttDefs() const
virtual const std::map< G4String, G4AttDef > * GetAttDefs() const