/[PAMELA software]/PamVMC/aux/PamVMCDigitizer/vmc_dig.cpp
ViewVC logotype

Contents of /PamVMC/aux/PamVMCDigitizer/vmc_dig.cpp

Parent Directory Parent Directory | Revision Log Revision Log


Revision 1.2 - (show annotations) (download)
Mon Jul 13 15:05:46 2009 UTC (15 years, 4 months ago) by pam-rm2
Branch: MAIN
CVS Tags: HEAD
Changes since 1.1: +1 -1 lines
Error occurred while calculating annotation data.
----------------------------------------------------------------------

1 // C/C++ headers
2 //
3 #include <fstream>
4 #include <string.h>
5 #include <iostream>
6 #include <stdio.h>
7 #include <ctype.h>
8 //
9 // ROOT headers
10 //
11
12 #include <TTree.h>
13
14 #include <TSystem.h>
15 #include <TStopwatch.h>
16 #include <TSystemDirectory.h>
17 #include <TGraphAsymmErrors.h>
18 #include <TString.h>
19 #include <TFile.h>
20 #include <TClass.h>
21 #include <TGraph.h>
22 #include <TLine.h>
23 #include <TH1F.h>
24 #include <TMath.h>
25 #include <TRandom3.h>
26 #include <TF1.h>
27 //
28 // RunInfo header
29 //
30 #include <RunInfo.h>
31 #include <PamLevel2.h>
32 #include <TClonesArray.h>
33
34 //PamVMC headers
35 #include "PamVMCDigMgr.h"
36 #include "PamVMCTrkHit.h"
37 #include "PamVMCTrkF77GpsSpe.h"
38 #include "PamVMCPrimaryGenerator.h"
39 using namespace std;
40 //This program is standalone Digitizer for PamVMC
41
42
43 void Usage ();
44
45 TClonesArray* FindColl(TTree*t, TString brname, TString type = "PamVMCDetectorHit");
46 void RegisterPDG();
47 //Here user should define all nuclei and particles that are not in PDG tabe of ROOT
48 //otherwise application will crash on TOF digitization whrer charge information needed;
49
50 void UnCompressCalo(TClonesArray* cal_comp, TClonesArray & cal_u);
51
52 int main(int argc, char * argv[] ){
53
54 if(argc!=4){
55 Usage();
56 return 0;
57 }
58
59 TString DataPATH = gSystem->Getenv("PWD");
60 TString inputfilename(argv[1]);
61 if(!inputfilename.Contains(".root")){
62 cerr<<"Error! Not root file on input!"<<endl;
63 return 0;
64 }
65
66 TString outputfilename(argv[2]);
67
68 Int_t seed = atoi(argv[3]);
69 TRandom3* rnd = new TRandom3(seed);
70 gRandom = rnd;
71
72 TStopwatch timer;
73 timer.Start();
74 TFile* rfin = new TFile(inputfilename, "UPDATE");
75 if(!rfin){
76 cerr<<"Error! Can't read root file!"<<endl;
77 return 0;
78 }
79
80 TTree::SetMaxTreeSize(1000*Long64_t(2000000000));
81
82 TTree * hits_tree = (TTree*)rfin->Get("hits");
83 if(!hits_tree){
84 cerr<<"Tree hits are not found in file "<<rfin->GetName()<<endl;
85 TString treename = rfin->GetName();
86 treename.Replace(treename.Index(".root",5),5,"");
87 hits_tree = (TTree*)rfin->Get(treename);
88 if(!hits_tree) return 0;
89 }
90
91
92
93 PamVMCDigMgr* d_mgr = PamVMCDigMgr::Instance();
94 d_mgr->Reset();
95 //AC CARD: C1D1, C2D1, CAT: TOP1, CAS: SID1
96 PamVMCAcDig* ac = new PamVMCAcDig();
97 ac->RegisterCollections("C1D1", FindColl(hits_tree,"C1D1"));
98 ac->RegisterCollections("C2D1", FindColl(hits_tree,"C2D1"));
99 ac->RegisterCollections("TOP1", FindColl(hits_tree,"TOP1"));
100 ac->RegisterCollections("SID1", FindColl(hits_tree,"SID1"));
101 d_mgr->SetDIG("AC",ac);
102 //TOF: S11Y, S12X, S21X, S22Y, S31Y, S32X
103 PamVMCTofDig* tof = new PamVMCTofDig();
104 tof->RegisterCollections("S11Y", FindColl(hits_tree,"S11Y"));
105 tof->RegisterCollections("S12X", FindColl(hits_tree,"S12X"));
106 tof->RegisterCollections("S21X", FindColl(hits_tree,"S21X"));
107 tof->RegisterCollections("S22Y", FindColl(hits_tree,"S22Y"));
108 tof->RegisterCollections("S31Y", FindColl(hits_tree,"S31Y"));
109 tof->RegisterCollections("S32X", FindColl(hits_tree,"S32X"));
110 d_mgr->SetDIG("Tof",tof);
111
112
113 //Calorimeter CAST
114 #define NSTRIPS 4224
115 PamVMCCaloDig* cal = new PamVMCCaloDig();
116 TClonesArray *cal_u = new TClonesArray("PamVMCDetectorHit",NSTRIPS); //uncompressed array
117 TClonesArray* cal_comp = FindColl(hits_tree,"CAST"); // compressed array from file
118 cal->RegisterCollections("CAST", cal_u);
119 d_mgr->SetDIG("CAST",cal);
120
121 //Tracker GPSSPE
122 PamVMCTrkDig* trk = new PamVMCTrkDig();
123 TBranch *gps = hits_tree->GetBranch("GPSSPE");
124 trk->RegisterCollections("TSPA", FindColl(hits_tree,"TSPA"));
125 pGPSSPEHits * hits = new pGPSSPEHits();
126 if(gps){
127 hits_tree->SetBranchAddress("GPSSPE",&hits);
128 hits_tree->SetBranchStatus("GPSSPE",1);
129 cout<<"TRKCOLL="<<hits<<endl;
130 trk->RegisterTrkData(hits);
131 }
132 d_mgr->SetDIG("TSPA",trk);
133
134 //S4
135 PamVMCS4Dig* s4 = new PamVMCS4Dig();
136 s4->RegisterCollections("S4", FindColl(hits_tree,"S4"));
137 d_mgr->SetDIG("S4",s4);
138 //ND: NDTI
139 PamVMCNDDig* nd = new PamVMCNDDig();
140 nd->RegisterCollections("NDTI", FindColl(hits_tree,"NDTI"));
141 d_mgr->SetDIG("NDTI",nd);
142
143 //Primaries
144 TClonesArray *prims = FindColl(hits_tree,"PRIM","PamVMCPrimary");
145
146 //Initialization
147 PamVMCRawMgr* r_mgr = PamVMCRawMgr::Instance();
148 r_mgr->CreateOutFile(outputfilename.Data());
149 d_mgr->Initialize(gRandom);
150 r_mgr->WriteToFile();
151
152 Int_t nev = hits_tree->GetEntries();
153 cout<<"TOTAL NUMBER OF EVENTS:"<<nev<<endl;
154 RegisterPDG();
155 //Event loop
156 d_mgr->PrintCollections();
157 for(Int_t i = 0; i<nev; i++){
158 hits_tree->GetEntry(i);
159 cout<<"--> Digitizing evnt"<<endl;
160 UnCompressCalo(cal_comp, *cal_u);
161 d_mgr->Digitize(i,((PamVMCPrimary*)prims->At(0))->fPDG);
162 cout<<"--> Writing evnt"<<endl;
163 r_mgr->WriteEvent();
164 }
165
166 d_mgr->FinishRun();
167 r_mgr->FinishRun();
168 d_mgr->PrintCollections();
169 timer.Stop();
170 gDirectory->Delete("tmp.root");
171 rfin->Close();
172 return 0;
173 }
174
175 TClonesArray* FindColl(TTree*t, TString brname, TString type){
176 if(t->GetBranch(brname)){
177 TClonesArray* coll = 0;//new TClonesArray(type,4224);
178 t->GetBranch(brname)->SetAddress(&coll);
179 cout<<"coll="<<coll<<endl;
180 return coll;
181 }
182 return 0;
183 }
184
185 void RegisterPDG(){
186 TDatabasePDG::Instance()->AddParticle("He-4","He-4",3.7274,kTRUE,0., 6. /* |e|/3 charge*/,"Ion",1000020040);
187 }
188
189 void UnCompressCalo(TClonesArray* cal_comp, TClonesArray & cal_u){
190 if(!cal_comp) return;
191 PamVMCDetectorHit* hit = NULL;
192 PamVMCDetectorHit* newhit = NULL;
193 cal_u.Clear("C");
194 for(Int_t j = 0; j<cal_comp->GetEntries(); j++){
195 hit=((PamVMCDetectorHit*)cal_comp->At(j));
196 newhit = (PamVMCDetectorHit*)cal_u.New(hit->GetPOS());
197 *newhit = *hit;
198 }
199 }
200
201 void Usage(){
202 cout<<"NAME"<<endl;
203 cout<<" vmc_dig - a standalone digitizer for PamVMC application"<<endl;
204 cout<<"SYNOPSIS"<<endl;
205 cout<<" vmc_dig INPROOTFILE OUTRAWFILE SEED"<<endl;
206 cout<<"DISCRIPTION"<<endl;
207 cout<<" This tool convert hits collections of PamVMC in raw format "<<endl;
208 cout<<" There are only 3 parameters: input root/output raw filenames and integer random seed."<<endl;
209 cout<<" If seed is 0, then TRandom3(0) object will be created."<<endl;
210 }

  ViewVC Help
Powered by ViewVC 1.1.23