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

* Fixing angle on Tiara Hyball

        - The ring center was located on the upper edge of the strip
          instead of the actual center
parent 4ef7f503
No related branches found
No related tags found
No related merge requests found
...@@ -632,7 +632,7 @@ void TTiaraHyballPhysics::AddWedgeDetector( double R,double Phi,double Z){ ...@@ -632,7 +632,7 @@ void TTiaraHyballPhysics::AddWedgeDetector( double R,double Phi,double Z){
for(int b = 0 ; b < Wedge_Sector_NumberOfStrip ; b++){ for(int b = 0 ; b < Wedge_Sector_NumberOfStrip ; b++){
StripCenter = Strip_1_1; StripCenter = Strip_1_1;
StripCenter.SetY(Wedge_R_Max-f*StripPitchRing); StripCenter.SetY(Wedge_R_Max-f*StripPitchRing-0.5*StripPitchRing);
StripCenter.SetZ(Z); StripCenter.SetZ(Z);
StripCenter.RotateZ(Phi+Wedge_Phi_Min+b*StripPitchSector); StripCenter.RotateZ(Phi+Wedge_Phi_Min+b*StripPitchSector);
lineX.push_back( StripCenter.X() ); lineX.push_back( StripCenter.X() );
......
...@@ -48,9 +48,9 @@ using namespace CLHEP; ...@@ -48,9 +48,9 @@ using namespace CLHEP;
namespace TIARA{ namespace TIARA{
// Energy and time Resolution // Energy and time Resolution
const G4double ResoTime = 0 ; const G4double ResoTime = 0 ;
const G4double ResoEnergyInnerBarrel = 0.058*MeV ;// = 136keV FWHM const G4double ResoEnergyInnerBarrel = 0.017*MeV ;// = 136keV FWHM
const G4double ResoEnergyOuterBarrel = 0.058*MeV ;// = 136keV FWHM const G4double ResoEnergyOuterBarrel = 0.017*MeV ;// = 136keV FWHM
const G4double ResoEnergyHyball = 0.029*MeV ;// = 70keV FWHM const G4double ResoEnergyHyball = 0.017*MeV ;// = 70keV FWHM
const G4double EnergyThreshold = 200*keV; const G4double EnergyThreshold = 200*keV;
......
...@@ -74,6 +74,9 @@ void Analysis::Init(){ ...@@ -74,6 +74,9 @@ void Analysis::Init(){
Si_E_TB = 0 ; Si_E_TB = 0 ;
Energy = 0; Energy = 0;
Original_ELab=0;
Original_ThetaLab=0;
XTarget = 0; XTarget = 0;
YTarget =0; YTarget =0;
BeamDirection = TVector3(0,0,1); BeamDirection = TVector3(0,0,1);
...@@ -86,7 +89,9 @@ void Analysis::Init(){ ...@@ -86,7 +89,9 @@ void Analysis::Init(){
void Analysis::TreatEvent(){ void Analysis::TreatEvent(){
// Reinitiate calculated variable // Reinitiate calculated variable
ReInitValue(); ReInitValue();
Original_ELab = Initial->GetKineticEnergy(0);
Original_ThetaLab = Initial->GetParticleDirection(0).Angle(Initial->GetBeamDirection())/deg;
////////////////////////////////////////// LOOP on TiaraHyball + SSSD Hit ////////////////////////////////////////// ////////////////////////////////////////// LOOP on TiaraHyball + SSSD Hit //////////////////////////////////////////
for(unsigned int countTiaraHyball = 0 ; countTiaraHyball < TH->Strip_E.size() ; countTiaraHyball++){ for(unsigned int countTiaraHyball = 0 ; countTiaraHyball < TH->Strip_E.size() ; countTiaraHyball++){
...@@ -102,6 +107,8 @@ void Analysis::TreatEvent(){ ...@@ -102,6 +107,8 @@ void Analysis::TreatEvent(){
if(XTarget>-1000 && YTarget>-1000){ if(XTarget>-1000 && YTarget>-1000){
TVector3 BeamImpact(XTarget,YTarget,0); TVector3 BeamImpact(XTarget,YTarget,0);
TVector3 HitDirection = TH -> GetRandomisedPositionOfInteraction(countTiaraHyball) - BeamImpact ; TVector3 HitDirection = TH -> GetRandomisedPositionOfInteraction(countTiaraHyball) - BeamImpact ;
// TVector3 HitDirection = TH -> GetPositionOfInteraction(countTiaraHyball) - BeamImpact ;
ThetaLab = HitDirection.Angle( BeamDirection ); ThetaLab = HitDirection.Angle( BeamDirection );
ThetaTHSurface = HitDirection.Angle(TVector3(0,0,-1) ); ThetaTHSurface = HitDirection.Angle(TVector3(0,0,-1) );
ThetaNormalTarget = HitDirection.Angle( TVector3(0,0,1) ) ; ThetaNormalTarget = HitDirection.Angle( TVector3(0,0,1) ) ;
...@@ -235,6 +242,8 @@ void Analysis::ReInitValue(){ ...@@ -235,6 +242,8 @@ void Analysis::ReInitValue(){
ELab = -1000; ELab = -1000;
ThetaLab = -1000; ThetaLab = -1000;
ThetaCM = -1000; ThetaCM = -1000;
Original_ELab = -1000;
Original_ThetaLab = -1000;
} }
///////////////////////////////////////////////////////////////////////////// /////////////////////////////////////////////////////////////////////////////
//////////////////////////////////////////////////////////////////////////// ////////////////////////////////////////////////////////////////////////////
...@@ -246,6 +255,8 @@ void Analysis::InitOutputBranch() { ...@@ -246,6 +255,8 @@ void Analysis::InitOutputBranch() {
RootOutput::getInstance()->GetTree()->Branch("TiaraImpactMatrixX",&TiaraIMX,"TiaraImpactMatrixX/D"); RootOutput::getInstance()->GetTree()->Branch("TiaraImpactMatrixX",&TiaraIMX,"TiaraImpactMatrixX/D");
RootOutput::getInstance()->GetTree()->Branch("TiaraImpactMatrixY",&TiaraIMY,"TiaraImpactMatrixY/D"); RootOutput::getInstance()->GetTree()->Branch("TiaraImpactMatrixY",&TiaraIMY,"TiaraImpactMatrixY/D");
RootOutput::getInstance()->GetTree()->Branch("TiaraImpactMatrixZ",&TiaraIMZ,"TiaraImpactMatrixZ/D"); RootOutput::getInstance()->GetTree()->Branch("TiaraImpactMatrixZ",&TiaraIMZ,"TiaraImpactMatrixZ/D");
RootOutput::getInstance()->GetTree()->Branch("Original_ELab",&Original_ELab,"Original_ELab/D");
RootOutput::getInstance()->GetTree()->Branch("Original_ThetaLab",&Original_ThetaLab,"Original_ThetaLab/D");
} }
///////////////////////////////////////////////////////////////////////////// /////////////////////////////////////////////////////////////////////////////
void Analysis::InitInputBranch(){ void Analysis::InitInputBranch(){
......
...@@ -83,6 +83,9 @@ class Analysis: public NPL::VAnalysis{ ...@@ -83,6 +83,9 @@ class Analysis: public NPL::VAnalysis{
double Si_E_TB ; double Si_E_TB ;
double Energy ; double Energy ;
double TargetThickness ; double TargetThickness ;
double Original_ThetaLab;
double Original_ELab;
double XTarget ; double XTarget ;
double YTarget ; double YTarget ;
......
void Results(){ void Results(){
TFile* File = new TFile("../../Outputs/Analysis/PhysicsTree.root"); TFile* File = new TFile("../../Outputs/Analysis/PhysicsTree.root");
TTree* Tree = (TTree*) File->FindObjectAny("PhysicsTree"); TTree* Tree = (TTree*) File->FindObjectAny("PhysicsTree");
...@@ -7,5 +5,20 @@ void Results(){ ...@@ -7,5 +5,20 @@ void Results(){
Tree->Draw("ELab:ThetaLab>>h2(1000,0,180,1000,0,20)","ELab>0","colz"); Tree->Draw("ELab:ThetaLab>>h2(1000,0,180,1000,0,20)","ELab>0","colz");
new TCanvas(); new TCanvas();
Tree->Draw("Ex>>h4(1000)","ELab>0 && Ex>-1000"); Tree->Draw("Ex>>h4(1000)","ELab>0 && Ex>-1000");
TCanvas* c = new TCanvas("CAS","Control Analyse-Simu");
c->Divide(2,2);
c->cd(1);
Tree->Draw("ELab:Original_ELab>>hcE(1000,0,30,1000,0,30)","","colz");
TH2* h = (TH2*) gDirectory->FindObjectAny("hcE");
h->GetXaxis()->SetTitle("Simulated_{} E_{lab} ");
h->GetYaxis()->SetTitle("Reconstructed_{} E_{lab} ");
TLine* lineE = new TLine (0,0,30,30);
lineE->Draw();
c->cd(2);
Tree->Draw("ThetaLab:Original_ThetaLab>>hcT(1000,0,180,1000,0,180)","","colz");
h = (TH2*) gDirectory->FindObjectAny("hcT");
h->GetXaxis()->SetTitle("Simulated #theta_{lab} ");
h->GetYaxis()->SetTitle("Reconstructed #theta_{lab} ");
TLine* lineT = new TLine (0,0,180,180);
lineT->Draw();
} }
...@@ -4,7 +4,7 @@ GeneralTarget ...@@ -4,7 +4,7 @@ GeneralTarget
%0.5mg/cm2 %0.5mg/cm2
Target Target
THICKNESS= 0.48 THICKNESS= 0.001
RADIUS= 5 RADIUS= 5
MATERIAL= CD2 MATERIAL= CD2
ANGLE= 0 ANGLE= 0
......
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