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