QwAnalysis
QwDriftChamber.h
Go to the documentation of this file.
1 /**********************************************************\
2 * File: QwDriftChamber.h *
3 * *
4 * Author: P. M. King, Rakitha Beminiwattha *
5 * J. H. Lee *
6 * Time-stamp: Wednesday, July 28 11:30:21 EDT 2010 *
7 \**********************************************************/
8 
9 
10 #ifndef __QWDRIFTCHAMBER__
11 #define __QWDRIFTCHAMBER__
12 
13 #include "TH1F.h"
14 #include "TH2F.h"
15 #include "TTree.h"
16 
17 
18 #include "QwHit.h"
19 #include "QwHitContainer.h"
20 
21 #include "QwTypes.h"
22 #include "QwDetectorInfo.h"
23 
24 #include <exception>
25 #include <iostream>
26 #include <fstream>
27 #include <string>
28 #include <iomanip>
29 
30 #include "VQwSubsystemTracking.h"
31 #include "QwF1TDContainer.h"
32 
33 ///
34 /// \ingroup QwTracking
36  /******************************************************************
37  * Class: QwDriftChamber
38  *
39  *
40  ******************************************************************/
41  private:
42  /// Private default constructor (not implemented, will throw linker error on use)
44 
45  public:
46  /// Constructor with name
47  QwDriftChamber(const TString& region_tmp);
48  /// Constructor with name and hit list
49  QwDriftChamber(const TString& region_tmp,std::vector< QwHit > &fWireHits_TEMP);
50  /// Virtual destructor
51  virtual ~QwDriftChamber();
52 
53  /// Define options
54  static void DefineOptions(QwOptions& options) {
55  options.AddOptions("Tracking options")("QwDriftChamber.print-f1tdc-configuration",
56  po::value<bool>(&fPrintF1TDCConfiguration)->default_bool_value(false),
57  "print F1TDC configuration");
58  }
59 
60  /* Member functions derived from VQwSubsystem. */
61 
62  virtual Int_t LoadInputParameters(TString mapfile){return 0;};
63 
64 
65  Int_t ProcessEvBuffer(const UInt_t roc_id, const UInt_t bank_id, UInt_t* buffer, UInt_t num_words);
66 
67  Int_t ProcessConfigurationBuffer(const UInt_t roc_id, const UInt_t bank_id, UInt_t* buffer, UInt_t num_words);
68  void PrintConfigurationBuffer(UInt_t *buffer,UInt_t num_words);
69  void ReportConfiguration();
70 
71  //separate meanings in VDC and HDC
72  virtual void SubtractReferenceTimes() = 0;
73  virtual void ProcessEvent() = 0;
74 
75  virtual void ClearEventData() = 0;
76 
77 
78 
79 
80  /* Unique member functions */
81 
82 
84 
85 
87 
88 
89 
90  virtual void GetHitList(QwHitContainer & grandHitContainer)
91  {
92  grandHitContainer.Append(fWireHits);
93  };
94 
95  void GetTDCHitList(QwHitContainer & grandHitContainer)
96  {
97  grandHitContainer.Append(fTDCHits);
98  };
99 
100  protected:
101 
102  static bool fPrintF1TDCConfiguration;
103 
104  virtual void FillRawTDCWord(Int_t bank_index, Int_t slot_num, Int_t chan, UInt_t data) = 0;
105  virtual Int_t AddChannelDefinition() = 0;
106  virtual Int_t BuildWireDataStructure(const UInt_t chan, const EQwDetectorPackage package, const Int_t octant, const Int_t plane, const Int_t wire)=0;
107  virtual Double_t CalculateDriftDistance(Double_t drifttime, QwDetectorID detector)=0;
108 
110  virtual void ConstructHistograms(TDirectory *folder, TString &prefix) = 0;
111  virtual void FillHistograms() = 0;
112 
113  virtual Int_t LoadTimeWireOffset(TString t0_map) = 0;
114  virtual void SubtractWireTimeOffset() = 0;
115  // virtual void ApplyTimeCalibration() = 0;
116 
117 
118  Int_t LinkReferenceChannel(const UInt_t chan, const Int_t plane, const Int_t wire);
120  Int_t RegisterROCNumber(const UInt_t roc_id, const UInt_t bank_id);
121  Int_t RegisterSubbank(const UInt_t bank_id);
122  Int_t RegisterSlotNumber(const UInt_t slot_id); // Tells this object that it will decode data from the current bank
123  Int_t GetTDCIndex(size_t bank_index, size_t slot_num) const;
124  Bool_t IsSlotRegistered(Int_t bank_index, Int_t slot_num) const
125  {
126  return (GetTDCIndex(bank_index,slot_num) != -1);
127  };
128 
129 
130  TString fRegion; /// Name of this subsystem (the region).
134 
135  //static const UInt_t kMaxNumberOfTDCsPerROC;
136  static const UInt_t kMaxNumberOfSlotsPerROC;
137  static const Int_t kReferenceChannelPlaneNumber; // plane is Int_t
138  static const Int_t kCodaMasterPlaneNumber;
139 
141 
142 
144 
145  std::vector< std::vector<Int_t> > fTDC_Index; // TDC index, indexed by bank_index and slot_number
146 
147  std::vector< std::pair<Int_t, Int_t> > fReferenceChannels;
148  // reference chans number <first:tdc_index, second:channel_number>
149  // fReferenceChannels[tdc_index,channel_number][ num of [tdc,chan] set]
150  std::vector< std::vector<Double_t> > fReferenceData;
151  std::vector< std::vector<Double_t> > fReferenceMaster;
152  // wire number < reference time >
153  // we use a wire number of QwHit to save a bank id of a reference time.
154  // thus, for fReferenceData, the wire (fElement) is the same as
155  // its bankid. Therefore we can use fReferenceData[bankid][reference time]
156 
157 
158  std::vector< QwHit > fTDCHits;
159  std::vector< QwHit > &fWireHits;
160  std::vector< Int_t > fWiresPerPlane;
161 
165 
166  // NOTE: The plane and wire indices count from "1" instead
167  // of from "0".
168  // When you're creating loops, just be careful that
169  // you don't try to use the first (zero-th) element
170  // of either index.
171 
172 
173  /*=====
174  * Histograms should be listed below here.
175  * They should be pointers to histograms which will be created
176  * inside the ConstructHistograms()
177  */
178  std::vector <TH1F*> TotHits;
179  std::vector <TH1F*> TOFP;
180  std::vector <TH1F*> TOFP_raw;
181  std::vector <TH1F*> WiresHit;
182  std::vector <TH2F*> TOFW;
183  std::vector <TH2F*> TOFW_raw;
184  std::vector <TH2F*> HitsWire;
185 
187  TotHits.clear();
188  TOFP.clear();
189  TOFP_raw.clear();
190  WiresHit.clear();
191  TOFW.clear();
192  TOFW_raw.clear();
193  HitsWire.clear();
194  }
195 
196 
197  //below are the data structures that are used in HDC/VDC
198 
199 
200  std::vector< std::vector< QwDetectorID > > fTDCPtrs;
201  // Indexed by TDC_index and Channel; gives the package, plane and wire assignment.
202  std::vector< std::vector< QwDetectorInfo > > fWireData;
203  // Indexed by package, plane this contains detector geometry information for each region;
204 
205  std::vector< std::vector< UInt_t > > fDirectionData;
206  //Indexed by pckg and plane each element represent the wire direction ( a value from 0 to 6)- Rakitha(10/23/2008)
207 
208 
209 
210 };
211 
212 #endif
std::vector< TH1F * > TOFP
std::vector< std::vector< Int_t > > fTDC_Index
virtual void SubtractReferenceTimes()=0
std::vector< TH1F * > WiresHit
virtual void GetHitList(QwHitContainer &grandHitContainer)
Get the hit list.
std::vector< std::vector< Double_t > > fReferenceMaster
size_t fCurrentBankIndex
Name of this subsystem (the region).
F1TDCs configuration and reference siganls container.
void GetTDCHitList(QwHitContainer &grandHitContainer)
#define default_bool_value(b)
Definition: QwOptions.h:51
void PrintConfigurationBuffer(UInt_t *buffer, UInt_t num_words)
std::vector< TH2F * > HitsWire
void InitHistogramPointers()
std::vector< TH2F * > TOFW_raw
void ClearAllBankRegistrations()
static const UInt_t kMaxNumberOfSlotsPerROC
An options class.
Definition: QwOptions.h:133
std::vector< std::vector< UInt_t > > fDirectionData
static const Int_t kCodaMasterPlaneNumber
Int_t ProcessEvBuffer(const UInt_t roc_id, const UInt_t bank_id, UInt_t *buffer, UInt_t num_words)
TODO: The non-event-type-aware ProcessEvBuffer routine should be replaced with the event-type-aware v...
std::vector< TH2F * > TOFW
EQwDetectorPackage
Definition: QwTypes.h:70
virtual void ProcessEvent()=0
Int_t ProcessConfigurationBuffer(const UInt_t roc_id, const UInt_t bank_id, UInt_t *buffer, UInt_t num_words)
virtual void ConstructHistograms()
Construct the histograms for this subsystem.
Definition: VQwSubsystem.h:209
std::vector< TH1F * > TOFP_raw
po::options_description_easy_init AddOptions(const std::string &blockname="Specialized options")
Add an option to a named block or create new block.
Definition: QwOptions.h:164
virtual Int_t LoadTimeWireOffset(TString t0_map)=0
static bool fPrintF1TDCConfiguration
static void DefineOptions(QwOptions &options)
Define options.
Draft skeleton for the decoding-to-QTR interface class.
Int_t LinkReferenceChannel(const UInt_t chan, const Int_t plane, const Int_t wire)
std::vector< std::vector< QwDetectorID > > fTDCPtrs
std::vector< std::pair< Int_t, Int_t > > fReferenceChannels
virtual void SubtractWireTimeOffset()=0
QwDriftChamber()
Private default constructor (not implemented, will throw linker error on use)
UInt_t kMaxNumberOfChannelsPerTDC
QwF1TDContainer * fF1TDContainer
std::vector< TH1F * > TotHits
virtual Double_t CalculateDriftDistance(Double_t drifttime, QwDetectorID detector)=0
std::vector< std::vector< Double_t > > fReferenceData
Int_t RegisterSubbank(const UInt_t bank_id)
void ReportConfiguration()
Int_t GetTDCIndex(size_t bank_index, size_t slot_num) const
Int_t RegisterSlotNumber(const UInt_t slot_id)
static const Int_t kReferenceChannelPlaneNumber
virtual ~QwDriftChamber()
Virtual destructor.
virtual Int_t LoadInputParameters(TString mapfile)
Mandatory parameter file definition.
Bool_t IsSlotRegistered(Int_t bank_index, Int_t slot_num) const
virtual Int_t BuildWireDataStructure(const UInt_t chan, const EQwDetectorPackage package, const Int_t octant, const Int_t plane, const Int_t wire)=0
virtual void FillHistograms()=0
Fill the histograms for this subsystem.
virtual Int_t AddChannelDefinition()=0
std::vector< QwHit > & fWireHits
MQwF1TDC fF1TDCDecoder
std::vector< QwHit > fTDCHits
void Append(const QwHitContainer &mylist)
virtual void ClearEventData()=0
void FillDriftDistanceToHits()
virtual void FillRawTDCWord(Int_t bank_index, Int_t slot_num, Int_t chan, UInt_t data)=0
void FillHardwareErrorSummary()
Hardware error summary.
std::vector< Int_t > fWiresPerPlane
Int_t RegisterROCNumber(const UInt_t roc_id, const UInt_t bank_id)
Tell the object that it will decode data from this ROC and sub-bank.
F1TDCReferenceContainer * fF1RefContainer
std::vector< std::vector< QwDetectorInfo > > fWireData