/[PAMELA software]/trieste/pamVMC/src/PamVMCApplication.cxx
ViewVC logotype

Contents of /trieste/pamVMC/src/PamVMCApplication.cxx

Parent Directory Parent Directory | Revision Log Revision Log


Revision 1.1 - (show annotations) (download)
Wed Mar 4 12:51:18 2009 UTC (15 years, 11 months ago) by pamelats
Branch: MAIN
Branch point for: pamVMC
Initial revision

1 // $Id: PamVMCApplication.cxx,v 1.0 2007/06/01
2 //
3 // Geant4 PAMELA VMC
4 #include <cstdlib>
5 #include <string>
6 #include <iostream>
7 #include <TROOT.h>
8 #include <TInterpreter.h>
9 #include <TVirtualMC.h>
10 #include <Riostream.h>
11 #include <TGeoManager.h>
12 #include <TVirtualGeoTrack.h>
13 #include <TLorentzVector.h>
14 #include "PamVMCApplication.h"
15 #include "PamVMCStack.h"
16 #include "PamVMCPrimaryGenerator.h"
17
18 #include "PamVMCFieldMgr.h"
19
20 #include "PamVMCSDMgr.h"
21
22 #include "PamVMCVolCrossMgr.h"
23
24 #include "PamVMCDetPamela.h"
25
26 //#include "TFluka.h"
27
28
29
30 ClassImp(PamVMCApplication)
31
32 PamVMCApplication::PamVMCApplication(const char *name, const char *title)
33 : TVirtualMCApplication(name,title),
34 fEventNo(0),
35 fVerbose(2),
36 fStack(0),
37 fdetector(0),
38 fPrimaryGenerator(0),
39 fRootManager("pamtest",kWrite)
40 {
41 // Standard constructor
42
43 // Create a user stack
44 fStack = new PamVMCStack(1000);
45
46
47 // Create a primary generator
48 fPrimaryGenerator = new PamVMCPrimaryGenerator(fStack);
49
50
51 fdetector = new PamVMCDetPamela();
52
53 PamVMCFieldMgr::Instance()->LoadField();
54 }
55
56
57 PamVMCApplication::PamVMCApplication()
58 : TVirtualMCApplication(),
59 fEventNo(0),
60 fVerbose(0),
61 fStack(0),
62 fdetector(0),
63 fPrimaryGenerator(0),
64 fRootManager()
65 {
66 // Default constructor
67 }
68
69 PamVMCApplication::~PamVMCApplication()
70 {
71 // Destructor
72
73 delete fStack;
74 delete fdetector;
75 delete fPrimaryGenerator;
76 delete gMC;
77 gMC = 0;
78 }
79
80
81 void PamVMCApplication::InitMC(const char* setup)
82 {
83 // Initialize MC.
84
85 fVerbose.InitMC();
86
87 gROOT->LoadMacro(setup);
88 gInterpreter->ProcessLine("Config()");
89 gMC->SetStack(fStack);
90 std::cout<<"!!! BEFORE INIT GMC"<<std::endl;
91
92 gMC->Init();
93 cout<<"INIT MC"<<endl;
94 gMC->BuildPhysics();
95 cout<<"PHYSICS BUILDED"<<endl;
96
97
98 fdetector->InitMC();
99
100 // fdetector->Print();
101
102 PamVMCSDMgr::Instance()->Register();
103
104 PamVMCDigMgr::Instance()->LoadCalib();
105 PamVMCRawMgr::Instance()->WriteToFile();
106
107 fPrimaryGenerator->Register();
108
109 std::cout<<"LIST OF ACTIVE SD'S:"<<std::endl;
110
111 PamVMCSDMgr::Instance()->Print();
112 }
113
114 void PamVMCApplication::RunMC(Int_t nofEvents)
115 {
116 // MC run.
117
118 fVerbose.RunMC(nofEvents);
119
120 gMC->ProcessRun(nofEvents);
121
122 fVerbose.FinishRun();
123 }
124
125 void PamVMCApplication::AddIons(){
126
127 fVerbose.AddIons();
128 gMC->DefineIon("He4",2,4,1,0.);
129 gMC->DefineIon("Boron",5,10,1,0.);
130 gMC->DefineIon("Carbon",6,12,1,0.);
131 }
132
133
134 void PamVMCApplication::ConstructGeometry()
135 {
136 // Construct geometry using detector contruction class
137
138 fVerbose.ConstructGeometry();
139 new TGeoManager("PAM_GEOM","PAMELA GEOMETRY");
140 fdetector->ConstructGeometry();
141 //TGeoManager::Import("/home/pamela/PamCAL210508/gp2tgeo.root");
142 gGeoManager->CloseGeometry();
143 gMC->SetRootGeometry();
144
145 cout<<"!!CONSTRUCT GEOMETRY!!"<<endl;
146 fdetector->ConstructCuts();
147
148 }
149
150 void PamVMCApplication::InitGeometry()
151 {
152 // Initialize geometry
153
154 cout<<"!!SETUP CUTS"<<endl;
155 //fdetector->InitGeometry();
156
157 fVerbose.InitGeometry();
158 }
159
160
161 void PamVMCApplication::GeneratePrimaries()
162 {
163 // Fill the user stack (derived from TVirtualMCStack) with primary particles.
164
165 fPrimaryGenerator->GeneratePrimary();
166 //fPrimaryGenerator->Print();
167 }
168
169 void PamVMCApplication::BeginEvent()
170 {
171 // User actions at beginning of event
172 PamVMCSDMgr::Instance()->PreEvent();
173 fEventNo++;
174 }
175
176 void PamVMCApplication::BeginPrimary()
177 {
178 // User actions at beginning of a primary track
179
180 fVerbose.BeginPrimary();
181 }
182
183 void PamVMCApplication::PreTrack()
184 {
185 // User actions at beginning of each track
186
187 }
188
189 void PamVMCApplication::Stepping()
190 {
191 if (fStack->GetCurrentTrack()->IsPrimary()){
192 PamVMCVolCross* tmp=(PamVMCVolCross*)PamVMCVolCrossMgr::Instance()->GetSD(gMC->CurrentVolName());
193 if (tmp){
194 tmp->CrossVol();
195 }
196 }
197 // User actions at each step
198 fVerbose.Stepping();
199 PamVMCDetectorSD* temp=(PamVMCDetectorSD*) PamVMCSDMgr::Instance()->GetSD(gMC->CurrentVolName());
200 if(temp){
201 fdstatus=INSIDE;
202 if(gMC->IsTrackEntering()){ fdstatus=ENTERING; }
203 if(gMC->IsTrackExiting() || gMC->IsTrackOut() || gMC->IsTrackStop() ) { fdstatus=EXITING;}
204 temp->FillHit(fdstatus,gMC);
205 }
206
207 }
208
209 void PamVMCApplication::PostTrack()
210 {
211 // User actions after finishing of each track
212
213 }
214
215 void PamVMCApplication::FinishPrimary()
216 {
217 // User actions after finishing of a primary track
218 fPrimaryGenerator->SetGood(PamVMCVolCrossMgr::Instance()->IsTrackGood());
219 PamVMCVolCrossMgr::Instance()->Reset();
220 fVerbose.FinishPrimary();
221 }
222
223 void PamVMCApplication::FinishEvent()
224 {
225 // User actions after finishing of an event
226
227
228
229 PamVMCSDMgr::Instance()->Digitize();
230 //EventNo needs for trigger, particle PDG needs for TOF..
231 PamVMCDigMgr::Instance()->Digitize(fEventNo,fPrimaryGenerator->GetParticle());
232
233 PamVMCSDMgr::Instance()->Compress();
234
235 PamVMCRawMgr::Instance()->WriteEvent();
236
237 fRootManager.Fill();
238
239 fPrimaryGenerator->ClearPrimCol();
240
241 PamVMCSDMgr::Instance()->Clear();
242
243 fStack->Reset();
244
245 //getchar();
246
247 }
248
249 void PamVMCApplication::FinishRun()
250 {
251 fVerbose.FinishRun();
252
253 fRootManager.WriteAll();
254
255 PamVMCDigMgr::Instance()->FinishRun(); //write RunTrailer to buffer
256
257 PamVMCRawMgr::Instance()->FinishRun(); //write RunTrailer to file and close;
258
259 PamVMCDigMgr::Instance()->PrintCollections();
260
261 }
262
263 void PamVMCApplication::Field(const Double_t* x, Double_t* b) const
264 {
265 PamVMCFieldMgr::Instance()->Field(x,b);
266 }

  ViewVC Help
Powered by ViewVC 1.1.23