FairRoot/PandaRoot
NeutronAnalysis_COSY.C
Go to the documentation of this file.
1 #include "TH1D.h"
2 #include "TH2D.h"
3 #include "TString.h"
4 #include "TFile.h"
5 #include "TTree.h"
6 #include "TStopwatch.h"
7 #include "TVector3.h"
8 #include "TROOT.h"
9 #include "TMath.h"
10 #include "TClonesArray.h"
11 #include "TSystem.h"
12 #include "TStyle.h"
13 //#include "/home/steinen/work/pandaroot/hypGe/PndHypGePoint.h"
14 
15 #include <stdio.h>
16 #include <iostream>
17 
18 
19 using namespace std;
20 int NeutronAnalysis_COSY(TString Filename_ext, Int_t StartEvent=0, Int_t NoOfEvents=0)
21 {
22  bool verbose = false;
23  // ----- Load libraries ------------------------------------------------
24  gROOT->Macro("$VMCWORKDIR/gconfig/rootlogon.C");
25  gSystem->Load("libHypGe");
26  // ----- Timer --------------------------------------------------------
27  TStopwatch timer;
28  timer.Start();
29  // ------------------------------------------------------------------------
30 
31 
32  //Files
33  //TString Path = getenv("SIMDATADIR");
34  TString Filename = "$SIMDATADIR/COSY/"+Filename_ext;
35  //TString Filename = "$SIMDATADIR/TripleBall40Offset10_urqmd_pbarC_1_5000Evts.root";
36  TString outfile= "$SIMDATADIR/COSY/Ana/ToCombine/" + Filename_ext;
37  outfile += "/Ana";
38  outfile += Filename_ext;
39  outfile.ReplaceAll(".root","");
40  outfile += "_Start_";
41  outfile += StartEvent;
42  outfile += "_nEvents_";
43  outfile += NoOfEvents;
44  outfile += ".root";
45  cout << outfile << endl;
46  TFile* InputFile = new TFile(Filename);
47  TFile* OutputFile = new TFile(outfile,"RECREATE");
48 
49  //getting simulation branches from input file
50  TTree *b=(TTree *) InputFile->Get("pndsim") ;
51  TClonesArray* hit_bar=new TClonesArray("PndHypGePoint");
52  b->SetBranchAddress("HypGePoint",&hit_bar);//Branch names
53  TClonesArray* mc_bar=new TClonesArray("PndMCTrack");
54  b->SetBranchAddress("MCTrack",&mc_bar);//Branch names
55 
56  //declaring histograms
57  //histogram with polar angle of detector interactions of primary neutrons
58  TH1D* hNHits = new TH1D("hNHits","Polar angle of primary neutrons interacting with the crystals",180, 90,180);
59  hNHits->SetXTitle("#Theta [#circ]");
60  hNHits->SetYTitle("Counts / 0.5 #circ");
61  //histogram with the polar momentum of primary neutrons interacting with the crystals
62  TH1D* hNHits_p = new TH1D("hNHits_p","Direction of polar momentum of primary neutrons interacting with the crystals",180, 90,180);
63  hNHits_p->SetXTitle("#Theta [#circ]");
64  hNHits_p->SetYTitle("Counts / 0.5 #circ");
65  //histogram with the polar momentum of primary neutrons in the solid angle of the crystal (not necessarily interacting with the crystal)
66  TH1D* hNAll = new TH1D("hNAll","All Neutrons",360, 0,180);
67  hNAll->SetXTitle("#Theta [#circ]");
68  hNAll->SetYTitle("Counts / 0.5 #circ");
69 
70  //histogram with the polar momentum of primary neutrons in the solid angle of the crystal (not necessarily interacting with the crystal)
71  TH1D* hPAll = new TH1D("hPAll","All Protons",360, 0,180);
72  hPAll->SetXTitle("#Theta [#circ]");
73  hPAll->SetYTitle("Counts / 0.5 #circ");
74 
75  //histogram with E_kin of primary neutrons
76  TH1D* hNMom = new TH1D("hNMom","E_{kin} of neutrons",1000, 0,0.5);
77  hNMom->SetXTitle("E_{kin} of neutrons [GeV]");
78  hNMom->SetYTitle("Counts / 0.5 MeV");
79  //histogram to see which detector is hit
80  TH1D* hCrystalHit = new TH1D("hCrystalHit","Hits per Crystal",1700,1,1700);
81  hCrystalHit->SetXTitle("Crystal number");
82  hCrystalHit->SetYTitle("Counts per crystal");
83  //histogram to see the energy deposited by the neutron
84  TH1D* hNeutronEnergyDeposit = new TH1D("hNeutronEnergyDeposit","Energy deposited by Neutrons per Hit",500,0,0.005);
85  hNeutronEnergyDeposit->SetXTitle("Energy [GeV]");
86  hNeutronEnergyDeposit->SetYTitle("Counts / 10 keV");
87 
88  //histogram to visualize the x-y-distribution of neutrons
89  TH2D* hNeutronXYDistribution = new TH2D("hNeutronXYDistribution","Distribution of Neutrons",40, -40,40,40,-40,40);
90  hNeutronXYDistribution->SetXTitle("X-postion [cm]");
91  hNeutronXYDistribution->SetYTitle("Y-postion [cm]");
92  gStyle->SetPalette(1);
93 
94 
95 
96  Int_t *ActualTrackID ;
97  //Int_t nEvents = b->GetEntriesFast();
98  //cout<< "Number of Simulated Events: "<<nEvents<<endl;
99 
100  Int_t Hits=0;
101  Int_t DetID=-10;
102  Int_t TrackID=-10;
103  Int_t iNeutrons = 0;
104 
105 
106 
107  TVector3 NHit;
108  TVector3 NMom;
109  TVector3 NAllMom;
110  TVector3 StartVertex;
111 
112 
113  Int_t EndEvent;
114  if (NoOfEvents ==0)
115  EndEvent = b->GetEntriesFast();
116  else
117  EndEvent = StartEvent + NoOfEvents;
118 
119  //event loop
120  for (Int_t k=StartEvent; k< EndEvent; k++)
121  //for (Int_t k=569494; k<569498; k++)
122  {
123  b->GetEntry(k);
124  if (!((k*100)% NoOfEvents)) // mark every % of proceeded events
125  {
126  cout << "Event number :\t" << k << endl;
127  }
128  //if(verbose) cout<<"Event No "<<j<<endl;
129 
130  DetID=-10;
131  TrackID=-10;
132 
133  // loop over the hits of an event
134  for (Int_t i=0; i<hit_bar->GetEntriesFast(); i++)
135  {
136  //if (i==0)
137  //Hits++;
138  //cout <<"Hit "<< i << endl;
139 
140  PndHypGePoint *hitgam=(PndHypGePoint*)hit_bar->At(i);
141  PndMCTrack *mcgam = (PndMCTrack*)mc_bar->At(hitgam->GetTrackID());
142 
143  //NMom = mcgam->GetMomentum();
144  if (hitgam->GetpdgCode() == 2112)
145  {
146  StartVertex = mcgam->GetStartVertex();
147  cout <<"Event " <<k<< " \t\tStartVertex Radius: " <<StartVertex.Mag()<< endl;
148 
149  if (StartVertex.Mag() <20) // prevent counting of neutrons created by nuclear reactions inside the detector 20 (cm) is random value, must be smaller than the distance from (0,0,0) to detector (here 30 cm)
150  {
151  if (!(DetID == hitgam->GetDetectorID() && TrackID==hitgam->GetTrackID())) // used to prevent counting of neutrons that reenter a crystal (neutrons left crystal through the central hole skewly and reentered the same crystal afterwards) and of neutrons that enter multiple crystals by scattering
152  {
153  DetID = hitgam->GetDetectorID();
154  TrackID = hitgam->GetTrackID();
155  //cout << "Event :\t" << k <<"\tHit : \t" << i << "\tDetector :\t" << DetID << "\tTrackID :\t"<< TrackID<<"\tMother :\t"<< mcgam->GetMotherID()<<endl;
156  iNeutrons++;
157  hCrystalHit->Fill(hitgam->GetDetectorID()); //fill histogram to see which detector is hit
158 
159 
160 
161  // fill vector with hit coordinates to get theta of the hit
162  NHit.SetX(hitgam->GetX());
163  NHit.SetY(hitgam->GetY());
164  NHit.SetZ(hitgam->GetZ()); // no shift like in PANDA!!!!
165  hNHits->Fill(180/TMath::Pi()*NHit.Theta()); //fill histogram with polar angle of detector interactions of primary neutrons
166  cout << "Event " << k <<"Crystal" << hitgam->GetDetectorID()<<" Theta " <<180/TMath::Pi()*NHit.Theta()<<endl;
167 
168  TVector3 NHit_p;
169  // fill vector with hit momentum to get theta of the hit momentum
170  NHit_p.SetX(hitgam->GetPx());
171  NHit_p.SetY(hitgam->GetPy());
172  NHit_p.SetZ(hitgam->GetPz());
173  hNHits_p->Fill(180/TMath::Pi()*NHit_p.Theta()); //fill histogram with the polar momentum of primary neutrons interacting with the crystals
174 
175  hNeutronXYDistribution->Fill(hitgam->GetX(),hitgam->GetY());
176  //hNMom->Fill(1./2.* NMom.Mag2()/0.939); //fill histogram with E_kin of primary neutrons
177  hNeutronEnergyDeposit->Fill(hitgam->GetEnergyLoss());
178 
179  }
180  }
181  }
182  }
183 
184  // loop over all entries in mctrack
185  //cout << "mcbar length :\t" <<mc_bar->GetEntriesFast() << endl;
186 
187 
188  for (Int_t m = 0; m<mc_bar->GetEntriesFast();m++)
189  {
190  //cout << "mcbar length :\t" <<mc_bar->GetEntriesFast() << endl;
191  PndMCTrack *mcgam = (PndMCTrack*)mc_bar->At(m);
192 
193  StartVertex = mcgam->GetStartVertex();
194  if (StartVertex.Mag())
195  {
196  NAllMom = mcgam->GetMomentum();
197  if (mcgam->GetPdgCode() == 2112 )
198  {
199  //cout << "event :\t"<<k<<"\tmcpdg :\t" << mcgam->GetPdgCode() << "\tmcmotherid :\t"<< mcgam->GetMotherID() <<endl;
200 
201  hNAll->Fill(180/TMath::Pi()*NAllMom.Theta()); //fill histogram with the polar momentum of primary neutrons in the solid angle of the crystal (not necessarily interacting with the crystal)
202  }
203  if (mcgam->GetPdgCode() == 2212 && mcgam->GetMotherID() != -1 )
204  {
205  //cout << "event :\t"<<k<<"\tmcpdg :\t" << mcgam->GetPdgCode() << "\tmcmotherid :\t"<< mcgam->GetMotherID() <<endl;
206 
207  hPAll->Fill(180/TMath::Pi()*NAllMom.Theta()); //fill histogram with the polar momentum of primary neutrons in the solid angle of the crystal (not necessarily interacting with the crystal)
208  }
209  }
210  }
211 
212 
213  }// end for k (events)
214 
215 
216  hNeutronXYDistribution->DrawCopy("colz");
217  //hNAll->DrawCopy();
218  //hNHits->SetLineColor(kRed);
219  //hNHits->DrawCopy("same");
220 
221  //Analysis of spectrum
222  hNHits->Write();
223  hNHits_p->Write();
224  hNAll->Write();
225  hPAll->Write();
226  hNMom->Write();
227  hCrystalHit->Write();
228  hNeutronEnergyDeposit->Write();
229  hNeutronXYDistribution->Write();
230 
231 
232 
233  //hNeutronEnergyDeposit->DrawCopy();
234 
235  OutputFile->Close();
236 
237  //hNeutronEnergyDeposit->Draw();
238  //hNHits->Draw();
239 
240 
241  // ----- Finish -------------------------------------------------------
242  timer.Stop();
243  Double_t rtime = timer.RealTime();
244  Double_t ctime = timer.CpuTime();
245  cout << endl << endl;
246  cout << "Macro finished succesfully." << endl;
247 
248  cout << "Real time " << rtime << " s, CPU time " << ctime << " s" << endl;
249  cout << endl;
250 
251  cout << "Hits\t" << Hits << endl;
252  // ------------------------------------------------------------------------
253  return iNeutrons;
254 }
Int_t i
Definition: run_full.C:25
__m128 m
Definition: P4_F32vec4.h:28
TTree * b
Int_t GetPdgCode() const
Definition: PndMCTrack.h:73
#define verbose
Double_t GetPz() const
Definition: PndHypGePoint.h:59
Int_t GetTrackID() const
Definition: PndHypGePoint.h:51
TVector3 GetMomentum() const
Definition: PndMCTrack.h:78
int NeutronAnalysis_COSY(TString Filename_ext, Int_t StartEvent=0, Int_t NoOfEvents=0)
TClonesArray * hit_bar
Double_t
Double_t GetPy() const
Definition: PndHypGePoint.h:58
TStopwatch timer
Definition: hit_dirc.C:51
Double_t GetPx() const
Definition: PndHypGePoint.h:57
Double_t ctime
Definition: hit_dirc.C:114
Double_t GetZ() const
Definition: PndHypGePoint.h:56
Double_t GetX() const
Definition: PndHypGePoint.h:54
Double_t GetY() const
Definition: PndHypGePoint.h:55
TClonesArray * mc_bar
TVector3 GetStartVertex() const
Definition: PndMCTrack.h:76
Int_t GetMotherID() const
Definition: PndMCTrack.h:74
Double_t rtime
Definition: hit_dirc.C:113
Double_t Pi
Int_t GetDetectorID() const
Definition: PndHypGePoint.h:53
Double_t GetEnergyLoss() const
Definition: PndHypGePoint.h:62
TString outfile
Double_t GetpdgCode() const
Definition: PndHypGePoint.h:63