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

* Adding SamuraiFDC2

parent 5abb3b07
No related branches found
No related tags found
No related merge requests found
Pipeline #86132 passed
add_custom_command(OUTPUT TSamuraiFDC2DataDict.cxx COMMAND ${CMAKE_BINARY_DIR}/scripts/build_dict.sh TSamuraiFDC2Data.h TSamuraiFDC2DataDict.cxx TSamuraiFDC2Data.rootmap libNPSamurai.dylib DEPENDS TSamuraiFDC2Data.h)
add_custom_command(OUTPUT TSamuraiFDC2PhysicsDict.cxx COMMAND ${CMAKE_BINARY_DIR}/scripts/build_dict.sh TSamuraiFDC2Physics.h TSamuraiFDC2PhysicsDict.cxx TSamuraiFDC2Physics.rootmap libNPSamurai.dylib DEPENDS TSamuraiFDC2Physics.h)
add_library(NPSamurai SHARED TSamuraiFDC2Data.cxx TSamuraiFDC2DataDict.cxx TSamuraiFDC2Physics.cxx TSamuraiFDC2PhysicsDict.cxx)
target_link_libraries(NPSamurai ${ROOT_LIBRARIES} NPCore)
install(FILES TSamuraiFDC2Data.h TSamuraiFDC2Physics.h DESTINATION ${CMAKE_INCLUDE_OUTPUT_DIRECTORY})
#include "TSamuraiFDC2Data.h"
#include <iostream>
TSamuraiFDC2Data::TSamuraiFDC2Data(){};
TSamuraiFDC2Data::~TSamuraiFDC2Data(){};
////////////////////////////////////////////////////////////////////////////////
void TSamuraiFDC2Data::SetData(const int& Det, const int& Layer, const int& Wire, const double& Time, const int& Edge){
fDC_DetectorNbr.push_back(Det);
fDC_LayerNbr.push_back(Layer);
fDC_WireNbr.push_back(Wire);
fDC_Time.push_back(Time);
fDC_Edge.push_back(Edge);
}
////////////////////////////////////////////////////////////////////////////////
void TSamuraiFDC2Data::Clear(){
fDC_DetectorNbr.clear();
fDC_LayerNbr.clear();
fDC_WireNbr.clear();
fDC_Time.clear();
fDC_Edge.clear();
}
////////////////////////////////////////////////////////////////////////////////
void TSamuraiFDC2Data::Print(){
using namespace std;
cout << " -- Event:" << endl;
cout << " - Multiplicity: " << Mult() << endl;
}
////////////////////////////////////////////////////////////////////////////////
unsigned int TSamuraiFDC2Data::MultLayer(unsigned int det , unsigned int layer, int edge){
unsigned int mult=0;
unsigned int size = fDC_DetectorNbr.size();
for(unsigned int i = 0 ; i< size ; i++ ){
if(fDC_DetectorNbr[i]==det)
if(fDC_LayerNbr[i]==layer)
if(fDC_Edge[i]==edge && edge!=-1) // edge type is specified (0 or 1)
mult++;
else if(edge==-1)// edge type is not specified
mult++;
}
return mult;
}
////////////////////////////////////////////////////////////////////////////////
std::vector<int> TSamuraiFDC2Data::GetWire(unsigned int det , unsigned int layer){
std::vector<int> wires;
unsigned int size = fDC_DetectorNbr.size();
for(unsigned int i = 0 ; i< size ; i++ ){
if(fDC_DetectorNbr[i]==det)
if(fDC_LayerNbr[i]==layer)
wires.push_back(fDC_WireNbr[i]);
}
return wires;
}
ClassImp(TSamuraiFDC2Data);
#ifndef TDCDATA_H
#define TDCDATA_H
#include "TObject.h"
#include <vector>
class TSamuraiFDC2Data: public TObject{
public:
TSamuraiFDC2Data();
~TSamuraiFDC2Data();
private:
std::vector<int> fDC_DetectorNbr;
std::vector<int> fDC_LayerNbr;
std::vector<int> fDC_WireNbr;
std::vector<double> fDC_Time;
std::vector<int> fDC_Edge;
public:
void Clear();
void Print();
void Clear(const Option_t*) {};
void Dump() const{};
public:
void SetData(const int& Det, const int& Layer, const int& Wire, const double& Time, const int& Edge);
unsigned int Mult(){return fDC_DetectorNbr.size();};
unsigned int MultLayer(unsigned int det , unsigned int layer, int edge=-1);
std::vector<int> GetWire(unsigned int det , unsigned int layer);
int const GetDetectorNbr(const unsigned int& i){return fDC_DetectorNbr[i];};
int const GetLayerNbr(const unsigned int& i){return fDC_LayerNbr[i];};
int const GetWireNbr(const unsigned int& i){return fDC_WireNbr[i];};
double const GetTime(const unsigned int& i){return fDC_Time[i];};
int const GetEdge(const unsigned int& i){return fDC_Edge[i];};
ClassDef(TSamuraiFDC2Data,1);
};
#endif
/*****************************************************************************
* 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 : October 2020 *
* Last update : *
*---------------------------------------------------------------------------*
* Decription: *
* This class hold SamuraiFDC2 treated data *
* *
*---------------------------------------------------------------------------*
* Comment: *
* *
*****************************************************************************/
#include "TSamuraiFDC2Physics.h"
// STL
#include <sstream>
#include <iostream>
#include <cmath>
#include <stdlib.h>
#include <limits>
// NPL
#include "RootInput.h"
#include "RootOutput.h"
#include "TAsciiFile.h"
#include "NPOptionManager.h"
#include "NPDetectorFactory.h"
// ROOT
#include "TChain.h"
///////////////////////////////////////////////////////////////////////////
ClassImp(TSamuraiFDC2Physics)
///////////////////////////////////////////////////////////////////////////
TSamuraiFDC2Physics::TSamuraiFDC2Physics(){
m_EventData = new TSamuraiFDC2Data ;
m_PreTreatedData = new TSamuraiFDC2Data ;
m_EventPhysics = this ;
//m_Spectra = NULL;
}
///////////////////////////////////////////////////////////////////////////
void TSamuraiFDC2Physics::BuildSimplePhysicalEvent(){
BuildPhysicalEvent();
}
///////////////////////////////////////////////////////////////////////////
void TSamuraiFDC2Physics::BuildPhysicalEvent(){
PreTreat();
return;
}
///////////////////////////////////////////////////////////////////////////
void TSamuraiFDC2Physics::PreTreat(){
ClearPreTreatedData();
unsigned int size = m_EventData->Mult();
for(unsigned int i = 0 ; i < size ; i++){
// EDGE=0 is the leading edge, IE, the real time.
// EDGE=1 is the trailing edge, so it helps build Tot
if(m_EventData->GetEdge(i)==0){
int det = m_EventData->GetDetectorNbr(i);
int layer = m_EventData->GetLayerNbr(i);
int wire = m_EventData->GetWireNbr(i);
double time = m_EventData->GetTime(i);
double etime = 0;
// look for matching trailing edge
for(unsigned int j = 0 ; j < size ; j++){
if(m_EventData->GetEdge(j)==1){
int edet = m_EventData->GetDetectorNbr(j);
int elayer = m_EventData->GetLayerNbr(j);
int ewire = m_EventData->GetWireNbr(j);
// same wire
if(wire==ewire && layer==elayer && det==edet){
etime = m_EventData->GetTime(j);
}
//if(etime<time)
// break;
//else
// etime=0;
}
}
// a valid wire must have an edge
if(etime && time){
Detector.push_back(det);
Layer.push_back(layer);
Wire.push_back(wire);
Time.push_back(time);
ToT.push_back(time-etime);
}
}
}
return;
}
///////////////////////////////////////////////////////////////////////////
void TSamuraiFDC2Physics::Clear(){
Detector.clear();
Layer.clear();
Wire.clear();
Time.clear();
ToT.clear();
ParticleDirection.clear();
MiddlePosition.clear();
}
///////////////////////////////////////////////////////////////////////////
//// Innherited from VDetector Class ////
///////////////////////////////////////////////////////////////////////////
void TSamuraiFDC2Physics::ReadConfiguration(NPL::InputParser parser){
}
///////////////////////////////////////////////////////////////////////////
void TSamuraiFDC2Physics::InitSpectra(){
//m_Spectra = new TSamuraiFDC2Spectra(m_NumberOfDetector);
}
///////////////////////////////////////////////////////////////////////////
void TSamuraiFDC2Physics::FillSpectra(){
// m_Spectra -> FillRawSpectra(m_EventData);
// m_Spectra -> FillPreTreatedSpectra(m_PreTreatedData);
// m_Spectra -> FillPhysicsSpectra(m_EventPhysics);
}
///////////////////////////////////////////////////////////////////////////
void TSamuraiFDC2Physics::CheckSpectra(){
// m_Spectra->CheckSpectra();
}
///////////////////////////////////////////////////////////////////////////
void TSamuraiFDC2Physics::ClearSpectra(){
// To be done
}
///////////////////////////////////////////////////////////////////////////
map< string , TH1*> TSamuraiFDC2Physics::GetSpectra() {
/* if(m_Spectra)
return m_Spectra->GetMapHisto();
else{
map< string , TH1*> empty;
return empty;
}*/
map< string , TH1*> empty;
return empty;
}
///////////////////////////////////////////////////////////////////////////
void TSamuraiFDC2Physics::WriteSpectra(){
// m_Spectra->WriteSpectra();
}
///////////////////////////////////////////////////////////////////////////
void TSamuraiFDC2Physics::AddParameterToCalibrationManager(){
}
///////////////////////////////////////////////////////////////////////////
void TSamuraiFDC2Physics::InitializeRootInputRaw(){
TChain* inputChain = RootInput::getInstance()->GetChain() ;
inputChain->SetBranchStatus( "SamuraiFDC2" , true );
// The following line is necessary only for system were the tree is splitted
// (older root version). The found argument silenced the Branches not found
// warning for non splitted tree.
if(inputChain->FindBranch("fDC_*"))
inputChain->SetBranchStatus( "fDC_*",true);
inputChain->SetBranchAddress( "SamuraiFDC2" , &m_EventData );
}
///////////////////////////////////////////////////////////////////////////
void TSamuraiFDC2Physics::InitializeRootInputPhysics(){
}
///////////////////////////////////////////////////////////////////////////
void TSamuraiFDC2Physics::InitializeRootOutput(){
TTree* outputTree = RootOutput::getInstance()->GetTree();
outputTree->Branch( "SamuraiFDC2" , "TSamuraiFDC2Physics" , &m_EventPhysics );
}
////////////////////////////////////////////////////////////////////////////////
// Construct Method to be pass to the DetectorFactory //
////////////////////////////////////////////////////////////////////////////////
NPL::VDetector* TSamuraiFDC2Physics::Construct(){
return (NPL::VDetector*) new TSamuraiFDC2Physics();
}
////////////////////////////////////////////////////////////////////////////////
// Registering the construct method to the factory //
////////////////////////////////////////////////////////////////////////////////
extern "C"{
class proxy_samuraiFDC2{
public:
proxy_samuraiFDC2(){
NPL::DetectorFactory::getInstance()->AddToken("SAMURAIFDC2","Samurai");
NPL::DetectorFactory::getInstance()->AddDetector("SAMURAIFDC2",TSamuraiFDC2Physics::Construct);
}
};
proxy_samuraiFDC2 p_samuraiFDC2;
}
#ifndef TSAMURAIFDC2PHYSICS_H
#define TSAMURAIFDC2PHYSICS_H
/*****************************************************************************
* 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 : October 2020 *
* Last update : *
*---------------------------------------------------------------------------*
* Decription: *
* This class hold SamuraiFDC2 treated data *
* *
*---------------------------------------------------------------------------*
* Comment: *
* *
* *
* *
*****************************************************************************/
// STL
#include <vector>
// NPL
#include "TSamuraiFDC2Data.h"
//#include "TSamuraiFDC2Spectra.h"
#include "NPCalibrationManager.h"
#include "NPVDetector.h"
#include "NPInputParser.h"
// ROOT
#include "TVector2.h"
#include "TVector3.h"
#include "TObject.h"
#include "TCanvas.h"
#include "TRandom3.h"
// Forward declaration
//class TSamuraiFDC2Spectra;
using namespace std ;
class TSamuraiFDC2Physics : public TObject, public NPL::VDetector{
public:
TSamuraiFDC2Physics();
~TSamuraiFDC2Physics() {};
public:
void Clear();
void Clear(const Option_t*) {};
public:
// Provide Physical Multiplicity
vector<int> Detector;
vector<int> Layer;
vector<int> Wire;
vector<double> Time;
vector<double> ToT;
// Computed variable
vector<TVector3> ParticleDirection;
vector<TVector3> MiddlePosition;
public:
// Projected position at given Z plan
TVector3 ProjectedPosition(double Z);
public: // Innherited from VDetector Class
// Read stream at ConfigFile to pick-up parameters of detector (Position,...) using Token
void ReadConfiguration(NPL::InputParser) ;
// Add Parameter to the CalibrationManger
void AddParameterToCalibrationManager() ;
// Activated associated Branches and link it to the private member DetectorData address
// In this method mother Branches (Detector) AND daughter leaf (fDetector_parameter) have to be activated
void InitializeRootInputRaw() ;
// Activated associated Branches and link it to the private member DetectorPhysics address
// In this method mother Branches (Detector) AND daughter leaf (parameter) have to be activated
void InitializeRootInputPhysics() ;
// Create associated branches and associated private member DetectorPhysics address
void InitializeRootOutput() ;
// This method is called at each event read from the Input Tree. Aime is to build treat Raw dat in order to extract physical parameter.
void BuildPhysicalEvent() ;
// Same as above, but only the simplest event and/or simple method are used (low multiplicity, faster algorythm but less efficient ...).
// This method aimed to be used for analysis performed during experiment, when speed is requiered.
// NB: This method can eventually be the same as BuildPhysicalEvent.
void BuildSimplePhysicalEvent() ;
// Same as above but for online analysis
void BuildOnlinePhysicalEvent() {BuildPhysicalEvent();};
// Those two method all to clear the Event Physics or Data
void ClearEventPhysics() {Clear();}
void ClearEventData() {m_EventData->Clear();}
// Method related to the TSpectra classes, aimed at providing a framework for online applications
// Instantiate the Spectra class and the histogramm throught it
void InitSpectra();
// Fill the spectra hold by the spectra class
void FillSpectra();
// Used for Online mainly, perform check on the histo and for example change their color if issues are found
void CheckSpectra();
// Used for Online only, clear all the spectra hold by the Spectra class
void ClearSpectra();
// Write Spectra to file
void WriteSpectra();
public: // Specific to SamuraiFDC2 Array
// Clear The PreTeated object
void ClearPreTreatedData() {m_PreTreatedData->Clear();}
// Remove bad channel, calibrate the data and apply threshold
void PreTreat();
// Retrieve raw and pre-treated data
TSamuraiFDC2Data* GetRawData() const {return m_EventData;}
TSamuraiFDC2Data* GetPreTreatedData() const {return m_PreTreatedData;}
private: // Root Input and Output tree classes
TSamuraiFDC2Data* m_EventData;//!
TSamuraiFDC2Data* m_PreTreatedData;//!
TSamuraiFDC2Physics* m_EventPhysics;//!
private: // Spectra Class
// TSamuraiFDC2Spectra* m_Spectra; // !
public: // Spectra Getter
map< string , TH1*> GetSpectra();
public: // Static constructor to be passed to the Detector Factory
static NPL::VDetector* Construct();
ClassDef(TSamuraiFDC2Physics,1) // SamuraiFDC2Physics structure
};
#endif
...@@ -8,4 +8,4 @@ StoppingPhysics 0 ...@@ -8,4 +8,4 @@ StoppingPhysics 0
OpticalPhysics 0 OpticalPhysics 0
HadronPhysicsINCLXX 0 HadronPhysicsINCLXX 0
HadronPhysicsQGSP_BIC_HP 0 HadronPhysicsQGSP_BIC_HP 0
Decay 0 Decay 1
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