Go to the documentation of this file.00001 #include <iostream>
00002 #include <vector>
00003 #include <string>
00004 #include <cassert>
00005 #include <exception>
00006 #include <boost/tuple/tuple.hpp>
00007
00008 #include "FWCore/Framework/interface/Event.h"
00009 #include "FWCore/Framework/interface/EventSetup.h"
00010 #include "FWCore/Framework/interface/ESHandle.h"
00011 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00012
00013 #include "SimG4Core/Notification/interface/BeginOfJob.h"
00014 #include "SimG4Core/Notification/interface/BeginOfEvent.h"
00015 #include "SimG4Core/Notification/interface/BeginOfTrack.h"
00016 #include "SimG4Core/Notification/interface/EndOfTrack.h"
00017
00018 #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
00019
00020
00021 #include "G4Step.hh"
00022 #include "G4Track.hh"
00023 #include "G4VSolid.hh"
00024 #include "G4LogicalVolumeStore.hh"
00025 #include "G4TouchableHistory.hh"
00026 #include "G4VPhysicalVolume.hh"
00027 #include "G4AffineTransform.hh"
00028
00029 #include "TrackingMaterialProducer.h"
00030
00031
00032 static
00033 const G4LogicalVolume* GetVolume(const std::string& name) {
00034 const G4LogicalVolumeStore* lvs = G4LogicalVolumeStore::GetInstance();
00035
00036 for (G4LogicalVolumeStore::const_iterator volume = lvs->begin(); volume != lvs->end(); ++volume) {
00037 if ((const std::string&) (*volume)->GetName() == name)
00038 return (*volume);
00039 }
00040
00041 return 0;
00042 }
00043
00044
00045 static inline
00046 const G4AffineTransform& GetTransform(const G4TouchableHistory* touchable, int depth)
00047 {
00048 return touchable->GetHistory()->GetTransform( touchable->GetHistory()->GetDepth() - depth );
00049 }
00050
00051
00052
00053
00054
00055
00056
00057 boost::tuple<const G4VPhysicalVolume*, int> GetSensitiveVolume( const G4VTouchable* touchable )
00058 {
00059 int depth = touchable->GetHistoryDepth();
00060 for (int level = 0; level < depth; ++level) {
00061 const G4VPhysicalVolume* volume = touchable->GetVolume(level);
00062 if (volume->GetLogicalVolume()->GetSensitiveDetector() != 0) {
00063 return boost::make_tuple(volume, level);
00064 }
00065 }
00066 return boost::tuple<const G4VPhysicalVolume*, int>(0, 0);
00067 }
00068
00069
00070 TrackingMaterialProducer::TrackingMaterialProducer(const edm::ParameterSet& iPSet)
00071 {
00072 edm::ParameterSet config = iPSet.getParameter<edm::ParameterSet>("TrackingMaterialProducer");
00073 m_selectedNames = config.getParameter< std::vector<std::string> >("SelectedVolumes");
00074 m_primaryTracks = config.getParameter<bool>("PrimaryTracksOnly");
00075 m_tracks = 0;
00076
00077 produces< std::vector<MaterialAccountingTrack> >();
00078 }
00079
00080
00081 TrackingMaterialProducer::~TrackingMaterialProducer(void)
00082 {
00083 }
00084
00085
00086 void TrackingMaterialProducer::update(const BeginOfJob* event)
00087 {
00088
00089 std::cout << "TrackingMaterialProducer: List of the selected volumes: " << std::endl;
00090 for (std::vector<std::string>::const_iterator volume_name = m_selectedNames.begin(); volume_name != m_selectedNames.end(); ++volume_name) {
00091 const G4LogicalVolume* volume = GetVolume(*volume_name);
00092 if (volume) {
00093 std::cout << '\t' << *volume_name << std::endl;
00094 m_selectedVolumes.push_back( volume );
00095 } else {
00096
00097 std::cerr << "TrackingMaterialProducer::update(const BeginOfJob*): WARNING: selected volume \"" << *volume_name << "\" not found in geometry " << std::endl;
00098 }
00099 }
00100 std::cout << std::endl;
00101 }
00102
00103
00104
00105 void TrackingMaterialProducer::update(const BeginOfEvent* event)
00106 {
00107 m_tracks = new std::vector<MaterialAccountingTrack>();
00108 }
00109
00110
00111
00112 void TrackingMaterialProducer::update(const BeginOfTrack* event)
00113 {
00114 m_track.reset();
00115
00116
00117 G4Track* track = const_cast<G4Track*>((*event)());
00118 if (m_primaryTracks and track->GetParentID() != 0) {
00119 track->SetTrackStatus(fStopAndKill);
00120 }
00121 }
00122
00123
00124
00125 void TrackingMaterialProducer::update(const G4Step* step)
00126 {
00127 const G4TouchableHistory* touchable = (G4TouchableHistory*)(step->GetTrack()->GetTouchable());
00128 if (not isSelected( touchable )) {
00129
00130 return;
00131 }
00132
00133
00134 const G4Material* material = touchable->GetVolume()->GetLogicalVolume()->GetMaterial();
00135 double length = step->GetStepLength() / cm;
00136 double X0 = material->GetRadlen() / cm;
00137 double Ne = material->GetElectronDensity() * cm3;
00138 double Xi = Ne / 6.0221415e23 * 0.307075 / 2.;
00139 double radiationLengths = length / X0;
00140 double energyLoss = length * Xi;
00141
00142
00143 G4ThreeVector globalPosPre = step->GetPreStepPoint()->GetPosition();
00144 G4ThreeVector globalPosPost = step->GetPostStepPoint()->GetPosition();
00145 GlobalPoint globalPositionIn( globalPosPre.x() / cm, globalPosPre.y() / cm, globalPosPre.z() / cm );
00146 GlobalPoint globalPositionOut( globalPosPost.x() / cm, globalPosPost.y() / cm, globalPosPost.z() / cm );
00147
00148
00149 bool enter_sensitive = false;
00150 bool leave_sensitive = false;
00151 double cosThetaPre = 0.0;
00152 double cosThetaPost = 0.0;
00153 int level = 0;
00154 const G4VPhysicalVolume* sensitive = 0;
00155 GlobalPoint position;
00156 boost::tuples::tie(sensitive, level) = GetSensitiveVolume(touchable);
00157 if (sensitive) {
00158 const G4VSolid & solid = *touchable->GetSolid( level );
00159 const G4AffineTransform & transform = GetTransform( touchable, level );
00160 G4ThreeVector pos = transform.Inverse().TransformPoint( G4ThreeVector( 0., 0., 0. ) );
00161 position = GlobalPoint( pos.x() / cm, pos.y() / cm, pos.z() / cm );
00162
00163 G4ThreeVector localPosPre = transform.TransformPoint( globalPosPre );
00164 EInside statusPre = solid.Inside( localPosPre );
00165 if (statusPre == kSurface) {
00166 enter_sensitive = true;
00167 G4ThreeVector globalDirPre = step->GetPreStepPoint()->GetMomentumDirection();
00168 G4ThreeVector localDirPre = transform.TransformAxis( globalDirPre );
00169 G4ThreeVector normalPre = solid.SurfaceNormal( localPosPre );
00170 cosThetaPre = normalPre.cosTheta( -localDirPre );
00171 }
00172
00173 G4ThreeVector localPosPost = transform.TransformPoint( globalPosPost );
00174 EInside statusPost = solid.Inside( localPosPost );
00175 if (statusPost == kSurface) {
00176 leave_sensitive = true;
00177 G4ThreeVector globalDirPost = step->GetPostStepPoint()->GetMomentumDirection();
00178 G4ThreeVector localDirPost = transform.TransformAxis( globalDirPost );
00179 G4ThreeVector normalPost = solid.SurfaceNormal( localPosPost );
00180 cosThetaPost = normalPost.cosTheta( localDirPost );
00181 }
00182 }
00183
00184
00185 if (enter_sensitive)
00186 m_track.enterDetector( sensitive, position, cosThetaPre );
00187 m_track.step(MaterialAccountingStep( length, radiationLengths, energyLoss, globalPositionIn, globalPositionOut ));
00188 if (leave_sensitive)
00189 m_track.leaveDetector( sensitive, cosThetaPost );
00190
00191
00192
00193
00194
00195
00196
00197
00198
00199
00200
00201
00202
00203
00204
00205
00206
00207
00208
00209 }
00210
00211
00212
00213 void TrackingMaterialProducer::update(const EndOfTrack* event)
00214 {
00215 const G4Track * track = (*event)();
00216 if (m_primaryTracks and track->GetParentID() != 0)
00217 return;
00218
00219 m_tracks->push_back(m_track);
00220
00221
00222 std::cout << "this track took " << m_track.steps().size() << " steps, and passed through " << m_track.detectors().size() << " sensitive detectors" << std::endl;
00223 std::cout << "\ttrack length: " << m_track.summary().length() << " cm" << std::endl;
00224 std::cout << "\tradiation lengths: " << m_track.summary().radiationLengths() << std::endl;
00225 std::cout << "\tenergy loss: " << m_track.summary().energyLoss() << " MeV" << std::endl;
00226
00227
00228
00229
00230
00231
00232
00233
00234
00235
00236
00237 std::cout << std::endl;
00238 }
00239
00240
00241 void TrackingMaterialProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup)
00242 {
00243
00244 std::auto_ptr<std::vector<MaterialAccountingTrack> > tracks( m_tracks );
00245 iEvent.put( tracks );
00246 m_tracks = 0;
00247 }
00248
00249
00250 bool TrackingMaterialProducer::isSelected( const G4VTouchable* touchable )
00251 {
00252 for (size_t i = 0; i < m_selectedVolumes.size(); ++i)
00253 if (m_selectedVolumes[i]->IsAncestor( touchable->GetVolume() ) or m_selectedVolumes[i] == touchable->GetVolume()->GetLogicalVolume())
00254 return true;
00255
00256 return false;
00257 }
00258
00259
00260
00261 #include "SimG4Core/Watcher/interface/SimWatcherFactory.h"
00262 #include "FWCore/Framework/interface/MakerMacros.h"
00263 DEFINE_SIMWATCHER(TrackingMaterialProducer);