Geant4  10.00.p02
WLSTrajectory.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 // $Id: WLSTrajectory.cc 72065 2013-07-05 09:54:59Z gcosmo $
27 //
30 //
31 //
32 #include "G4AttDef.hh"
33 #include "G4AttValue.hh"
34 #include "G4AttDefStore.hh"
35 
36 #include "G4UIcommand.hh"
37 #include "G4UnitsTable.hh"
38 
39 #include "WLSTrajectory.hh"
40 #include "WLSTrajectoryPoint.hh"
41 #include "G4ParticleTable.hh"
42 #include "G4ParticleTypes.hh"
43 
44 #include "G4Polyline.hh"
45 #include "G4Colour.hh"
46 #include "G4VisAttributes.hh"
47 #include "G4VVisManager.hh"
48 #include "G4Polymarker.hh"
49 
50 //#define G4ATTDEBUG
51 #ifdef G4ATTDEBUG
52 #include "G4AttCheck.hh"
53 #endif
54 
56 
57 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
58 
60  : fpPointsContainer(0), fTrackID(0), fParentID(0),
61  fPDGCharge(0.0), fPDGEncoding(0), fParticleName(""),
62  fInitialMomentum(G4ThreeVector())
63 {
64  fWLS = false;
65  fDrawIt = false;
66  fForceNoDraw = false;
67  fForceDraw = false;
68 
69  fParticleDefinition = NULL;
70 }
71 
72 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
73 
75 {
80  fTrackID = aTrack->GetTrackID();
81  fParentID = aTrack->GetParentID();
82  fInitialMomentum = aTrack->GetMomentum();
84  // Following is for the first trajectory point
85  fpPointsContainer->push_back(new WLSTrajectoryPoint(aTrack));
86 
87  fWLS = false;
88  fDrawIt = false;
89  fForceNoDraw = false;
90  fForceDraw = false;
91 }
92 
93 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
94 
96 {
99  fPDGCharge = right.fPDGCharge;
100  fPDGEncoding = right.fPDGEncoding;
101  fTrackID = right.fTrackID;
102  fParentID = right.fParentID;
105 
106  for(size_t i=0;i<right.fpPointsContainer->size();++i) {
107  WLSTrajectoryPoint* rightPoint
108  = (WLSTrajectoryPoint*)((*(right.fpPointsContainer))[i]);
109  fpPointsContainer->push_back(new WLSTrajectoryPoint(*rightPoint));
110  }
111 
112  fWLS = right.fWLS;
113  fDrawIt = right.fDrawIt;
114  fForceNoDraw = right.fForceNoDraw;
115  fForceDraw = right.fForceDraw;
116 }
117 
118 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
119 
121 {
122  for(size_t i=0;i<fpPointsContainer->size();++i){
123  delete (*fpPointsContainer)[i];
124  }
125  fpPointsContainer->clear();
126 
127  delete fpPointsContainer;
128 }
129 
130 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
131 
132 void WLSTrajectory::ShowTrajectory(std::ostream& os) const
133 {
134  // Invoke the default implementation in G4VTrajectory...
136  // ... or override with your own code here.
137 }
138 
139 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
140 
142 {
143  // i_mode is no longer available as an argument of G4VTrajectory.
144  // In this exampple it was always called with an argument of 50.
145  const G4int i_mode = 50;
146  // Consider using commands /vis/modeling/trajectories.
147 
148  // Invoke the default implementation in G4VTrajectory...
149  // G4VTrajectory::DrawTrajectory(i_mode);
150  // ... or override with your own code here.
151 
152  //Taken from G4VTrajectory and modified to select colours based on particle
153  //type and to selectively eliminate drawing of certain trajectories.
154 
155  if (!fForceDraw && (!fDrawIt || fForceNoDraw)) return;
156 
157  // If i_mode>=0, draws a trajectory as a polyline and, if i_mode!=0,
158  // adds markers - yellow circles for step points and magenta squares
159  // for auxiliary points, if any - whose screen size in pixels is
160  // given by std::abs(i_mode)/1000. E.g: i_mode = 5000 gives easily
161  // visible markers.
162 
164  if (!pVVisManager) return;
165 
166  const G4double markerSize = std::abs(i_mode)/1000;
167  G4bool lineRequired (i_mode >= 0);
168  G4bool markersRequired (markerSize > 0.);
169 
170  G4Polyline trajectoryLine;
171  G4Polymarker stepPoints;
172  G4Polymarker auxiliaryPoints;
173 
174  for (G4int i = 0; i < GetPointEntries() ; i++) {
175  G4VTrajectoryPoint* aTrajectoryPoint = GetPoint(i);
176  const std::vector<G4ThreeVector>* auxiliaries
177  = aTrajectoryPoint->GetAuxiliaryPoints();
178  if (auxiliaries) {
179  for (size_t iAux = 0; iAux < auxiliaries->size(); ++iAux) {
180  const G4ThreeVector pos((*auxiliaries)[iAux]);
181  if (lineRequired) {
182  trajectoryLine.push_back(pos);
183  }
184  if (markersRequired) {
185  auxiliaryPoints.push_back(pos);
186  }
187  }
188  }
189  const G4ThreeVector pos(aTrajectoryPoint->GetPosition());
190  if (lineRequired) {
191  trajectoryLine.push_back(pos);
192  }
193  if (markersRequired) {
194  stepPoints.push_back(pos);
195  }
196  }
197 
198  if (lineRequired) {
199  G4Colour colour;
200 
202  if(fWLS) //WLS photons are red
203  colour = G4Colour(1.,0.,0.);
204  else{ //Scintillation and Cerenkov photons are green
205  colour = G4Colour(0.,1.,0.);
206  }
207  }
208  else //All other particles are blue
209  colour = G4Colour(0.,0.,1.);
210 
211  G4VisAttributes trajectoryLineAttribs(colour);
212  trajectoryLine.SetVisAttributes(&trajectoryLineAttribs);
213  pVVisManager->Draw(trajectoryLine);
214  }
215  if (markersRequired) {
216  auxiliaryPoints.SetMarkerType(G4Polymarker::squares);
217  auxiliaryPoints.SetScreenSize(markerSize);
218  auxiliaryPoints.SetFillStyle(G4VMarker::filled);
219  G4VisAttributes auxiliaryPointsAttribs(G4Colour(0.,1.,1.)); // Magenta
220  auxiliaryPoints.SetVisAttributes(&auxiliaryPointsAttribs);
221  pVVisManager->Draw(auxiliaryPoints);
222 
224  stepPoints.SetScreenSize(markerSize);
225  stepPoints.SetFillStyle(G4VMarker::filled);
226  G4VisAttributes stepPointsAttribs(G4Colour(1.,1.,0.)); // Yellow.
227  stepPoints.SetVisAttributes(&stepPointsAttribs);
228  pVVisManager->Draw(stepPoints);
229  }
230 }
231 
232 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
233 
235 {
236  fpPointsContainer->push_back(new WLSTrajectoryPoint(aStep));
237 }
238 
239 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
240 
242 {
243  return (G4ParticleTable::GetParticleTable()->FindParticle(fParticleName));
244 }
245 
246 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
247 
249 {
250  if(!secondTrajectory) return;
251 
252  WLSTrajectory* second = (WLSTrajectory*)secondTrajectory;
253  G4int ent = second->GetPointEntries();
254  // initial point of the second trajectory should not be merged
255  for(G4int i=1; i<ent; ++i) {
256  fpPointsContainer->push_back((*(second->fpPointsContainer))[i]);
257  }
258  delete (*second->fpPointsContainer)[0];
259  second->fpPointsContainer->clear();
260 }
261 
262 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
263 
264 const std::map<G4String,G4AttDef>* WLSTrajectory::GetAttDefs() const
265 {
266  G4bool isNew;
267  std::map<G4String,G4AttDef>* store
268  = G4AttDefStore::GetInstance("Trajectory",isNew);
269 
270  if (isNew) {
271 
272  G4String ID("ID");
273  (*store)[ID] = G4AttDef(ID,"Track ID","Bookkeeping","","G4int");
274 
275  G4String PID("PID");
276  (*store)[PID] = G4AttDef(PID,"Parent ID","Bookkeeping","","G4int");
277 
278  G4String PN("PN");
279  (*store)[PN] = G4AttDef(PN,"Particle Name","Physics","","G4String");
280 
281  G4String Ch("Ch");
282  (*store)[Ch] = G4AttDef(Ch,"Charge","Physics","e+","G4double");
283 
284  G4String PDG("PDG");
285  (*store)[PDG] = G4AttDef(PDG,"PDG Encoding","Physics","","G4int");
286 
287  G4String IMom("IMom");
288  (*store)[IMom] = G4AttDef(IMom,
289  "Momentum of track at start of trajectory",
290  "Physics","G4BestUnit","G4ThreeVector");
291 
292  G4String IMag("IMag");
293  (*store)[IMag] = G4AttDef(IMag,
294  "Magnitude of momentum of track at start of trajectory",
295  "Physics","G4BestUnit","G4double");
296 
297  G4String NTP("NTP");
298  (*store)[NTP] = G4AttDef(NTP,"No. of points","Bookkeeping","","G4int");
299 
300  }
301  return store;
302 }
303 
304 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
305 
306 std::vector<G4AttValue>* WLSTrajectory::CreateAttValues() const
307 {
308  std::vector<G4AttValue>* values = new std::vector<G4AttValue>;
309 
310  values->push_back
312 
313  values->push_back
315 
316  values->push_back(G4AttValue("PN",fParticleName,""));
317 
318  values->push_back
320 
321  values->push_back
323 
324  values->push_back
325  (G4AttValue("IMom",G4BestUnit(fInitialMomentum,"Energy"),""));
326 
327  values->push_back
328  (G4AttValue("IMag",G4BestUnit(fInitialMomentum.mag(),"Energy"),""));
329 
330  values->push_back
332 
333 #ifdef G4ATTDEBUG
334  G4cout << G4AttCheck(values,GetAttDefs());
335 #endif
336  return values;
337 }
virtual void Draw(const G4Circle &, const G4Transform3D &objectTransformation=G4Transform3D())=0
G4ParticleDefinition * GetDefinition() const
virtual const std::map< G4String, G4AttDef > * GetAttDefs() const
G4ThreadLocal G4Allocator< WLSTrajectory > * WLSTrajectoryAllocator
G4int GetParentID() const
virtual void AppendStep(const G4Step *aStep)
CLHEP::Hep3Vector G4ThreeVector
virtual const std::vector< G4ThreeVector > * GetAuxiliaryPoints() const
static G4VVisManager * GetConcreteInstance()
virtual G4VTrajectoryPoint * GetPoint(G4int i) const
virtual void MergeTrajectory(G4VTrajectory *secondTrajectory)
void SetMarkerType(MarkerType)
void SetFillStyle(FillStyle)
static G4String ConvertToString(G4bool boolVal)
Definition: G4UIcommand.cc:357
Definition of the WLSTrajectoryPoint class.
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
G4String fParticleName
#define G4ThreadLocal
Definition: tls.hh:52
int G4int
Definition: G4Types.hh:78
G4ParticleDefinition * GetParticleDefinition()
const G4String & GetParticleName() const
virtual int GetPointEntries() const
G4GLOB_DLL std::ostream G4cout
bool G4bool
Definition: G4Types.hh:79
static const double second
Definition: G4SIunits.hh:138
virtual std::vector< G4AttValue > * CreateAttValues() const
Definition: G4Step.hh:76
G4int GetTrackID() const
virtual void ShowTrajectory(std::ostream &os=G4cout) const
void SetVisAttributes(const G4VisAttributes *)
Definition: G4Visible.cc:80
Definition of the WLSTrajectory class.
G4ThreeVector GetMomentum() const
virtual const G4ThreeVector GetPosition() const =0
WLSTrajectoryPointContainer * fpPointsContainer
static G4ParticleTable * GetParticleTable()
G4ParticleDefinition * fParticleDefinition
std::vector< G4VTrajectoryPoint * > WLSTrajectoryPointContainer
G4double fPDGCharge
std::map< G4String, G4AttDef > * GetInstance(const G4String &storeKey, G4bool &isNew)
static G4OpticalPhoton * OpticalPhotonDefinition()
double G4double
Definition: G4Types.hh:76
virtual void ShowTrajectory(std::ostream &os=G4cout) const
virtual ~WLSTrajectory()
G4double GetPDGCharge() const
virtual void DrawTrajectory() const
G4ThreeVector fInitialMomentum
static const G4double pos
void SetScreenSize(G4double)