SM_SD.cc

Go to the documentation of this file.
00001 //! \file
00002 //!
00003 //! Source file for Symmetric_Moeller sensitive detector class.
00004 //!
00005 //! The GT_SD class defines the member functions for processing the hit
00006 //! information whenever a hit is simulated in the detector.  These
00007 //! routines initialise the hit collection, process hits, and perform
00008 //! the end of event action which fills the hit information arrays.
00009 //!
00010 //! \author Y. Ma based on code for olympus simulation package by D.K. Hasell
00011 //! \version 1.0
00012 //! \date 2010-10-11
00013 //!
00014 //! \ingroup detector
00015 
00016 // *+****1****+****2****+****3****+****4****+****5****+****6****+****7****+****
00017 
00018 // Include the header file and the user header files referenced here.
00019 
00020 #include "SM_SD.h"
00021 #include "SM_Hit.h"
00022 #include "SM_Data.h"
00023 #include "SM_Messenger.h"
00024 
00025 #include "EventAction.h"
00026 
00027 // Include the GEANT4 header files referenced in this file.
00028 
00029 #include "G4VSensitiveDetector.hh"
00030 #include "G4HCofThisEvent.hh"
00031 #include "G4SDManager.hh"
00032 #include "G4Step.hh"
00033 #include "G4TouchableHistory.hh"
00034 #include "G4StepPoint.hh"
00035 #include "G4TouchableHandle.hh"
00036 #include "G4ThreeVector.hh"
00037 
00038 #include "Randomize.hh"
00039 
00040 // Use the standard namespace.
00041 
00042 using namespace std;
00043 
00044 // *+****1****+****2****+****3****+****4****+****5****+****6****+****7****+****
00045 
00046 // SM_SD source code.
00047 
00048 SM_SD::SM_SD( G4String name ) : G4VSensitiveDetector( name ) {
00049 
00050    // Insert hits collection name.
00051 
00052    collectionName.insert( "SM_HC" );
00053 
00054    // Create pointer to GT SD Messenger.
00055 
00056    SM_Messenger::Instance()->setSM_SDptr( this );
00057 
00058    // Set defaults.
00059 
00060    SM_threshold = 0. * eV;
00061    SM_Eresol = 0. * eV;
00062 
00063    // Initially transformations are not defined so Boolean is false.
00064 
00065    for( int i = 0; i < N_SM; ++i ) SM_Transform[i] = false;
00066 }
00067 
00068 // Destructor.
00069 
00070 SM_SD::~SM_SD(){}
00071 
00072 // Initialize hits collection.
00073 
00074 void SM_SD::Initialize( G4HCofThisEvent * HCE ) {
00075 
00076    SM_HC = new
00077       SM_HitsCollection( SensitiveDetectorName, collectionName[0] ); 
00078 
00079    static G4int HCID = -1;   
00080    if( HCID < 0 ) {
00081       HCID = G4SDManager::GetSDMpointer()->GetCollectionID(collectionName[0]);
00082    }
00083 
00084    HCE->AddHitsCollection( HCID, SM_HC );
00085 }
00086 
00087 // Process hits.
00088 
00089 G4bool SM_SD::ProcessHits( G4Step * step, G4TouchableHistory * ROhist ) {
00090 
00091    // Get the total energy deposited in this step.
00092 
00093    G4double edep = step->GetTotalEnergyDeposit();
00094 
00095    // Skip this hit if below threshold.
00096 
00097    if( edep < SM_threshold ) return false;
00098 
00099    // Get the PreStepPoint and the TouchableHandle.
00100 
00101    G4StepPoint * prestep = step->GetPreStepPoint();
00102 
00103    G4TouchableHandle touchable = prestep->GetTouchableHandle();
00104 
00105    // Get copy number (used to identify which SM).
00106 
00107    G4int copyno = touchable->GetCopyNumber();
00108 
00109    // Get the track.
00110 
00111    G4Track * track = step->GetTrack();
00112 
00113    // Get the global time and track ID.
00114 
00115    G4double time = prestep->GetGlobalTime();
00116    G4int id = track->GetTrackID();
00117 
00118    // Get the true world and local coordinate positions.
00119 
00120    G4ThreeVector trueworld = prestep->GetPosition();
00121 
00122    // If this first time for this detector store transformation and inverse.
00123 
00124    if( !SM_Transform[copyno] ) {
00125       SM_WorldtoLocal[copyno] = touchable->GetHistory()->GetTopTransform();
00126       SM_LocaltoWorld[copyno] = SM_WorldtoLocal[copyno].Inverse();
00127       SM_Transform[copyno] = true;
00128    }
00129 
00130    G4ThreeVector truelocal = SM_WorldtoLocal[copyno].TransformPoint(trueworld);
00131 
00132    // Create local and world coordinates smeared by resolution.
00133 
00134    G4ThreeVector local = truelocal;
00135 
00136    G4ThreeVector world = SM_LocaltoWorld[copyno].TransformPoint(local);
00137 
00138    // Create a pointer to a new Hit.
00139 
00140    SM_Hit * hit = new SM_Hit();
00141 
00142    // Fill the Hit information.
00143 
00144    hit->copyno  = copyno;
00145    hit->trackid = id;
00146    hit->edep    = edep + CLHEP::RandGauss::shoot( 0.0, SM_Eresol );
00147    hit->time    = time;
00148    hit->tworld  = trueworld;
00149    hit->tlocal  = truelocal;
00150    hit->world   = world;
00151    hit->local   = local;
00152 
00153    // Add the Hit to the Hit Collection.
00154   
00155  SM_HC->insert( hit );
00156    
00157    //hit->Print();
00158 
00159    return true;
00160 }
00161 
00162 // Routine to process the hits at the end of an event.
00163 
00164 void SM_SD::EndOfEvent( G4HCofThisEvent * HCE ){
00165 
00166    SM_Data * smdata = EventAction::smdata;
00167 
00168    smdata->Reset();
00169 
00170    G4int N_Hits = SM_HC->entries();
00171    
00172    smdata->nSM = N_Hits;
00173 
00174    for( G4int i = 0; i < N_Hits; ++i ) {
00175 
00176       smdata->id.push_back( ( *SM_HC )[i]->copyno );
00177       smdata->tr.push_back( ( *SM_HC )[i]->trackid );
00178 
00179       smdata->e.push_back( ( *SM_HC )[i]->edep / MeV );
00180       smdata->t.push_back( ( *SM_HC )[i]->time / ns );
00181 
00182       smdata->tx.push_back( ( *SM_HC )[i]->tworld.x() / cm );
00183       smdata->ty.push_back( ( *SM_HC )[i]->tworld.y() / cm );
00184       smdata->tz.push_back( ( *SM_HC )[i]->tworld.z() / cm );
00185 
00186       smdata->txl.push_back( ( *SM_HC )[i]->tlocal.x() / cm );
00187       smdata->tyl.push_back( ( *SM_HC )[i]->tlocal.y() / cm );
00188       smdata->tzl.push_back( ( *SM_HC )[i]->tlocal.z() / cm );
00189 
00190       smdata->x.push_back( ( *SM_HC )[i]->world.x() / cm );
00191       smdata->y.push_back( ( *SM_HC )[i]->world.y() / cm );
00192       smdata->z.push_back( ( *SM_HC )[i]->world.z() / cm );
00193 
00194       smdata->xl.push_back( ( *SM_HC )[i]->local.x() / cm );
00195       smdata->yl.push_back( ( *SM_HC )[i]->local.y() / cm );
00196       smdata->zl.push_back( ( *SM_HC )[i]->local.z() / cm );
00197 
00198    }
00199 
00200 }
00201 
00202 void SM_SD::clear(){} 
00203 
00204 void SM_SD::DrawAll() {} 
00205 
00206 // Print all the hits.
00207 
00208 void SM_SD::PrintAll() {
00209    G4int N_Hits = SM_HC->entries();
00210    G4cout << "\nSM Hits Collection N_Hits = " << N_Hits << "\n" << G4endl; 
00211    for ( G4int i = 0; i < N_Hits; ++i ) (*SM_HC)[i]->Print();
00212 } 
00213 
00214 // Set/Get energy threshold.
00215 
00216 G4double SM_SD::setThreshold( G4double thres ) { return SM_threshold = thres; }
00217 G4double SM_SD::getThreshold() { return SM_threshold; }
00218 
00219 // Set/Get Energy resolution.
00220 
00221 G4double SM_SD::setEresol( G4double res ) { return SM_Eresol = res; }
00222 G4double SM_SD::getEresol() { return SM_Eresol; }