Skip to content
Snippets Groups Projects
Commit adc59c43 authored by Adrien Matta's avatar Adrien Matta :skull_crossbones:
Browse files

* Adding missing InteractionScorers file

parent 48e56f6a
No related branches found
No related tags found
No related merge requests found
/*****************************************************************************
* Copyright (C) 2009-2016 this file is part of the NPTool Project *
* *
* For the licensing terms see $NPTOOL/Licence/NPTool_Licence *
* For the list of contributors see $NPTOOL/Licence/Contributors *
*****************************************************************************/
/*****************************************************************************
* Original Author: Adrien MATTA contact address: matta@lpccaen.in2p3.fr *
* *
* Creation Date : February 2013 *
* Last update : *
*---------------------------------------------------------------------------*
* Decription: *
* File old the scorer to record Hit energy,time and position *
* *
*---------------------------------------------------------------------------*
* Comment: *
* *
*****************************************************************************/
#include "InteractionScorers.hh"
#include "G4UnitsTable.hh"
using namespace InteractionScorers ;
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
vector<InteractionData>::iterator InteractionDataVector::find(const unsigned int& index){
for(vector<InteractionData>::iterator it= m_Data.begin() ; it !=m_Data.end() ; it++){
if((*it).GetIndex()==index)
return it;
}
return m_Data.end();
}
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
PS_Interactions::PS_Interactions(G4String name,TInteractionCoordinates* Inter, int depth) :G4VPrimitiveScorer(name, depth){
m_Level = depth;
m_InterractionCoordinates=Inter;
}
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
G4bool PS_Interactions::ProcessHits(G4Step* aStep, G4TouchableHistory*){
static G4StepPoint* point;
point = aStep->GetPreStepPoint();
t_Energy = aStep->GetTotalEnergyDeposit();
t_Time = point->GetGlobalTime();
t_Position = point->GetPosition();
t_Index = aStep->GetTrack()->GetTrackID();
vector<InteractionData>::iterator it;
it = m_DataVector.find(t_Index);
if(it!=m_DataVector.end())
it->Add(t_Energy);
else
m_DataVector.Set(t_Index,t_Energy,t_Time,t_Position.x(),t_Position.y(),t_Position.z(),t_Position.theta(),t_Position.phi());
return TRUE;
}
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
void PS_Interactions::Initialize(G4HCofThisEvent*){
// Clear is called by EventAction
}
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
void PS_Interactions::EndOfEvent(G4HCofThisEvent*){
unsigned int size = m_DataVector.size();
for(unsigned int i = 0 ; i < size ; i++)
m_InterractionCoordinates->SetInteraction(m_DataVector[i]->GetEnergy(),m_DataVector[i]->GetTime(),m_DataVector[i]->GetPositionX(),m_DataVector[i]->GetPositionY(),m_DataVector[i]->GetPositionZ(),m_DataVector[i]->GetTheta()/deg,m_DataVector[i]->GetPhi()/deg);
m_DataVector.clear();
}
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
void PS_Interactions::clear(){
m_DataVector.clear();
}
#ifndef InteractionScorers_h
#define InteractionScorers_h 1
/*****************************************************************************
* Copyright (C) 2009-2016 this file is part of the NPTool Project *
* *
* For the licensing terms see $NPTOOL/Licence/NPTool_Licence *
* For the list of contributors see $NPTOOL/Licence/Contributors *
*****************************************************************************/
/*****************************************************************************
* Original Author: Adrien MATTA contact address: matta@lpccaen.in2p3.fr *
* *
* Creation Date : February 2013 *
* Last update : *
*---------------------------------------------------------------------------*
* Decription: *
* File old the scorer to record Hit energy,time and position *
* *
*---------------------------------------------------------------------------*
* Comment: *
* *
*****************************************************************************/
#include "G4VPrimitiveScorer.hh"
#include "TInteractionCoordinates.h"
#include "NPImage.h"
#include <map>
using namespace std;
using namespace CLHEP;
namespace InteractionScorers {
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
class InteractionData{
public:
InteractionData(){m_Index=0;};
InteractionData(const unsigned int& Index ,const double& Energy, const double& Time , const double& PositionX, const double& PositionY, const double& PositionZ, const double& Theta, const double& Phi){
m_Index = Index;
m_Energy = Energy;
m_Time = Time;
m_PositionX = PositionX;;
m_PositionY = PositionY;
m_PositionZ = PositionZ;
m_Theta = Theta;
m_Phi = Phi;
}
~InteractionData(){};
private:
unsigned int m_Index;
double m_Energy;
double m_Time;
double m_PositionX;
double m_PositionY;
double m_PositionZ;
double m_Theta;
double m_Phi;
public:
unsigned int GetIndex() const{return m_Index;};
double GetEnergy() const{return m_Energy;};
double GetTime() const{return m_Time;};
double GetPositionX() const{return m_PositionX;};
double GetPositionY() const{return m_PositionY;};
double GetPositionZ() const{return m_PositionZ;};
double GetTheta() const{return m_Theta;};
double GetPhi() const{return m_Phi;};
public:
void Set(const unsigned int& Index, const double& Energy, const double& Time , const double& PositionX, const double& PositionY, const double& PositionZ, const double& Theta, const double& Phi){
m_Index = Index;
m_Energy = Energy;
m_Time = Time;
m_PositionX = PositionX;;
m_PositionY = PositionY;
m_PositionZ = PositionZ;
m_Theta = Theta;
m_Phi = Phi;
}
void Add(const double& Energy){m_Energy+=Energy;};
unsigned int GetIndex(){return m_Index;};
};
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
// Manage a vector of Interaction hit
class InteractionDataVector{
public:
InteractionDataVector(){};
~InteractionDataVector(){};
private:
vector<InteractionData> m_Data;
public:
vector<InteractionData>::iterator find(const unsigned int& index) ;
void clear(){m_Data.clear();} ;
vector<InteractionData>::iterator end() {return m_Data.end();};
vector<InteractionData>::iterator begin() {return m_Data.begin();};
unsigned int size() {return m_Data.size();};
void Add(const unsigned int& index,const double& Energy) {find(index)->Add(Energy);};
void Set(const unsigned int& index,const double& Energy, const double& Time , const double& PositionX, const double& PositionY, const double& PositionZ, const double& Theta, const double& Phi) {m_Data.push_back(InteractionData(index,Energy,Time,PositionX,PositionY,PositionZ,Theta,Phi));};
InteractionData* operator[](const unsigned int& i){return &m_Data[i];};
};
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
class PS_Interactions : public G4VPrimitiveScorer{
public: // with description
PS_Interactions(G4String name, TInteractionCoordinates* Inter,G4int depth=0);
~PS_Interactions(){};
protected: // with description
G4bool ProcessHits(G4Step*, G4TouchableHistory*);
public:
void Initialize(G4HCofThisEvent*);
void EndOfEvent(G4HCofThisEvent*);
void clear();
void DrawAll(){};
void PrintAll(){m_InterractionCoordinates->Dump();};
// Level at which to find the copy number linked to the detector number
G4int m_Level;
private:
InteractionDataVector m_DataVector;
TInteractionCoordinates* m_InterractionCoordinates;
unsigned int t_Index;
double t_Energy;
double t_Time;
G4ThreeVector t_Position;
};
}
#endif
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment