GammaWare  Head Version for release 0.9
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
GeantLMOF.h
Go to the documentation of this file.
1 #ifndef GW_GEANTLMOF_H
2 #define GW_GEANTLMOF_H
3 
4 #ifndef ROOT_Rtypes
5 #include "Rtypes.h"
6 #endif
7 
8 #ifndef GW_GAMMATRACKER_H
9 #include <GammaTracker.h>
10 #endif
11 
12 #include <iostream>
13 #include <fstream>
14 #include <string>
15 #include <sstream>
16 #include <vector>
17 
18 #define MAXCHARACTER 500
19 
20 class TGeoVolume;
21 
22 namespace Gw {
23 
59 class GeantLMOF
60 {
61  friend GeantLMOF &fastout(GeantLMOF &reader) ;
62  friend GeantLMOF &checkout(GeantLMOF &reader) ;
63 
64 private:
65  bool fFastOut;
66 
68  GeantLMOF *fGamma;
69 
70 protected:
72  UInt_t OnlyType(const GeantLMOF &input, Int_t fTypemin = 1, Int_t fTypemax = 1);
73 
74 public:
76  enum Status { kEmpty, kGood, kEoE, kFail };
77 
78 private:
79  Double_t fRecoilBeta;
80  Double_t fRecoilSigma;
81  Double_t fRecoilDx;
82  Double_t fRecoilDy;
83  Double_t fRecoilDz;
84  Double_t fRecoilAngle;
85  Double_t fTime;
86 
87  UInt_t fNbEmitted;
88  UInt_t fNbImpacts;
89 
90  std::vector<Int_t> fType;
91  std::vector<Int_t> fNdet;
92  std::vector<Int_t> fEvNum;
93  std::vector<Int_t> fNInt;
94  std::vector<Int_t> fPrimary;
95  std::vector<Double_t> fEnergy;
96  std::vector<Double_t> fDx;
97  std::vector<Double_t> fDy;
98  std::vector<Double_t> fDz;
99  std::vector<Double_t> fEdep;
100  std::vector<Double_t> fPx;
101  std::vector<Double_t> fPy;
102  std::vector<Double_t> fPz;
103  std::vector<Double_t> fNseg;
104 
105  Double_t fDistFactor;
106  Double_t fEnerFactor;
107 
108 private:
109  std::ifstream fInfile; // input stream
110 
111  Int_t fCurrent, fMaxFiles, fStatus; // current file number, max file number and current status
112 
113  std::vector<Double_t> fEventsRead; // total number of events read for each files
114  std::vector<std::string> fFileName; // name of each file
115 
116  std::string fPath; // path to reach input files
117 
118  Char_t buf[MAXCHARACTER]; // buffer to get data from files
119 
120 private:
121  void ReadHeader();
122 
123 public:
124  GeantLMOF();
125  ~GeantLMOF();
126 
128 
132  static TGeoVolume *ImportAgataG4(const char *asolid = "asolid", const char *aclust = "aclust", const char *aeuler = "aeuler", const char *top = "World");
133 
134 
135  void SetRecoil(Float_t b=0.0, Float_t dx = 0.0, Float_t dy = 0.0, Float_t dz = 1.0, Float_t angle = 0.0) {fRecoilBeta=b; fRecoilDx=dx; fRecoilDy=dy; fRecoilDz=dz; fRecoilAngle=angle;}
136 
137  Double_t GetDistFactor() const { return fDistFactor; }
138  Double_t GetEnerFactor() const { return fEnerFactor; }
139 
140  Double_t GetfRecoilBeta() const { return fRecoilBeta; }
141  void GetRecoil(Double_t &b, Double_t &dx, Double_t &dy, Double_t &dz) const { b = fRecoilBeta; dx = fRecoilDx; dy = fRecoilDy; dz = fRecoilDz;}
142 
143  UInt_t NbEmitted() const { return fNbEmitted; }
144  UInt_t NbImpacts() const { return fNbImpacts; }
145 
146  const std::vector<Double_t> &GetEnergy() { return fEnergy;}
147  const std::vector<Double_t> &GetDx() { return fDx;}
148  const std::vector<Double_t> &GetDy() { return fDy;}
149  const std::vector<Double_t> &GetDz() { return fDz;}
150 
151  const std::vector<Double_t> &GetEdep() { return fEdep;}
152  const std::vector<Double_t> &GetPx() { return fPx;}
153  const std::vector<Double_t> &GetPy() { return fPy;}
154  const std::vector<Double_t> &GetPz() { return fPz;}
155  const std::vector<Double_t> &GetNseg() { return fNseg;}
156 
157  const std::vector<Int_t> &GetType() { return fType;}
158  const std::vector<Int_t> &GetNdet() { return fNdet;}
159  const std::vector<Int_t> &GetEvNum() { return fEvNum;}
160  const std::vector<Int_t> &GetNInt() { return fNInt;}
161  const std::vector<Int_t> &GetPrimary() { return fPrimary;}
162 
163  const std::vector<Double_t> &GetEventsRead() { return fEventsRead; }
164 
166  void operator>> (GammaTracker *) ;
167 
168  GeantLMOF &operator<< (GeantLMOF &(*pf)(GeantLMOF &)) { return (*pf)(*this);}
169 
171 
178  void AddFile(const char *fname) { std::string s = fname; fFileName.push_back(s); Double_t zero = 0.0; fEventsRead.push_back(zero); fMaxFiles++; }
179  void SetPath(const char *path ="./") { fPath = path; }
180 
182 
186  void NextFile();
187 
189 
197  Bool_t NextEvent(UInt_t mult = kMaxUInt);
198 
199  void ShowCurrentEvent();
201 
205  void ShowCurrentConditions();
206 };
207 
208 GeantLMOF &fastout(GeantLMOF &reader) { reader.fFastOut = true ; return reader; }
209 GeantLMOF &checkout(GeantLMOF &reader) { reader.fFastOut = false ; return reader; }
210 
211 }
212 #endif
213 
TBrowser * b
#define MAXCHARACTER
Definition: GeantLMOF.h:18
const std::vector< Int_t > & GetEvNum()
Definition: GeantLMOF.h:159
Base class to build tracker families.
Definition: GammaTracker.h:29
void AddFile(const char *fname)
to add an input file
Definition: GeantLMOF.h:178
void NextFile()
switch to next file.
Definition: GeantLMOF.cpp:129
GeantLMOF & fastout(GeantLMOF &reader)
Definition: GeantLMOF.h:208
const std::vector< Double_t > & GetDy()
Definition: GeantLMOF.h:148
Double_t GetfRecoilBeta() const
Definition: GeantLMOF.h:140
const std::vector< Int_t > & GetNdet()
Definition: GeantLMOF.h:158
Bool_t NextEvent(UInt_t mult=kMaxUInt)
To load the next event.
Definition: GeantLMOF.cpp:205
void GetRecoil(Double_t &b, Double_t &dx, Double_t &dy, Double_t &dz) const
Definition: GeantLMOF.h:141
friend GeantLMOF & fastout(GeantLMOF &reader)
Definition: GeantLMOF.h:208
Double_t GetEnerFactor() const
Definition: GeantLMOF.h:138
void ShowCurrentEvent()
Definition: GeantLMOF.cpp:71
friend GeantLMOF & checkout(GeantLMOF &reader)
Definition: GeantLMOF.h:209
Status
Reading status.
Definition: GeantLMOF.h:76
void ShowCurrentConditions()
The current status is shown on the standard ouptut.
Definition: GeantLMOF.cpp:58
const std::vector< Int_t > & GetNInt()
Definition: GeantLMOF.h:160
GeantLMOF & checkout(GeantLMOF &reader)
Definition: GeantLMOF.h:209
static TGeoVolume * ImportAgataG4(const char *asolid="asolid", const char *aclust="aclust", const char *aeuler="aeuler", const char *top="World")
to import the agata geometry from the Agata-Geant4 world
Definition: GeantLMOF.cpp:53
Double_t GetDistFactor() const
Definition: GeantLMOF.h:137
const std::vector< Double_t > & GetPz()
Definition: GeantLMOF.h:154
const std::vector< Double_t > & GetPx()
Definition: GeantLMOF.h:152
const std::vector< Double_t > & GetEnergy()
Definition: GeantLMOF.h:146
const std::vector< Int_t > & GetPrimary()
Definition: GeantLMOF.h:161
Class to read/decode Geant List Mode Output Files.
Definition: GeantLMOF.h:59
const std::vector< Double_t > & GetEdep()
Definition: GeantLMOF.h:151
UInt_t NbEmitted() const
Definition: GeantLMOF.h:143
const std::vector< Double_t > & GetNseg()
Definition: GeantLMOF.h:155
void SetPath(const char *path="./")
Definition: GeantLMOF.h:179
const std::vector< Double_t > & GetPy()
Definition: GeantLMOF.h:153
const std::vector< Double_t > & GetDx()
Definition: GeantLMOF.h:147
void SetRecoil(Float_t b=0.0, Float_t dx=0.0, Float_t dy=0.0, Float_t dz=1.0, Float_t angle=0.0)
Definition: GeantLMOF.h:135
const std::vector< Double_t > & GetDz()
Definition: GeantLMOF.h:149
void operator>>(GammaTracker *)
fill the tracker with the current event
Definition: GeantLMOF.cpp:183
UInt_t NbImpacts() const
Definition: GeantLMOF.h:144
const std::vector< Double_t > & GetEventsRead()
Definition: GeantLMOF.h:163
GeantLMOF & operator<<(GeantLMOF &(*pf)(GeantLMOF &))
Definition: GeantLMOF.h:168
UInt_t OnlyType(const GeantLMOF &input, Int_t fTypemin=1, Int_t fTypemax=1)
Fill this using another GeantLMOF by keeping only emmitted particles with fType between fTypemin & fT...
Definition: GeantLMOF.cpp:149
const std::vector< Int_t > & GetType()
Definition: GeantLMOF.h:157