Geant4  10.03.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4RTRun.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: G4RTRun.cc 66264 2012-12-14 10:17:44Z allison $
28 //
29 //
30 //
31 
33 //G4RTRun.cc
35 
36 #include "G4RTRun.hh"
37 #include "G4RTRun.hh"
38 #include "G4TheMTRayTracer.hh"
39 
40 #include "G4Colour.hh"
41 #include "G4VisAttributes.hh"
42 #include "G4Event.hh"
43 #include "G4TrajectoryContainer.hh"
44 
45 #include "G4RayTrajectory.hh"
46 #include "G4RayTrajectoryPoint.hh"
47 
49 {
50  colorMap = new G4THitsMap<G4Colour>("G4RTRun","ColorMap");
51 
52  backgroundColour = G4TheMTRayTracer::theInstance->backgroundColour;
53  lightDirection = G4TheMTRayTracer::theInstance->lightDirection;
54  attenuationLength = G4TheMTRayTracer::theInstance->attenuationLength;
55 }
56 
58 {
59  colorMap->clear();
60  delete colorMap;
61 }
62 
63 void G4RTRun::RecordEvent(const G4Event* evt)
64 {
65  G4TrajectoryContainer * trajectoryContainer = evt->GetTrajectoryContainer();
66  if(!trajectoryContainer) return;
67  G4RayTrajectory* trajectory = static_cast<G4RayTrajectory*>( (*trajectoryContainer)[0] );
68  if(!trajectory) return;
69 
70  G4int nPoint = trajectory->GetPointEntries();
71  if(nPoint==0) return;
72 
73  G4int evId = evt->GetEventID();
74  G4Colour initialCol(backgroundColour);
75  if( trajectory->GetPointC(nPoint-1)->GetPostStepAtt() )
76  { initialCol = GetSurfaceColour(trajectory->GetPointC(nPoint-1)); }
77  G4Colour rayColour = Attenuate(trajectory->GetPointC(nPoint-1),initialCol);
78 
79  for(int i=nPoint-2;i>=0;i--)
80  {
81  G4Colour surfaceCol = GetSurfaceColour(trajectory->GetPointC(i));
82  G4double weight = 1.0 - surfaceCol.GetAlpha();
83  G4Colour mixedCol = GetMixedColour(rayColour,surfaceCol,weight);
84  rayColour = Attenuate(trajectory->GetPointC(i),mixedCol);
85  }
86 
87  colorMap->add(evId,rayColour);
88 }
89 
90 void G4RTRun::Merge(const G4Run* aLocalRun)
91 {
92  const G4RTRun* theLocalRun = static_cast<const G4RTRun*>(aLocalRun);
93  if(theLocalRun) *(colorMap) += *(theLocalRun->colorMap);
94  G4Run::Merge(aLocalRun);
95 }
96 
97 G4Colour G4RTRun::GetSurfaceColour(G4RayTrajectoryPoint* point)
98 {
99  const G4VisAttributes* preAtt = point->GetPreStepAtt();
100  const G4VisAttributes* postAtt = point->GetPostStepAtt();
101 
102  G4bool preVis = ValidColour(preAtt);
103  G4bool postVis = ValidColour(postAtt);
104 
105  G4Colour transparent(1.,1.,1.,0.);
106 
107  if(!preVis&&!postVis) return transparent;
108 
110 
111  G4Colour preCol(1.,1.,1.);
112  G4Colour postCol(1.,1.,1.);
113 
114  if(preVis)
115  {
116  G4double brill = (1.0-(-lightDirection).dot(normal))/2.0;
117  G4double red = preAtt->GetColour().GetRed();
118  G4double green = preAtt->GetColour().GetGreen();
119  G4double blue = preAtt->GetColour().GetBlue();
120  preCol = G4Colour
121  (red*brill,green*brill,blue*brill,preAtt->GetColour().GetAlpha());
122  }
123  else
124  { preCol = transparent; }
125 
126  if(postVis)
127  {
128  G4double brill = (1.0-(-lightDirection).dot(-normal))/2.0;
129  G4double red = postAtt->GetColour().GetRed();
130  G4double green = postAtt->GetColour().GetGreen();
131  G4double blue = postAtt->GetColour().GetBlue();
132  postCol = G4Colour
133  (red*brill,green*brill,blue*brill,postAtt->GetColour().GetAlpha());
134  }
135  else
136  { postCol = transparent; }
137 
138  if(!preVis) return postCol;
139  if(!postVis) return preCol;
140 
141  G4double weight = 0.5;
142  return GetMixedColour(preCol,postCol,weight);
143 }
144 
145 G4Colour G4RTRun::GetMixedColour(G4Colour surfCol,G4Colour transCol,G4double weight)
146 {
147  G4double red = weight*surfCol.GetRed() + (1.-weight)*transCol.GetRed();
148  G4double green = weight*surfCol.GetGreen() + (1.-weight)*transCol.GetGreen();
149  G4double blue = weight*surfCol.GetBlue() + (1.-weight)*transCol.GetBlue();
150  G4double alpha = weight*surfCol.GetAlpha() + (1.-weight)*transCol.GetAlpha();
151  return G4Colour(red,green,blue,alpha);
152 }
153 
154 G4Colour G4RTRun::Attenuate(G4RayTrajectoryPoint* point, G4Colour sourceCol)
155 {
156  const G4VisAttributes* preAtt = point->GetPreStepAtt();
157 
158  G4bool visible = ValidColour(preAtt);
159  if(!visible) return sourceCol;
160 
161  G4Colour objCol = preAtt->GetColour();
162  G4double stepRed = objCol.GetRed();
163  G4double stepGreen = objCol.GetGreen();
164  G4double stepBlue = objCol.GetBlue();
165  G4double stepAlpha = objCol.GetAlpha();
166  G4double stepLength = point->GetStepLength();
167 
168  G4double attenuationFuctor;
169  if(stepAlpha > 0.9999999){ stepAlpha = 0.9999999; } // patch to the next line
170  attenuationFuctor = -stepAlpha/(1.0-stepAlpha)*stepLength/attenuationLength;
171 
172  G4double KtRed = std::exp((1.0-stepRed)*attenuationFuctor);
173  G4double KtGreen = std::exp((1.0-stepGreen)*attenuationFuctor);
174  G4double KtBlue = std::exp((1.0-stepBlue)*attenuationFuctor);
175  if(KtRed>1.0){KtRed=1.0;}
176  if(KtGreen>1.0){KtGreen=1.0;}
177  if(KtBlue>1.0){KtBlue=1.0;}
178  return G4Colour(sourceCol.GetRed()*KtRed,
179  sourceCol.GetGreen()*KtGreen,sourceCol.GetBlue()*KtBlue);
180 }
181 
182 G4bool G4RTRun::ValidColour(const G4VisAttributes* visAtt)
183 {
184  G4bool val = true;
185  if(!visAtt)
186  { val = false; }
187  else if(!(visAtt->IsVisible()))
188  { val = false; }
189  else if(visAtt->IsForceDrawingStyle()
191  { val = false; }
192  return val;
193 }
194 
195 
virtual void Merge(const G4Run *)
Definition: G4Run.cc:54
void clear()
Definition: G4THitsMap.hh:267
G4double GetAlpha() const
Definition: G4Colour.hh:142
Definition: test07.cc:36
virtual void Merge(const G4Run *)
Definition: G4RTRun.cc:90
virtual void RecordEvent(const G4Event *)
Definition: G4RTRun.cc:63
const G4VisAttributes * GetPreStepAtt() const
G4bool IsVisible() const
const G4Colour & GetColour() const
G4RTRun()
Definition: G4RTRun.cc:48
G4double GetBlue() const
Definition: G4Colour.hh:141
Definition: test07.cc:36
int G4int
Definition: G4Types.hh:78
G4TrajectoryContainer * GetTrajectoryContainer() const
Definition: G4Event.hh:189
virtual int GetPointEntries() const
static double normal(HepRandomEngine *eptr)
Definition: RandPoisson.cc:77
G4int GetEventID() const
Definition: G4Event.hh:151
G4ThreeVector lightDirection
G4double GetRed() const
Definition: G4Colour.hh:139
bool G4bool
Definition: G4Types.hh:79
G4double GetGreen() const
Definition: G4Colour.hh:140
Definition: G4Run.hh:46
G4ThreeVector GetSurfaceNormal() const
const G4VisAttributes * GetPostStepAtt() const
virtual ~G4RTRun()
Definition: G4RTRun.cc:57
G4int add(const G4int &key, T *&aHit) const
Definition: G4THitsMap.hh:106
G4RayTrajectoryPoint * GetPointC(G4int i) const
G4bool IsForceDrawingStyle() const
G4Colour backgroundColour
double G4double
Definition: G4Types.hh:76
G4double GetStepLength() const
static const G4double alpha
G4double attenuationLength
ForcedDrawingStyle GetForcedDrawingStyle() const