Forked from
Chad Lantz / zdcLG
This fork has diverged from the upstream repository.
-
Chad Lantz authoredChad Lantz authored
PMTSD.cc 4.91 KiB
// ********************************************************************
// * License and Disclaimer *
// * *
// * The Geant4 software is copyright of the Copyright Holders of *
// * the Geant4 Collaboration. It is provided under the terms and *
// * conditions of the Geant4 Software License, included in the file *
// * LICENSE and available at http://cern.ch/geant4/license . These *
// * include a list of copyright holders. *
// * *
// * Neither the authors of this software system, nor their employing *
// * institutes,nor the agencies providing financial support for this *
// * work make any representation or warranty, express or implied, *
// * regarding this software system or assume any liability for its *
// * use. Please see the license in the file LICENSE and URL above *
// * for the full disclaimer and the limitation of liability. *
// * *
// * This code implementation is the result of the scientific and *
// * technical work of the GEANT4 collaboration. *
// * By using, copying, modifying or distributing the software (or *
// * any work based on the software) you agree to acknowledge its *
// * use in resulting scientific publications, and indicate your *
// * acceptance of all terms of the Geant4 Software license. *
// ********************************************************************
//
// @Author Chad Lantz
#include "PMTSD.hh"
#include "PMTHit.hh"
#include "RunAction.hh"
#include "lgAnalysis.hh"
#include "G4HCofThisEvent.hh"
#include "G4Step.hh"
#include "G4ThreeVector.hh"
#include "G4SDManager.hh"
#include "G4RunManager.hh"
#include "G4UImanager.hh"
#include "G4ios.hh"
#include "G4Poisson.hh"
#include "G4ParticleTypes.hh"
#include "G4ParticleDefinition.hh"
#include "TString.h"
#include <string>
#include <iostream>
#include <cmath>
/*
*
*/
PMTSD::PMTSD(G4String sdName, G4bool isIntermediate)
: G4VSensitiveDetector(sdName),
INTERMEDIATE(isIntermediate)
{
collectionName.insert(sdName);
HCID = -1;
VISUALIZE = false;
int nVectors = (INTERMEDIATE) ? 3 : 9;
for(int i = 0; i < nVectors; i++) m_PtrVec.push_back( new std::vector<double> );
}
/*
*
*/
PMTSD::~PMTSD(){
for(auto vec : m_PtrVec ) delete vec;
}
/* Mandatory method
*
*/
void PMTSD::Initialize(G4HCofThisEvent* HCE){
hitCollection = new HitsCollection(GetName(), collectionName[0]);
std::string name = GetName();
if(HCID<0)
{ HCID = G4SDManager::GetSDMpointer()->GetCollectionID( name );}
HCE->AddHitsCollection( HCID, hitCollection );
// If the UI window exists, set visualize to true so we record and draw hits
if( G4UImanager::GetUIpointer()->GetG4UIWindow() ) VISUALIZE = true;
m_nHits = 0;
for(auto vec : m_PtrVec ) vec->clear();
}
/* Mandatory method
*
*/
G4bool PMTSD::ProcessHits(G4Step* aStep,G4TouchableHistory*){
G4Track* theTrack = aStep->GetTrack();
m_nHits++;
G4ThreeVector origin = theTrack->GetVertexPosition();
G4ThreeVector hitPos = theTrack->GetPosition();
G4ThreeVector momentum = aStep->GetPreStepPoint()->GetMomentum();
//Display hits if the UI is on
if(VISUALIZE){
PMTHit* newHit = new PMTHit();
newHit->setPos( theTrack->GetPosition() );
hitCollection->insert ( newHit );
}
// Fill ouput vectors
if(INTERMEDIATE){
m_PtrVec[0]->push_back( hitPos.x() ); //Where the photon landed on the PMT in x
m_PtrVec[1]->push_back( hitPos.y() ); //Where the photon landed on the PMT in y
m_PtrVec[2]->push_back( hitPos.z() ); //Where the photon landed on the PMT in z
}else{
m_PtrVec[0]->push_back( origin.x() ); //Origin of the photon in x
m_PtrVec[1]->push_back( origin.y() ); //Origin of the photon in y
m_PtrVec[2]->push_back( origin.z() ); //Origin of the photon in z
m_PtrVec[3]->push_back( hitPos.x() ); //Where the photon landed on the PMT in x
m_PtrVec[4]->push_back( hitPos.y() ); //Where the photon landed on the PMT in y
m_PtrVec[5]->push_back( hitPos.z() ); //Where the photon landed on the PMT in z
m_PtrVec[6]->push_back( theTrack->GetGlobalTime() ); //Time of arrival of the photon at the PMT
G4double pTheta = atan(sqrt(pow(momentum.x(),2) + pow(momentum.z(),2) )/fabs(momentum.y() ));
m_PtrVec[7]->push_back( pTheta ); //Incident angle at the PMT about the y-axis
m_PtrVec[8]->push_back( aStep->GetPreStepPoint()->GetTotalEnergy() ); //Photon energy
}
return true;
}
/* Mandatory method
*
*/
void PMTSD::EndOfEvent(G4HCofThisEvent*){
if(INTERMEDIATE){
G4AnalysisManager::Instance()->FillNtupleIColumn(m_NtupleNo, 0, m_nHits);
}
}