StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
StVpdSimMaker.h
1 /***************************************************************************
2  *
3  * $Id: StVpdSimMaker.h,v 1.2 2021/01/27 04:33:30 geurts Exp $
4  *
5  * Author: Nickolas Luttrell (Rice University)
6  ***************************************************************************
7  *
8  * Description: StVpdSimMaker.h - Virtual base class for Vertex Position
9  * detector simulations
10  *
11  ***************************************************************************
12  *
13  * $Log: StVpdSimMaker.h,v $
14  * Revision 1.2 2021/01/27 04:33:30 geurts
15  * move StVpdSimConfig.h to StBTofUtil area for shared used between StBTofCalibMaker and StVpdSimMaker
16  *
17  * Revision 1.1 2017/03/02 18:38:17 jeromel
18  * First version of Vpd simulations
19  *
20  *
21  ***************************************************************************/
22 
23 #ifndef StVpdSimMaker_HH
24 #define StVpdSimMaker_HH
25 
26 #include "StMaker.h"
27 #include "St_DataSet.h"
28 class TH1F;
29 class TH2F;
30 class TH3F;
31 class TNtuple;
32 class StEvent;
33 class StBTofCollection;
34 class StVpdSimConfig;
35 
36 #include "tables/St_g2t_vpd_hit_Table.h"
37 #include "StMcEvent/StMcEvent.hh"
38 #include "StMcEvent/StMcBTofHitCollection.hh"
39 #include "StBTofUtil/StVpdSimConfig.h"
40 
41 #include <vector>
42 #ifndef ST_NO_NAMESPACES
43 
44 #endif
45 
46 class StVpdSimMaker : public StMaker {
47 
48 public:
49  StVpdSimMaker(const char *name = "VpdSim");
50  virtual ~StVpdSimMaker();
51 
52  void Reset();
53  virtual int Init();
54  int InitRun(int);
55  int FinishRun(int);
56  virtual int Make();
57  virtual int Finish();
58 
59  // Define get and set functions
60 
65 
66  string pullHistFileName();
67  string getParamsFileName() { return mParamsFileName; }
68  void setParamsFile(string fileName = "db/vpdSimParams/vpdSimParams.dat") { mParamsFileName = fileName; }
69  void setBookHisto(bool bookHist) { mBookHisto = bookHist; }
70 
71  virtual const char *GetCVS() const{
72  static const char cvs[] = "Tag $Name: $ $Id: StVpdSimMaker.h,v 1.2 2021/01/27 04:33:30 geurts Exp $ built " __DATE__ " " __TIME__; return cvs;
73  }
74 
75 protected:
76 
77 
79  St_DataSet * mGeantData = nullptr;
80  StEvent * mEvent = nullptr;
81  StMcEvent * mMcEvent = nullptr;
83  std::map<int, StVpdSimConfig::SingleTubeParams> mSimParams;
84 
86  struct VpdSingleHit {
87  int tray;
88  int tubeId;
89  double tof;
90  double t0;
91  double de;
92  double pathL;
93  double q;
94  };
95 
96  // Various general use variables
97 
99  bool mBookHisto;
102 
103  int mNHits = 0;
104  double mVx = 0;
105  double mVy = 0;
106  double mVz = 0;
107  double mSumTubeTime = 0;
108  double mTubeTAvg = 0;
109 
110  double mTStart = 0;
111  double mTubeTAvgWest = 0;
112  double mTubeTAvgEast = 0;
113  float mVpdVertex = 0;
114 
115  // Histogram variables
116  string mHistoFileName = "";
121  TH1F* mNHitsWest;
122  TH1F* mNHitsEast;
123  TH1F* mLeTimeWest;
124  TH1F* mLeTimeEast;
125  TH1F* mTStartHist;
128  TH1F* mZVertexHist;
130  TH1F* mVpdResHist;
132 
133  // Functions
134 
136  int vpdResponse(VpdSingleHit &Hit, g2t_vpd_hit_st* vpd_hit, int vId);
138  double thresholdCut(std::vector<VpdSingleHit> Hits, std::vector<int> Tube_Counts, TH1F* TubeHits, TH1F* NHits);
140  int storeMcVpdHit(StMcBTofHit* mcVpdHit);
142  int fillEvent();
144  int bookHistograms();
145 
146  ClassDef(StVpdSimMaker, 2)
147 };
148 
149 #endif /* StVpdSimMaker_h */
150 
151 
152 
153 // end of StVpdSimMaker.h
StMcBTofHitCollection * mMcBTofHitCollection
The Mc hit collection.
Definition: StVpdSimMaker.h:78
string mParamsFileName
path/name of the calibration file to be passed
string mHistoFileName
histogram file name
int tubeId
The tube id of a given Vpd, with values [1,19].
Definition: StVpdSimMaker.h:88
double mTStart
Start time for an event.
TH1F * mZVertexHist
Histogram of the provided Mc Vertices for all events.
StBTofCollection * mVpdCollection
The StBTofCollection of StBTofHit&#39;s (in this case, vpd hits)
Definition: StVpdSimMaker.h:82
std::map< int, StVpdSimConfig::SingleTubeParams > mSimParams
Map of the calibration parameters to be applied.
Definition: StVpdSimMaker.h:83
bool mUseFileParameters
Default is kFALSE.
double t0
Time offset (currently unused)
Definition: StVpdSimMaker.h:90
TH1F * mVpdResHist
Histogram of the difference between Mc and calculated vertices.
StVpdSimConfig * mSimConfig
The calibration parameters for Vpd.
Definition: StVpdSimMaker.h:98
VpdSingleHit contains the parameters that describe a vpd hit.
Definition: StVpdSimMaker.h:86
int vpdResponse(VpdSingleHit &Hit, g2t_vpd_hit_st *vpd_hit, int vId)
Extracts relevant parameters from a Vpd hit.
TH3F * mResVsNumHits
Histogram fo the difference between Mc and calculated vertices vs. number of hits.
int tray
The vpd tray (121==West, 122==East)
Definition: StVpdSimMaker.h:87
virtual int Make()
int fillEvent()
Fill StEvent from the McBTofCollection.
StMcEvent * mMcEvent
The McEvent info.
Definition: StVpdSimMaker.h:81
StMcBTofHitCollection * GetMcBTofHitCollection() const
Returns the StMcBTofHitCollection of Mc Vpd hits.
Definition: StVpdSimMaker.h:64
float mVpdVertex
The calculated vertex as seen by the Vpd.
TH1F * mNRawHitsWest
Number of hits on each west Vpd tube before threshold cuts.
TH1F * mLeTimeEast
Leading edge times (currently equal to Time of Flight) for east Vpd.
TH1F * mTubeHitsWest
Number of hits on each west Vpd tube after threshold cuts.
TH1F * mTStartHist
mTStart times (in ns)
bool mBookHisto
Default is kFALSE.
Definition: StVpdSimMaker.h:99
double q
Charge (currently set to 0 for Vpd hits)
Definition: StVpdSimMaker.h:93
double mTubeTAvgWest
Corrected event time for the west Vpd.
StBTofCollection * GetVpdCollection() const
Returns the StBTofCollection of Vpd hits.
Definition: StVpdSimMaker.h:62
int storeMcVpdHit(StMcBTofHit *mcVpdHit)
Builds the McBTofCollection, insures no duplicate hits.
TH2F * mLeTubeHitsEast
Leading edge times and Number of tubes hit across events for East Vpd.
double mSumTubeTime
Tracks the time measured by each Vpd tube and sums them.
TH2F * mLeTubeHitsWest
Leading edge times and Number of tubes hit across events for West Vpd.
StVpdSimMaker(const char *name="VpdSim")
double mTubeTAvgEast
Corrected event time for the east Vpd.
int bookHistograms()
Creat the QA histograms.
double thresholdCut(std::vector< VpdSingleHit > Hits, std::vector< int > Tube_Counts, TH1F *TubeHits, TH1F *NHits)
Determines average time information from East and West Vpd, cuts zero-velocity particles.
double mTubeTAvg
Average time lapse seen by the east or west Vpd.
TH1F * mNHitsEast
Number of tubes hit across events for east Vpd.
virtual ~StVpdSimMaker()
StEvent * mEvent
The StEvent info.
Definition: StVpdSimMaker.h:80
TH1F * mNRawHitsEast
Number of hits on each east Vpd tube before threshold cuts.
TH1F * mTubeHitsEast
Number of hits on each east Vpd tube after threshold cuts.
virtual int Init()
TH1F * mLeTimeWest
Leading edge times (currently equal to Time of Flight) for west Vpd.
St_DataSet * mGeantData
Geant data passed into the StVpdSimMaker.
Definition: StVpdSimMaker.h:79
TH1F * mNHitsWest
Number of tubes hit across events for west Vpd.
double de
Energy deposition in GeV.
Definition: StVpdSimMaker.h:91
Event data structure to hold all information from a Monte Carlo simulation. This class is the interfa...
Definition: StMcEvent.hh:169
virtual int Finish()
TH1F * mVpdVertexHist
Histogram of the calculated Vpd vertices for all events.
double tof
Time of flight given in ns.
Definition: StVpdSimMaker.h:89
double pathL
Path length in cm (currently set to -9999)
Definition: StVpdSimMaker.h:92