QwAnalysis
QwDriftChamberHDC Class Reference

#include <QwDriftChamberHDC.h>

+ Inheritance diagram for QwDriftChamberHDC:
+ Collaboration diagram for QwDriftChamberHDC:

Public Member Functions

 QwDriftChamberHDC (TString name)
 
virtual ~QwDriftChamberHDC ()
 
void SubtractReferenceTimes ()
 
void ProcessEvent ()
 
Int_t LoadChannelMap (TString mapfile)
 Mandatory map file definition. More...
 
void ClearEventData ()
 
void ProcessOptions (QwOptions &options)
 Process the command line options. More...
 
- Public Member Functions inherited from QwDriftChamber
 QwDriftChamber (const TString &region_tmp)
 Constructor with name. More...
 
 QwDriftChamber (const TString &region_tmp, std::vector< QwHit > &fWireHits_TEMP)
 Constructor with name and hit list. More...
 
virtual ~QwDriftChamber ()
 Virtual destructor. More...
 
virtual Int_t LoadInputParameters (TString mapfile)
 Mandatory parameter file definition. More...
 
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 version. More...
 
Int_t ProcessConfigurationBuffer (const UInt_t roc_id, const UInt_t bank_id, UInt_t *buffer, UInt_t num_words)
 
void PrintConfigurationBuffer (UInt_t *buffer, UInt_t num_words)
 
void ReportConfiguration ()
 
void FillDriftDistanceToHits ()
 
void FillHardwareErrorSummary ()
 Hardware error summary. More...
 
virtual void GetHitList (QwHitContainer &grandHitContainer)
 Get the hit list. More...
 
void GetTDCHitList (QwHitContainer &grandHitContainer)
 
- Public Member Functions inherited from VQwSubsystemTracking
 VQwSubsystemTracking (TString name)
 Constructor with name. More...
 
virtual ~VQwSubsystemTracking ()
 Default destructor. More...
 
Int_t LoadGeometryDefinition (TString mapfile)
 Load geometry definition for tracking subsystems. More...
 
Int_t LoadCrosstalkDefinition (TString mapfile)
 Load crosstalk definition for tracking subsystems. More...
 
const QwGeometryGetDetectorInfo () const
 Get the detector geometry information. More...
 
virtual void ConstructBranchAndVector (TTree *tree, TString &prefix, std::vector< Double_t > &values)
 Construct the branch and tree vector. More...
 
virtual void ConstructBranch (TTree *tree, TString &prefix)
 Construct the branch and tree vector. More...
 
virtual void ConstructBranch (TTree *tree, TString &prefix, QwParameterFile &trim_file)
 Construct the branch and tree vector based on the trim file. More...
 
virtual void FillTreeVector (std::vector< Double_t > &values) const
 Fill the tree vector. More...
 
- Public Member Functions inherited from VQwSubsystem
 VQwSubsystem (const TString &name)
 Constructor with name. More...
 
 VQwSubsystem (const VQwSubsystem &orig)
 Copy constructor by object. More...
 
virtual ~VQwSubsystem ()
 Default destructor. More...
 
TString GetSubsystemName () const
 
Bool_t HasDataLoaded () const
 
void SetParent (QwSubsystemArray *parent)
 Set the parent of this subsystem to the specified array. More...
 
QwSubsystemArrayGetParent (const unsigned int parent=0) const
 Get the parent of this subsystem. More...
 
VQwSubsystemGetSibling (const std::string &name) const
 Get the sibling with specified name. More...
 
Bool_t PublishInternalValue (const TString &name, const TString &desc, const VQwHardwareChannel *value) const
 Publish a variable name to the parent subsystem array. More...
 
virtual Bool_t PublishInternalValues () const
 Publish all variables of the subsystem. More...
 
virtual Bool_t PublishByRequest (TString device_name)
 Try to publish an internal variable matching the submitted name. More...
 
Bool_t RequestExternalValue (const TString &name, VQwHardwareChannel *value) const
 Request a named value which is owned by an external subsystem; the request will be handled by the parent subsystem array. More...
 
virtual const VQwHardwareChannelReturnInternalValue (const TString &name) const
 Return a pointer to a varialbe to the parent subsystem array to be delivered to a different subsystem. More...
 
virtual Bool_t ReturnInternalValue (const TString &name, VQwHardwareChannel *value) const
 Return a named value to the parent subsystem array to be delivered to a different subsystem. More...
 
virtual std::vector< TString > GetParamFileNameList ()
 
virtual std::map< TString,
TString > 
GetDetectorMaps ()
 
virtual Int_t LoadDetectorMaps (QwParameterFile &file)
 Parse parameter file to find the map files. More...
 
virtual Int_t LoadEventCuts (TString mapfile)
 Optional event cut file. More...
 
void SetEventTypeMask (const UInt_t mask)
 Set event type mask. More...
 
UInt_t GetEventTypeMask () const
 Get event type mask. More...
 
virtual Int_t ProcessEvBuffer (const UInt_t event_type, const UInt_t roc_id, const UInt_t bank_id, UInt_t *buffer, UInt_t num_words)
 
virtual void ExchangeProcessedData ()
 Request processed data from other subsystems for internal use in the second event processing stage. Not all derived classes will require data from other subsystems. More...
 
virtual void ProcessEvent_2 ()
 Process the event data again, including data from other subsystems. Not all derived classes will require a second stage of event data processing. More...
 
virtual void AtEndOfEventLoop ()
 Perform actions at the end of the event loop. More...
 
virtual void RandomizeEventData (int helicity=0, double time=0.0)
 
virtual void EncodeEventData (std::vector< UInt_t > &buffer)
 
virtual void PrintInfo () const
 Print some information about the subsystem. More...
 
virtual VQwSubsystemoperator= (VQwSubsystem *value)
 Assignment Note: Must be called at the beginning of all subsystems routine call to operator=(VQwSubsystem *value) by VQwSubsystem::operator=(value) More...
 
virtual void PrintDetectorMaps (Bool_t status) const
 
virtual void ConstructHistograms ()
 Construct the histograms for this subsystem. More...
 
virtual void ConstructHistograms (TDirectory *folder)
 Construct the histograms for this subsystem in a folder. More...
 
virtual void ConstructHistograms (TString &prefix)
 Construct the histograms for this subsystem with a prefix. More...
 
virtual void ConstructBranchAndVector (TTree *tree, std::vector< Double_t > &values)
 Construct the branch and tree vector. More...
 
virtual void ConstructTree ()
 Construct the tree for this subsystem. More...
 
virtual void ConstructTree (TDirectory *folder)
 Construct the tree for this subsystem in a folder. More...
 
virtual void ConstructTree (TString &prefix)
 Construct the tree for this subsystem with a prefix. More...
 
virtual void ConstructTree (TDirectory *folder, TString &prefix)
 Construct the tree for this subsystem in a folder with a prefix. More...
 
virtual void FillTree ()
 Fill the tree for this subsystem. More...
 
virtual void DeleteTree ()
 Delete the tree for this subsystem. More...
 
- Public Member Functions inherited from MQwHistograms
void ShareHistograms (const MQwHistograms *source)
 Share histogram pointers between objects. More...
 
- Public Member Functions inherited from MQwCloneable< VQwSubsystem, QwDriftChamberHDC >
virtual ~MQwCloneable ()
 Virtual destructor. More...
 
virtual VQwSubsystemClone () const
 Concrete clone method. More...
 
const VQwFactory< VQwSubsystem > * Factory () const
 Factory getter. More...
 
- Public Member Functions inherited from VQwCloneable< VQwSubsystem >
virtual ~VQwCloneable ()
 Virtual destructor. More...
 
std::string GetClassName () const
 Get demangled name of this class. More...
 

Static Public Member Functions

static void DefineOptions (QwOptions &options)
 
- Static Public Member Functions inherited from QwDriftChamber
static void DefineOptions (QwOptions &options)
 Define options. More...
 
- Static Public Member Functions inherited from VQwSubsystem
static void DefineOptions ()
 Define options function (note: no virtual static functions in C++) More...
 
- Static Public Member Functions inherited from MQwCloneable< VQwSubsystem, QwDriftChamberHDC >
static VQwSubsystemCreate (const std::string &name)
 Object creation. More...
 
static QwDriftChamberHDCCast (QwDriftChamberHDC *type)
 Object dynamic cast. More...
 

Protected Member Functions

void FillRawTDCWord (Int_t bank_index, Int_t slot_num, Int_t chan, UInt_t data)
 
Int_t AddChannelDefinition ()
 
Int_t BuildWireDataStructure (const UInt_t chan, const EQwDetectorPackage package, const Int_t octant, const Int_t plane, const Int_t wire)
 
Double_t CalculateDriftDistance (Double_t drifttime, QwDetectorID detector)
 
void ConstructHistograms (TDirectory *folder, TString &prefix)
 Construct the histograms for this subsystem in a folder with a prefix. More...
 
void FillHistograms ()
 Fill the histograms for this subsystem. More...
 
Int_t LoadTimeWireOffset (TString t0_map)
 
void LoadTtoDParameters (TString ttod_map)
 
void SubtractWireTimeOffset ()
 
void UpdateHits ()
 
- Protected Member Functions inherited from QwDriftChamber
Int_t LinkReferenceChannel (const UInt_t chan, const Int_t plane, const Int_t wire)
 
void ClearAllBankRegistrations ()
 
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. More...
 
Int_t RegisterSubbank (const UInt_t bank_id)
 
Int_t RegisterSlotNumber (const UInt_t slot_id)
 
Int_t GetTDCIndex (size_t bank_index, size_t slot_num) const
 
Bool_t IsSlotRegistered (Int_t bank_index, Int_t slot_num) const
 
void InitHistogramPointers ()
 
- Protected Member Functions inherited from VQwSubsystem
void UpdatePublishedValue (const TString &name, VQwHardwareChannel *data_channel)
 
void ClearAllBankRegistrations ()
 Clear all registration of ROC and Bank IDs for this subsystem. More...
 
Int_t RegisterSubbank (const UInt_t bank_id)
 Tell the object that it will decode data from this sub-bank in the ROC currently open for registration. More...
 
Int_t GetSubbankIndex () const
 
Int_t GetSubbankIndex (const UInt_t roc_id, const UInt_t bank_id) const
 
void SetDataLoaded (Bool_t flag)
 
Int_t FindIndex (const std::vector< UInt_t > &myvec, const UInt_t value) const
 
Bool_t Compare (VQwSubsystem *source)
 
- Protected Member Functions inherited from MQwHistograms
 MQwHistograms ()
 Default constructor. More...
 
 MQwHistograms (const MQwHistograms &source)
 Copy constructor. More...
 
virtual ~MQwHistograms ()
 Virtual destructor. More...
 
virtual MQwHistogramsoperator= (const MQwHistograms &value)
 
void Fill_Pointer (TH1_ptr hist_ptr, Double_t value)
 
void AddHistogram (TH1 *h)
 Register a histogram. More...
 

Protected Attributes

Double_t trig_h1
 
std::vector< std::vector
< std::vector< Double_t > > > 
fTimeWireOffsets
 
std::vector< Double_t > fTtoDNumbers
 
Int_t fR2Octant
 
- Protected Attributes inherited from QwDriftChamber
TString fRegion
 
size_t fCurrentBankIndex
 Name of this subsystem (the region). More...
 
Int_t fCurrentSlot
 
Int_t fCurrentTDCIndex
 
UInt_t kMaxNumberOfChannelsPerTDC
 
Int_t fNumberOfTDCs
 
std::vector< std::vector< Int_t > > fTDC_Index
 
std::vector< std::pair< Int_t,
Int_t > > 
fReferenceChannels
 
std::vector< std::vector
< Double_t > > 
fReferenceData
 
std::vector< std::vector
< Double_t > > 
fReferenceMaster
 
std::vector< QwHitfTDCHits
 
std::vector< QwHit > & fWireHits
 
std::vector< Int_t > fWiresPerPlane
 
MQwF1TDC fF1TDCDecoder
 
QwF1TDContainerfF1TDContainer
 
F1TDCReferenceContainerfF1RefContainer
 
std::vector< TH1F * > TotHits
 
std::vector< TH1F * > TOFP
 
std::vector< TH1F * > TOFP_raw
 
std::vector< TH1F * > WiresHit
 
std::vector< TH2F * > TOFW
 
std::vector< TH2F * > TOFW_raw
 
std::vector< TH2F * > HitsWire
 
std::vector< std::vector
< QwDetectorID > > 
fTDCPtrs
 
std::vector< std::vector
< QwDetectorInfo > > 
fWireData
 
std::vector< std::vector
< UInt_t > > 
fDirectionData
 
- Protected Attributes inherited from VQwSubsystemTracking
QwGeometry fDetectorInfo
 Geometry information of this subsystem. More...
 
size_t fTreeArrayIndex
 Tree indices. More...
 
size_t fTreeArrayNumEntries
 
Double_t fF1TDCResolutionNS
 
- Protected Attributes inherited from VQwSubsystem
std::map< TString,
VQwHardwareChannel * > 
fPublishedInternalValues
 Map of published internal values. More...
 
std::vector< std::vector
< TString > > 
fPublishList
 List of parameters to be published (loaded at the channel map) More...
 
TString fSystemName
 Name of this subsystem. More...
 
UInt_t fEventTypeMask
 Mask of event types. More...
 
Bool_t fIsDataLoaded
 Has this subsystem gotten data to be processed? More...
 
std::vector< TString > fDetectorMapsNames
 
std::map< TString, TString > fDetectorMaps
 
Int_t fCurrentROC_ID
 ROC ID that is currently being processed. More...
 
Int_t fCurrentBank_ID
 Bank ID that is currently being processed. More...
 
std::vector< UInt_t > fROC_IDs
 Vector of ROC IDs associated with this subsystem. More...
 
std::vector< std::vector
< UInt_t > > 
fBank_IDs
 Vector of Bank IDs per ROC ID associated with this subsystem. More...
 
std::vector< QwSubsystemArray * > fArrays
 Vector of pointers to subsystem arrays that contain this subsystem. More...
 
- Protected Attributes inherited from MQwHistograms
std::vector< TH1_ptrfHistograms
 Histograms associated with this data element. More...
 

Additional Inherited Members

- Static Protected Attributes inherited from QwDriftChamber
static bool fPrintF1TDCConfiguration = true
 
static const UInt_t kMaxNumberOfSlotsPerROC = 21
 
static const Int_t kReferenceChannelPlaneNumber = 99
 
static const Int_t kCodaMasterPlaneNumber = 98
 

Detailed Description

Definition at line 19 of file QwDriftChamberHDC.h.

Constructor & Destructor Documentation

QwDriftChamberHDC::QwDriftChamberHDC ( TString  name)
inline

Definition at line 26 of file QwDriftChamberHDC.h.

26 : VQwSubsystem(name), QwDriftChamber(name,fTDCHits) { };
QwDriftChamber()
Private default constructor (not implemented, will throw linker error on use)
std::vector< QwHit > fTDCHits
virtual QwDriftChamberHDC::~QwDriftChamberHDC ( )
inlinevirtual

Definition at line 27 of file QwDriftChamberHDC.h.

27 { };

Member Function Documentation

Int_t QwDriftChamberHDC::AddChannelDefinition ( )
protectedvirtual

Implements QwDriftChamber.

Definition at line 403 of file QwDriftChamberHDC.cc.

References QwDriftChamber::fTDC_Index, QwDriftChamber::fTDCPtrs, QwDriftChamber::fWireData, QwDriftChamber::fWiresPerPlane, QwDriftChamber::kReferenceChannelPlaneNumber, and OK.

Referenced by LoadChannelMap().

404 {
405  bool temp_local_debug = false;
406 
407  if (temp_local_debug){
408  std::cout << " QwDriftChamberHDC::AddChannelDefinition"<<std::endl;
409  }
410 
411  std::size_t i =0;
412 
413  fWireData.resize(fWiresPerPlane.size());
414 
415  for (i=0; i<fWiresPerPlane.size(); i++) {
416 
417  if(temp_local_debug){
418  std::cout << "wire #" << i
419  << " " << fWiresPerPlane.at(i)
420  << std::endl;
421  }
422  fWireData.at(i).resize(fWiresPerPlane.at(i));
423  }
424 
425  Int_t mytdc = 0;
426  for ( i=0; i<fTDC_Index.size(); i++){
427  for (size_t j=0; j<fTDC_Index.at(i).size(); j++){
428  mytdc = fTDC_Index.at(i).at(j);
429  if (mytdc not_eq -1){
430  for (size_t k=0; k<fTDCPtrs.at(mytdc).size(); k++){
431  // Int_t package = fTDCPtrs.at(mytdc).at(k).fPackage;
432  Int_t plane = fTDCPtrs.at(mytdc).at(k).fPlane;
433  if (( plane>0) and (plane not_eq (Int_t) kReferenceChannelPlaneNumber) ){
434  Int_t wire = fTDCPtrs.at(mytdc).at(k).fElement;
435  fWireData.at(plane).at(wire).SetElectronics(i,j,k);
436  if(temp_local_debug) {
437  std::cout << "mytdc " << mytdc
438  << "wire #" << wire
439  << " plane " << plane
440  << " (i j k) (" << i <<" " << j << " "<< k << ")"
441  << std::endl;
442  }
443  }
444  }
445  }
446  }
447  }
448 
449  return OK;
450 }
std::vector< std::vector< Int_t > > fTDC_Index
std::vector< std::vector< QwDetectorID > > fTDCPtrs
static const Int_t kReferenceChannelPlaneNumber
#define OK
std::vector< Int_t > fWiresPerPlane
std::vector< std::vector< QwDetectorInfo > > fWireData

+ Here is the caller graph for this function:

Int_t QwDriftChamberHDC::BuildWireDataStructure ( const UInt_t  chan,
const EQwDetectorPackage  package,
const Int_t  octant,
const Int_t  plane,
const Int_t  wire 
)
protectedvirtual

Implements QwDriftChamber.

Definition at line 367 of file QwDriftChamberHDC.cc.

References QwDriftChamber::fCurrentTDCIndex, QwDriftChamber::fTDCPtrs, QwDriftChamber::fWiresPerPlane, QwDriftChamber::kReferenceChannelPlaneNumber, QwDriftChamber::LinkReferenceChannel(), and OK.

Referenced by LoadChannelMap().

372 {
373  if (plane == kReferenceChannelPlaneNumber ){
374  LinkReferenceChannel(chan, plane, wire);
375  }
376  else {
377  fTDCPtrs.at(fCurrentTDCIndex).at(chan).fOctant = octant;
378  fTDCPtrs.at(fCurrentTDCIndex).at(chan).fPackage = package;
379  fTDCPtrs.at(fCurrentTDCIndex).at(chan).fPlane = plane;
380  fTDCPtrs.at(fCurrentTDCIndex).at(chan).fElement = wire;
381 
382  // std::cout << "fWiresPerPlane.size() " << fWiresPerPlane.size()
383  // << " plane " << plane
384  // << std::endl;
385 
386  if (plane >= (Int_t) fWiresPerPlane.size()){ // plane is Int_t
387  fWiresPerPlane.resize(plane+1);
388  // size() is one more larger than last plane number
389  // For HDC, plane 1,2,....,12
390  // vector 0,1,2,....,12
391  // thus, vector.size() returns 13
392  // So the magic number "1" is.
393  // Wednesday, July 28 21:56:34 EDT 2010, jhlee
394  }
395  if (wire>=fWiresPerPlane.at(plane)){
396  fWiresPerPlane.at(plane) = wire+1;
397  }
398  }
399  return OK;
400 }
Int_t LinkReferenceChannel(const UInt_t chan, const Int_t plane, const Int_t wire)
std::vector< std::vector< QwDetectorID > > fTDCPtrs
static const Int_t kReferenceChannelPlaneNumber
#define OK
std::vector< Int_t > fWiresPerPlane

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

Double_t QwDriftChamberHDC::CalculateDriftDistance ( Double_t  drifttime,
QwDetectorID  detector 
)
protectedvirtual

Implements QwDriftChamber.

Definition at line 209 of file QwDriftChamberHDC.cc.

References fTtoDNumbers.

210 {
211  //0.00545449393 0.0668865488 0.000352462179 -2.00383196E-05 3.57577417E-07 -2.82802562E-09 7.89009965E-12
212  // Double_t dt_ = 0.12 * (drifttime-trig_h1 + 933.0);
213  //Double_t dt_ = 0.12 * (drifttime);
214 
215  Double_t dt_= drifttime;
216  //Double_t dd_ = 0.0,t0=0.0;
217  Double_t dd_ = 0.0;
218  Double_t resolution=1.0;
219  //dt_-=t0;
220  if(dt_>=0 && dt_<130){
221  Int_t index = (Int_t) (dt_/resolution);
222  dd_=( dt_-resolution*index ) /resolution * ( fTtoDNumbers.at ( index+1 )-fTtoDNumbers.at ( index ) ) +fTtoDNumbers.at ( index );}
223  else{
224  dd_=-50;
225  }
226  //dd_ = -0.0138896
227  // + 0.00987685 * dt_
228  // + 0.00100368 * dt_ * dt_
229  // + (-1.79785E-06 * dt_ * dt_ * dt_)
230  // + ( -8.96859E-08 * dt_ * dt_ * dt_ * dt_)
231  // + (6.11736E-10 * dt_ * dt_ * dt_ * dt_ * dt_)
232  // + ( -1.15889E-12 * dt_ * dt_ * dt_ * dt_ * dt_ * dt_);
233  /*
234  dd_ = 0.078067
235  + 0.0437269 * dt_
236  + 0.00117295 * dt_ * dt_
237  + (-2.25244E-05 * dt_ * dt_ * dt_ )
238  + (1.56369E-07 * dt_ * dt_ * dt_ * dt_ )
239  + (-4.93323E-10 * dt_ * dt_ * dt_ * dt_ * dt_ )
240  + (5.91555E-13 * dt_ * dt_ * dt_ * dt_ * dt_ * dt_ );
241  */
242  /* modified on 418
243  dd_ = 0.134955
244  + 0.0173593 * dt_
245  + 0.00254912 * dt_ * dt_
246  + (-4.26414E-05 * dt_ * dt_ * dt_ )
247  + (2.86094E-07 * dt_ * dt_ * dt_ * dt_ )
248  + (-8.85512E-10 * dt_ * dt_ * dt_ * dt_ * dt_ )
249  + (1.04669E-12 * dt_ * dt_ * dt_ * dt_ * dt_ * dt_ );
250  */
251  /*dd_ = -0.0382463
252  + 0.0704644 * dt_
253  - 0.00089724454912 * dt_ * dt_
254  + (4.8703E-05 * dt_ * dt_ * dt_ )
255  + (-8.40592E-07 * dt_ * dt_ * dt_ * dt_ )
256  + (5.61285E-09 * dt_ * dt_ * dt_ * dt_ * dt_ )
257  + (-1.31489E-11 * dt_ * dt_ * dt_ * dt_ * dt_ * dt_ );
258  */
259  //dd_ = 0.015993
260  // + 0.0466621 * dt_
261  // + 0.00180132 * dt_ * dt_
262  // + (-5.96573E-05 * dt_ * dt_ * dt_ )
263  // + (1.18497E-06 * dt_ * dt_ * dt_ * dt_ )
264  // + (-1.2071E-08 * dt_ * dt_ * dt_ * dt_ * dt_ )
265  // + (4.50584E-11 * dt_ * dt_ * dt_ * dt_ * dt_ * dt_ );
266  /* dd_ = -0.0181483
267  + 0.0867558 * dt_
268  -0.00115513 * dt_ * dt_
269  + (4.41846E-05 * dt_ * dt_ * dt_)
270  + (-6.9308E-07 * dt_ * dt_ * dt_ * dt_)
271  + (4.30973E-09 * dt_ * dt_ * dt_ * dt_ * dt_)
272  + (-9.32219E-12 * dt_ * dt_ * dt_ * dt_ * dt_ * dt_);
273  */
274  //dd_ = -0.119434
275  // + 0.0321642 * dt_
276  // + 0.0011334 * dt_ * dt_
277  // + (-1.08359E-05 * dt_ * dt_ * dt_)
278  // + (-2.78272E-08 * dt_ * dt_ * dt_ * dt_)
279  // + (5.97029E-10 * dt_ * dt_ * dt_ * dt_ * dt_)
280  // + (-1.67945E-12 * dt_ * dt_ * dt_ * dt_ * dt_ * dt_);
281 
282  //std::cout<<" Drift distance "<<dd_<<" Drift time "<<dt_<<" Original value "<<drifttime<<std::endl;
283 
284  // This formula returns a drift distance in [mm], but the tracking code
285  // expects a number in [cm]. Since the geometry definitions are in [cm],
286  // we should stick to [cm] as unit of distance. (wdconinc)
287  return dd_ / 10.0;
288 }
std::vector< Double_t > fTtoDNumbers
void QwDriftChamberHDC::ClearEventData ( )
virtual

Implements QwDriftChamber.

Definition at line 804 of file QwDriftChamberHDC.cc.

References F1TDCReferenceContainer::ClearEventData(), QwDetectorID::fElement, QwDriftChamber::fF1RefContainer, QwDetectorID::fPlane, QwDriftChamber::fReferenceData, QwDriftChamber::fTDCHits, QwDriftChamber::fWireData, and VQwSubsystem::SetDataLoaded().

805 {
806  SetDataLoaded(kFALSE);
807  QwDetectorID this_det;
808  // Loop through fTDCHits, to decide which detector elements need to be cleared.
809  std::vector<QwHit>::iterator end=fTDCHits.end();
810  for (std::vector<QwHit>::iterator hit1=fTDCHits.begin(); hit1!=end; hit1++) {
811  this_det = hit1->GetDetectorID();
812  fWireData.at(this_det.fPlane).at(this_det.fElement).ClearHits();
813  }
814  fTDCHits.clear();
815 
817  std::size_t i = 0;
818  for (i=0; i<fReferenceData.size(); i++) {
819  fReferenceData.at(i).clear();
820  }
821  return;
822 }
Int_t fElement
trace number for R1; wire number for R2 &amp; R3; PMT number for others
Definition: QwTypes.h:253
Int_t fPlane
R or theta index for R1; plane index for R2 &amp; R3.
Definition: QwTypes.h:251
std::vector< std::vector< Double_t > > fReferenceData
std::vector< QwHit > fTDCHits
F1TDCReferenceContainer * fF1RefContainer
std::vector< std::vector< QwDetectorInfo > > fWireData
void SetDataLoaded(Bool_t flag)
Definition: VQwSubsystem.h:305

+ Here is the call graph for this function:

void QwDriftChamberHDC::ConstructHistograms ( TDirectory *  folder,
TString &  prefix 
)
protectedvirtual

Construct the histograms for this subsystem in a folder with a prefix.

Implements QwDriftChamber.

Definition at line 622 of file QwDriftChamberHDC.cc.

References QwDriftChamber::fWiresPerPlane, VQwSubsystem::GetSubsystemName(), QwDriftChamber::HitsWire, QwDriftChamber::TOFP, QwDriftChamber::TOFP_raw, QwDriftChamber::TOFW, QwDriftChamber::TOFW_raw, QwDriftChamber::TotHits, and QwDriftChamber::WiresHit.

623 {
624  // If we have defined a subdirectory in the ROOT file, then change into it.
625  if ( folder ) folder->cd();
626  // Now create the histograms...
627  TString region = GetSubsystemName();
628  // Loop over the number of planes.
629 
630  const Short_t buffer_size = 2000;
631  Float_t bin_offset = -0.5;
632 
633  std::size_t total_plane_number = 0;
634  total_plane_number = fWiresPerPlane.size();
635 
636  TotHits.resize(total_plane_number);
637  TOFP.resize(total_plane_number);
638  TOFP_raw.resize(total_plane_number);
639  WiresHit.resize(total_plane_number);
640  TOFW.resize(total_plane_number);
641  TOFW_raw.resize(total_plane_number);
642  HitsWire.resize(total_plane_number);
643 
644  std::size_t iplane = 0;
645 
646  //std::cout << "QwDriftChamberHDC::ConstructHistograms, "
647  // << "we are contructing histograms with index from 0 to " <<total_plane_number
648  // << "\n"
649  // << "Thus, fWiresPerPlane.size() returns "
650  // << total_plane_number
651  // << " and its array definition is ["
652  // << total_plane_number
653  // << "]."
654  // << " And hist[i] <-> hist.at(i) <-> fWiresPerplane[i] <-> fWiresPerPlane.at(i)"
655  // << std::endl;
656 
657  // wire_per_plane is the number of wire per plane?
658  //
659  // we skip the first zero-th plane or wire histogram. thus
660  // i starts with '1'. hist[0] is NULL
661 
662  for ( iplane=1; iplane<total_plane_number; iplane++) {
663 
664  // push_back can "push" iplane = 1 into TotHits.at(0) ??
665  TotHits[iplane] = new TH1F(
666  Form("%s%sHitsOnEachWirePlane%d", prefix.Data(), region.Data(), (Int_t) iplane),
667  Form("Total hits on all wires in plane %d", (Int_t) iplane),
668  fWiresPerPlane[iplane], bin_offset, fWiresPerPlane[iplane]+bin_offset
669  );
670 
671  TotHits[iplane]->GetXaxis()->SetTitle("Wire #");
672  TotHits[iplane]->GetYaxis()->SetTitle("Events");
673 
674  WiresHit[iplane] = new TH1F(
675  Form("%s%sWiresHitPlane%d", prefix.Data(), region.Data(), (Int_t) iplane),
676  Form("Number of Wires Hit in plane %d",(Int_t) iplane),
677  20, bin_offset, 20+bin_offset
678  );
679  WiresHit[iplane]->GetXaxis()->SetTitle("Wires Hit per Event");
680  WiresHit[iplane]->GetYaxis()->SetTitle("Events");
681 
682  HitsWire[iplane] = new TH2F(
683  Form("%s%sHitsOnEachWirePerEventPlane%d", prefix.Data(), region.Data(), (Int_t) iplane),
684  Form("hits on all wires per event in plane %d", (Int_t) iplane),
685  fWiresPerPlane[iplane],bin_offset,fWiresPerPlane[iplane]+bin_offset,
686  7, -bin_offset, 7-bin_offset
687  );
688  HitsWire[iplane]->GetXaxis()->SetTitle("Wire Number");
689  HitsWire[iplane]->GetYaxis()->SetTitle("Hits");
690 
691  TOFP[iplane] = new TH1F(
692  Form("%s%sTimeofFlightPlane%d", prefix.Data(), region.Data(), (Int_t) iplane),
693  Form("Subtracted time of flight for events in plane %d", (Int_t) iplane),
694  400,0,0
695  );
696  TOFP[iplane] -> SetDefaultBufferSize(buffer_size);
697  TOFP[iplane] -> GetXaxis()->SetTitle("Time of Flight");
698  TOFP[iplane] -> GetYaxis()->SetTitle("Hits");
699 
700 
701  TOFP_raw[iplane] = new TH1F(
702  Form("%s%sRawTimeofFlightPlane%d", prefix.Data(), region.Data(), (Int_t) iplane),
703  Form("Raw time of flight for events in plane %d", (Int_t) iplane),
704  // 400,-65000,65000);
705  400, 0,0
706  );
707  TOFP_raw[iplane] -> SetDefaultBufferSize(buffer_size);
708  TOFP_raw[iplane]->GetXaxis()->SetTitle("Time of Flight");
709  TOFP_raw[iplane]->GetYaxis()->SetTitle("Hits");
710 
711  TOFW[iplane] = new TH2F(
712  Form("%s%sTimeofFlightperWirePlane%d", prefix.Data(), region.Data(), (Int_t) iplane),
713  Form("Subtracted time of flight for each wire in plane %d",(Int_t) iplane),
714  fWiresPerPlane[iplane], bin_offset, fWiresPerPlane[iplane]+bin_offset,
715  100,-40000,65000
716  );
717  // why this range is not -65000 ??
718  TOFW[iplane]->GetXaxis()->SetTitle("Wire Number");
719  TOFW[iplane]->GetYaxis()->SetTitle("Time of Flight");
720 
721  TOFW_raw[iplane] = new TH2F(
722  Form("%s%sRawTimeofFlightperWirePlane%d", prefix.Data() ,region.Data(),(Int_t)iplane),
723  Form("Raw time of flight for each wire in plane %d",(Int_t)iplane),
724  fWiresPerPlane[iplane], bin_offset, fWiresPerPlane[iplane]+bin_offset,
725  100,-40000,65000
726  );
727  // why this range is not -65000 ??
728  TOFW_raw[iplane]->GetXaxis()->SetTitle("Wire Number");
729  TOFW_raw[iplane]->GetYaxis()->SetTitle("Time of Flight");
730  }
731  return;
732 }
std::vector< TH1F * > TOFP
std::vector< TH1F * > WiresHit
std::vector< TH2F * > HitsWire
std::vector< TH2F * > TOFW_raw
std::vector< TH2F * > TOFW
std::vector< TH1F * > TOFP_raw
std::vector< TH1F * > TotHits
std::vector< Int_t > fWiresPerPlane
TString GetSubsystemName() const
Definition: VQwSubsystem.h:93

+ Here is the call graph for this function:

void QwDriftChamberHDC::DefineOptions ( QwOptions options)
static

Definition at line 610 of file QwDriftChamberHDC.cc.

References QwOptions::AddOptions().

611 {
612  options.AddOptions("Tracking options") ("R2-octant",
613  po::value<int>()->default_value(1),
614  "MD Package 2 of R2 is in front of" );
615 }
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

+ Here is the call graph for this function:

void QwDriftChamberHDC::FillHistograms ( )
protectedvirtual

Fill the histograms for this subsystem.

Implements QwDriftChamber.

Definition at line 736 of file QwDriftChamberHDC.cc.

References QwLog::endl(), QwDetectorID::fElement, QwDetectorID::fPlane, QwDriftChamber::fTDCHits, QwDriftChamber::fWireData, QwDriftChamber::fWiresPerPlane, QwDetectorInfo::GetNumHits(), VQwSubsystem::HasDataLoaded(), QwDriftChamber::HitsWire, QwMessage, QwDriftChamber::TOFP, QwDriftChamber::TOFP_raw, QwDriftChamber::TOFW, QwDriftChamber::TOFW_raw, QwDriftChamber::TotHits, and QwDriftChamber::WiresHit.

737 {
738  Bool_t local_debug = false;
739  if (not HasDataLoaded()) return;
740 
741  QwDetectorID this_detid;
742  QwDetectorInfo *this_det;
743 
744  // Fill all of the histograms.
745 
746  std::vector<Int_t> wireshitperplane(fWiresPerPlane.size(),0);
747 
748  UInt_t raw_time = 0;
749  Double_t time = 0.0;
750 
751  Int_t plane = 0;
752  Int_t element = 0;
753 
754  std::vector<QwHit>::iterator end=fTDCHits.end();
755  for(std::vector<QwHit>::iterator hit=fTDCHits.begin(); hit!=end; hit++) {
756 
757  this_detid = hit->GetDetectorID();
758  plane = this_detid.fPlane;
759  element = this_detid.fElement;
760 
761  if (plane<=0 or element<=0) {
762  if (local_debug) {
763  QwMessage << "QwDriftChamberHDC::FillHistograms: Bad plane or element index:"
764  << " fPlane = " << plane
765  << " fElement= " << element
766  << QwLog::endl;
767  }
768  continue;
769  }
770 
771 
772  this_det= &(fWireData.at(plane).at(element));
773 
774  if (hit->IsFirstDetectorHit()) {
775  // If this is the first hit for this detector, then let's plot the
776  // total number of hits this wire had.
777  HitsWire[plane]->Fill(element,this_det->GetNumHits());
778 
779  // Also increment the total number of events in whichthis wire was hit.
780  TotHits[plane]->Fill(element,1);
781  // Increment the number of wires hit in this plane
782  wireshitperplane.at(plane) += 1;
783  }
784 
785  raw_time = hit->GetRawTime();
786  time = hit->GetTimeNs();
787 
788  // Fill ToF histograms
789  TOFP_raw[plane]->Fill(raw_time);
790  TOFW_raw[plane]->Fill(element, raw_time);
791  TOFP [plane]->Fill(time);
792  TOFW [plane]->Fill(element, time);
793  }
794 
795  std::size_t iplane = 0;
796 
797  for (iplane=1; iplane<fWiresPerPlane.size(); iplane++) {
798  WiresHit[iplane]->Fill(wireshitperplane[iplane]);
799  }
800  return;
801 }
#define QwMessage
Predefined log drain for regular messages.
Definition: QwLog.h:50
std::vector< TH1F * > TOFP
std::vector< TH1F * > WiresHit
std::vector< TH2F * > HitsWire
std::vector< TH2F * > TOFW_raw
Bool_t HasDataLoaded() const
Definition: VQwSubsystem.h:94
std::vector< TH2F * > TOFW
std::vector< TH1F * > TOFP_raw
Int_t fElement
trace number for R1; wire number for R2 &amp; R3; PMT number for others
Definition: QwTypes.h:253
Int_t fPlane
R or theta index for R1; plane index for R2 &amp; R3.
Definition: QwTypes.h:251
std::vector< TH1F * > TotHits
static std::ostream & endl(std::ostream &)
End of the line.
Definition: QwLog.cc:299
std::vector< QwHit > fTDCHits
std::vector< Int_t > fWiresPerPlane
std::vector< std::vector< QwDetectorInfo > > fWireData

+ Here is the call graph for this function:

void QwDriftChamberHDC::FillRawTDCWord ( Int_t  bank_index,
Int_t  slot_num,
Int_t  chan,
UInt_t  data 
)
protectedvirtual

Implements QwDriftChamber.

Definition at line 291 of file QwDriftChamberHDC.cc.

References QwDriftChamber::fDirectionData, QwDriftChamber::fF1RefContainer, QwDriftChamber::fReferenceData, QwDriftChamber::fTDCHits, QwDriftChamber::fTDCPtrs, QwDriftChamber::GetTDCIndex(), kDirectionNull, QwDriftChamber::kReferenceChannelPlaneNumber, kRegionID2, F1TDCReferenceContainer::SetReferenceSignal(), and QwHit::WireMatches().

292 {
293  Bool_t local_debug = false;
294  Int_t tdcindex = 0;
295 
296  tdcindex = GetTDCIndex(bank_index,slot_num);
297 
298  if (tdcindex not_eq -1) {
299 
300  Int_t hitcnt = 0;
301  EQwDirectionID direction = kDirectionNull;
302 
303  fF1RefContainer->SetReferenceSignal(bank_index, slot_num, chan, data, local_debug);
304 
305 
306  Int_t plane = fTDCPtrs.at(tdcindex).at(chan).fPlane;
307  Int_t wire = fTDCPtrs.at(tdcindex).at(chan).fElement;
308  EQwDetectorPackage package = fTDCPtrs.at(tdcindex).at(chan).fPackage;
309  Int_t octant = fTDCPtrs.at(tdcindex).at(chan).fOctant;
310 
311 
312 
313  if (plane == -1 or wire == -1){
314  // This channel is not connected to anything. Do nothing.
315  }
316  else if (plane == kReferenceChannelPlaneNumber){
317  fReferenceData.at(wire).push_back(data);
318  //now wire contains the value fCurrentBankIndex so we can assign the ref timing data to it.
319  }
320  else {
321 
322  direction = (EQwDirectionID)fDirectionData.at(package-1).at(plane-1);
323  //Wire Direction is accessed from the vector -Rakitha (10/23/2008)
324  //hitCount gives the total number of hits on a given wire -Rakitha (10/23/2008)
325  hitcnt = std::count_if(fTDCHits.begin(), fTDCHits.end(),
326  boost::bind(
327  &QwHit::WireMatches, _1, kRegionID2, boost::ref(package), boost::ref(plane), boost::ref(wire)
328  )
329  );
330 
331  fTDCHits.push_back(
332  QwHit(
333  bank_index,
334  slot_num,
335  chan,
336  hitcnt,
337  kRegionID2,
338  package,
339  octant,
340  plane,
341  direction,
342  wire,
343  data
344  )
345  );
346  //in order-> bank index, slot num, chan, hitcount, region=2, package, plane,,direction, wire,wire hit time
347  if(local_debug) {
348  std::cout << "At QwDriftChamberHDC::FillRawTDCWord " << "\n";
349  std::cout << " hitcnt " << hitcnt
350  << " plane " << plane
351  << " wire " << wire
352  << " package " << package
353  << " fTDCHits.size() " << fTDCHits.size()
354  << std::endl;
355  }
356 
357  }
358 
359  }
360 
361  return;
362 }
std::vector< std::vector< UInt_t > > fDirectionData
EQwDetectorPackage
Definition: QwTypes.h:70
std::vector< std::vector< QwDetectorID > > fTDCPtrs
void SetReferenceSignal(Int_t bank_index, Int_t slot, Int_t chan, UInt_t data, Bool_t debug=false)
Bool_t WireMatches(EQwRegionID region, EQwDetectorPackage package, Int_t plane, Int_t wire)
Definition: QwHit.cc:357
std::vector< std::vector< Double_t > > fReferenceData
Int_t GetTDCIndex(size_t bank_index, size_t slot_num) const
EQwDirectionID
Definition: QwTypes.h:41
static const Int_t kReferenceChannelPlaneNumber
Hit structure uniquely defining each hit.
Definition: QwHit.h:43
std::vector< QwHit > fTDCHits
F1TDCReferenceContainer * fF1RefContainer

+ Here is the call graph for this function:

Int_t QwDriftChamberHDC::LoadChannelMap ( TString  mapfile)
virtual

Mandatory map file definition.

Implements VQwSubsystem.

Definition at line 466 of file QwDriftChamberHDC.cc.

References AddChannelDefinition(), BuildWireDataStructure(), QwDriftChamber::fCurrentBankIndex, QwDriftChamber::fCurrentSlot, VQwSubsystem::fDetectorMaps, QwDriftChamber::fDirectionData, QwDriftChamber::fF1RefContainer, fR2Octant, QwParameterFile::GetUInt(), QwOptions::GetValue(), gQwOptions, kDirectionNull, kPackage1, kPackage2, QwDriftChamber::kReferenceChannelPlaneNumber, LoadTimeWireOffset(), LoadTtoDParameters(), OK, QwDriftChamber::RegisterROCNumber(), QwDriftChamber::RegisterSlotNumber(), and QwDriftChamber::RegisterSubbank().

467 {
468  Bool_t local_debug = false;
469  LoadTtoDParameters ( "R2_TtoDTable.map" );
470  LoadTimeWireOffset ( "R2_timeoffset.map" );
471 
472  TString varname, varvalue;
473  UInt_t value = 0;
474  UInt_t chan = 0;
475  UInt_t DIRMODE = 0;
476  Int_t plane = 0;
477  Int_t wire = 0;
478  TString name = "";
479 
480  EQwDetectorPackage package = kPackageNull;
481  EQwDirectionID direction = kDirectionNull;
482 
483  fDirectionData.resize(2);//currently we have 2 package - Rakitha (10/23/2008)
484  fDirectionData.at(0).resize(12); //currently we have 12 wire planes in each package - Rakitha (10/23/2008)
485  fDirectionData.at(1).resize(12); //currently we have 12 wire planes in each package - Rakitha (10/23/2008)
486 
487  QwParameterFile mapstr(mapfile.Data()); //Open the file
488  fDetectorMaps.insert(mapstr.GetParamFileNameContents());
489 
490  while (mapstr.ReadNextLine()) {
491  mapstr.TrimComment('!'); // Remove everything after a '!' character.
492  mapstr.TrimWhitespace(); // Get rid of leading and trailing spaces.
493  if (mapstr.LineIsEmpty()) continue;
494 
495  if (mapstr.HasVariablePair("=",varname,varvalue)) {
496  // This is a declaration line. Decode it.
497  varname.ToLower();
498  value = QwParameterFile::GetUInt(varvalue);
499  if (value ==0){
500  value = atol(varvalue.Data());
501  }
502  if (varname=="roc") {
503  RegisterROCNumber(value,0);
504  DIRMODE=0;
505  }
506  else if (varname=="bank") {
507  RegisterSubbank(value);
508  DIRMODE=0;
509  }
510  else if (varname=="pkg") {
511  //this will identify the coming sequence is wire plane to direction mapping - Rakitha
512  DIRMODE=1;
513  package=(EQwDetectorPackage)value;
514  }
515  else if (varname=="slot") {
516  RegisterSlotNumber(value);
517  DIRMODE=0;
518  }
519 
520  }
521  else if (DIRMODE==0) {
522  // Break this line into tokens to process it.
523  chan = mapstr.GetTypedNextToken<Int_t>();
524  plane = mapstr.GetTypedNextToken<Int_t>();
525  // wire = mapstr.GetTypedNextToken<Int_t>();
526  name = mapstr.GetTypedNextToken<TString>();
527 
528  Int_t octant = 0;
529  fR2Octant = gQwOptions.GetValue<int>("R2-octant");
530  if (package == kPackage1) octant = (fR2Octant + 4) % 8;
531  if (package == kPackage2) octant = fR2Octant;
532 
533  if(local_debug) {
534  printf("chan %8d plan %4d wire %12s\n", chan, plane, name.Data());
535  }
536  if (plane==kReferenceChannelPlaneNumber) {
537  fF1RefContainer -> AddF1TDCReferenceSignal(new F1TDCReferenceSignal(fCurrentBankIndex, fCurrentSlot, chan, name));
538 
539  if (name=="MasterTrigger" ) {
540  wire = 0;
541  BuildWireDataStructure(chan, package, octant, plane, wire);
542  }
543  }
544  else {
545  wire = name.Atoi();
546  // VDC and HDC
547  BuildWireDataStructure(chan, package, octant, plane, wire);
548  }
549 
550  }
551  else if (DIRMODE==1) {
552  //this will decode the wire plane directions - Rakitha
553  plane = mapstr.GetTypedNextToken<Int_t>();
554  direction = (EQwDirectionID) mapstr.GetTypedNextToken<Int_t>();
555  fDirectionData.at(package-1).at(plane-1) = direction;
556  }
557 
558  }
559 
560 
561 
562 
563  // Construct the wire data structures.
565 
566  /*
567  for (size_t i=0; i<fDirectionData.at(0).size(); i++){
568  std::cout<<"Direction data Plane "<<i+1<<" "<<fDirectionData.at(0).at(i)<<std::endl;
569  }
570  */
571  //
572 
573  mapstr.Close(); // Close the file (ifstream)
574  // / ReportConfiguration();
575 
576  return OK;
577 }
std::map< TString, TString > fDetectorMaps
Definition: VQwSubsystem.h:322
size_t fCurrentBankIndex
Name of this subsystem (the region).
static UInt_t GetUInt(const TString &varvalue)
std::vector< std::vector< UInt_t > > fDirectionData
QwOptions gQwOptions
Definition: QwOptions.cc:27
EQwDetectorPackage
Definition: QwTypes.h:70
T GetValue(const std::string &key)
Get a templated value.
Definition: QwOptions.h:240
Int_t BuildWireDataStructure(const UInt_t chan, const EQwDetectorPackage package, const Int_t octant, const Int_t plane, const Int_t wire)
Int_t RegisterSubbank(const UInt_t bank_id)
EQwDirectionID
Definition: QwTypes.h:41
Int_t RegisterSlotNumber(const UInt_t slot_id)
static const Int_t kReferenceChannelPlaneNumber
void LoadTtoDParameters(TString ttod_map)
Int_t LoadTimeWireOffset(TString t0_map)
#define OK
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

+ Here is the call graph for this function:

Int_t QwDriftChamberHDC::LoadTimeWireOffset ( TString  t0_map)
protectedvirtual

Implements QwDriftChamber.

Definition at line 902 of file QwDriftChamberHDC.cc.

References VQwSubsystem::fDetectorMaps, fTimeWireOffsets, and OK.

Referenced by LoadChannelMap().

903 {
904  //std::cout << "beginning to load t0 file... " << std::endl;
905  //
906  QwParameterFile mapstr ( t0_map.Data() );
907  fDetectorMaps.insert(mapstr.GetParamFileNameContents());
908 
909  TString varname,varvalue;
910  Int_t plane=0,wire=0;
911  Double_t t0 = 0.0;
912  EQwDetectorPackage package = kPackageNull;
913 
914  while ( mapstr.ReadNextLine() )
915  {
916  mapstr.TrimComment ( '!' );
917  mapstr.TrimWhitespace();
918  if ( mapstr.LineIsEmpty() ) continue;
919  if ( mapstr.HasVariablePair ( "=",varname,varvalue ) )
920  {
921  varname.ToLower();
922  if ( varname=="package" )
923  {
924  package = ( EQwDetectorPackage ) atoi ( varvalue.Data() );
925  if ( package> ( Int_t ) fTimeWireOffsets.size() ) fTimeWireOffsets.resize ( package );
926  }
927  else if ( varname=="plane" )
928  {
929  //std::cout << "package: " << fTimeWireOffsets.at(package-1).size()<< std::endl;
930  plane = atoi ( varvalue.Data() );
931  if ( plane> ( Int_t ) fTimeWireOffsets.at ( package-1 ).size() ) fTimeWireOffsets.at ( package-1 ).resize ( plane );
932  //std::cout << "plane: " << fTimeWireOffsets.at(package-1).size()<< std::endl;
933 
934  // To Siyuan, * : can package be obtained before plane in while loop? if plane is before package
935  // we have at(-1), thus, if condition is always "false", I guess.
936  // * : if, else if then can we expect something more?
937  // from Han
938  }
939  continue;
940  }
941 
942  wire = mapstr.GetTypedNextToken<Int_t>();
943  t0 = mapstr.GetTypedNextToken<Int_t>();
944 
945  if ( wire > ( Int_t ) fTimeWireOffsets.at ( package-1 ).at ( plane-1 ).size() )
946  {
947  fTimeWireOffsets.at ( package-1 ).at ( plane-1 ).resize ( wire );
948 
949  fTimeWireOffsets.at ( package-1 ).at ( plane-1 ).at ( wire-1 ) = t0;
950  }
951 
952  }
953  //
954 
955  mapstr.Close(); // Close the file (ifstream)
956  return OK;
957 }
std::map< TString, TString > fDetectorMaps
Definition: VQwSubsystem.h:322
std::vector< std::vector< std::vector< Double_t > > > fTimeWireOffsets
EQwDetectorPackage
Definition: QwTypes.h:70
#define OK

+ Here is the caller graph for this function:

void QwDriftChamberHDC::LoadTtoDParameters ( TString  ttod_map)
protected

Definition at line 882 of file QwDriftChamberHDC.cc.

References VQwSubsystem::fDetectorMaps, and fTtoDNumbers.

Referenced by LoadChannelMap().

883 {
884  QwParameterFile mapstr(ttod_map.Data());
885  fDetectorMaps.insert(mapstr.GetParamFileNameContents());
886 
887  while (mapstr.ReadNextLine())
888  {
889  mapstr.TrimComment('!');
890  mapstr.TrimWhitespace();
891  if (mapstr.LineIsEmpty()) continue;
892 
893  /* Double_t t = */ mapstr.GetTypedNextToken<Double_t>();
894  Double_t d = mapstr.GetTypedNextToken<Double_t>();
895  fTtoDNumbers.push_back(d);
896  }
897 
898  mapstr.Close(); // Close the file (ifstream)
899 }
std::map< TString, TString > fDetectorMaps
Definition: VQwSubsystem.h:322
std::vector< Double_t > fTtoDNumbers

+ Here is the caller graph for this function:

void QwDriftChamberHDC::ProcessEvent ( )
virtual

Implements QwDriftChamber.

Definition at line 452 of file QwDriftChamberHDC.cc.

References QwDriftChamber::FillDriftDistanceToHits(), VQwSubsystem::HasDataLoaded(), SubtractReferenceTimes(), SubtractWireTimeOffset(), and UpdateHits().

453 {
454  if (not HasDataLoaded()) return;
455 
456  SubtractReferenceTimes(); // A.U. unit
457  UpdateHits(); // Fill QwDetectorInfo, fTimeRes (ns), and fTimeNs (ns) in QwHits
458 
459  SubtractWireTimeOffset(); // Subtract t0 offset (ns) and educated guesss t0 (a.u.) from fTime (a.u.) and fTimeNs (ns)
460  FillDriftDistanceToHits(); // must call GetTimeNs() instead of GetTime()
461 
462  return;
463 }
Bool_t HasDataLoaded() const
Definition: VQwSubsystem.h:94
void FillDriftDistanceToHits()

+ Here is the call graph for this function:

void QwDriftChamberHDC::ProcessOptions ( QwOptions options)
virtual

Process the command line options.

Reimplemented from VQwSubsystem.

Definition at line 617 of file QwDriftChamberHDC.cc.

References fR2Octant, and QwOptions::GetValue().

618 {
619  fR2Octant = options.GetValue<int>("R2-octant");
620 }
T GetValue(const std::string &key)
Get a templated value.
Definition: QwOptions.h:240

+ Here is the call graph for this function:

void QwDriftChamberHDC::SubtractReferenceTimes ( )
virtual

Implements QwDriftChamber.

Definition at line 17 of file QwDriftChamberHDC.cc.

References QwDriftChamber::fF1RefContainer, QwDriftChamber::fF1TDContainer, QwDriftChamber::fTDCHits, F1TDCReferenceContainer::GetReferenceTimeAU(), VQwSubsystem::GetSubsystemName(), QwMessage, and QwF1TDContainer::ReferenceSignalCorrection().

Referenced by ProcessEvent().

18 {
19  UInt_t bank_index = 0;
20  Double_t raw_time_arb_unit = 0.0;
21  Double_t ref_time_arb_unit = 0.0;
22  Double_t time_arb_unit = 0.0;
23  // EQwDetectorPackage package = kPackageNull;
24 
25  Bool_t local_debug = false;
26  Int_t slot_num = 0;
27 
28  TString reference_name1 = "MasterTrigger";
29  TString reference_name2 = "CopyMasterTrigger";
30 
31  std::vector<QwHit>::iterator end=fTDCHits.end();
32  for ( std::vector<QwHit>::iterator hit=fTDCHits.begin(); hit!=end; hit++ )
33  {
34 
35  bank_index = hit -> GetSubbankID();
36  slot_num = hit -> GetModule();
37  raw_time_arb_unit = (Double_t) hit -> GetRawTime();
38 
39 
40  // John L. introduced this package match reference time subtraction for only Region 3,
41  // I was told it would be better to extend this concept to Region2, so I did, and I
42  // got the worst result..... Q2 = +2.735199e-02 its RMS +1.192504e-02
43  // if I don't use this method in Region2, Q2 = +2.292274e-02 and RMS = +7.133603e-03
44  // if I don't use this method in Region3, I couldn't get "peak" structure
45  // in time distributions of four VDC chambers.
46  //
47  // Why does this method give us a good result for only Region 3? MUX?
48  // I don't know ..
49  // Friday, February 17 14:01:12 EST 2012, jhlee
50 
51  // package = hit -> GetPackage();
52  // if (package == kPackageUp) { // Up is 1
53  // reference_name1 = "TrigScintPkg1";
54  // }
55  // else if (package == kPackageDown) {// Down is 2
56  // reference_name1 = "TrigScintPkg2";
57  // }
58  // else {
59  // reference_name1 = "MasterTrigger";
60  // }
61 
62  ref_time_arb_unit = fF1RefContainer -> GetReferenceTimeAU(bank_index, reference_name1);
63  //
64  // if there is no reference time due to a channel error, try to use a copy of mater trigger
65  //
66  if(ref_time_arb_unit==0.0) {
67  ref_time_arb_unit = fF1RefContainer->GetReferenceTimeAU(bank_index, reference_name2);
68  }
69  // second time, it returns 0.0, we simply ignore this event ....
70  // set time zero. ReferenceSignalCorrection() will return zero, and increase RFM counter...
71  //
72  time_arb_unit = fF1TDContainer->ReferenceSignalCorrection(raw_time_arb_unit, ref_time_arb_unit, bank_index, slot_num);
73 
74  hit -> SetTime(time_arb_unit);
75  hit -> SetRawRefTime((UInt_t) ref_time_arb_unit);
76 
77  if(local_debug) {
78  QwMessage << this->GetSubsystemName()
79  << " BankIndex " << std::setw(2) << bank_index
80  << " Slot " << std::setw(2) << slot_num
81  << " RawTime : " << std::setw(6) << raw_time_arb_unit
82  << " RefTime : " << std::setw(6) << ref_time_arb_unit
83  << " time : " << std::setw(6) << time_arb_unit
84  << std::endl;
85 
86  }
87  }
88 
89 
90 
91  return;
92 }
#define QwMessage
Predefined log drain for regular messages.
Definition: QwLog.h:50
Double_t GetReferenceTimeAU(Int_t bank_index, TString name)
Double_t ReferenceSignalCorrection(Double_t raw_time, Double_t ref_time, Int_t bank_index, Int_t slot)
QwF1TDContainer * fF1TDContainer
std::vector< QwHit > fTDCHits
TString GetSubsystemName() const
Definition: VQwSubsystem.h:93
F1TDCReferenceContainer * fF1RefContainer

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void QwDriftChamberHDC::SubtractWireTimeOffset ( )
protectedvirtual

Implements QwDriftChamber.

Definition at line 854 of file QwDriftChamberHDC.cc.

References VQwSubsystemTracking::fF1TDCResolutionNS, QwDriftChamber::fTDCHits, and fTimeWireOffsets.

Referenced by ProcessEvent().

855 {
856  Int_t plane = 0;
857  Int_t wire = 0;
858  EQwDetectorPackage package = kPackageNull;
859  Double_t t0_NS = 0.0;
860 
861  Double_t educated_guess_t0_correction_AU = 11255.0;
862  Double_t educated_guess_t0_correction_NS = educated_guess_t0_correction_AU*fF1TDCResolutionNS;
863 
864  //hit -> SetTime(time+educated_guess_t0_correction); // an educated guess
865 
866  std::size_t nhits=fTDCHits.size();
867 
868  for(std::size_t i=0;i<nhits;i++)
869  {
870  package = fTDCHits.at(i).GetPackage();
871  plane = fTDCHits.at(i).GetPlane();
872  wire = fTDCHits.at(i).GetElement();
873  t0_NS = fTimeWireOffsets.at ( package-1 ).at ( plane-1 ).at ( wire-1 );
874  t0_NS = t0_NS - educated_guess_t0_correction_NS;
875 
876  fTDCHits.at(i).SubtractTimeNsOffset(t0_NS); // time unit is ns, Replace fTimeNs
877  fTDCHits.at(i).SubtractTimeAuOffset(t0_NS/fF1TDCResolutionNS); // time unit is a.u. Replace fTime
878  }
879  return ;
880 }
std::vector< std::vector< std::vector< Double_t > > > fTimeWireOffsets
EQwDetectorPackage
Definition: QwTypes.h:70
std::vector< QwHit > fTDCHits

+ Here is the caller graph for this function:

void QwDriftChamberHDC::UpdateHits ( )
protected

Definition at line 826 of file QwDriftChamberHDC.cc.

References VQwSubsystemTracking::fDetectorInfo, VQwSubsystemTracking::fF1TDCResolutionNS, QwDriftChamber::fTDCHits, and QwGeometry::in().

Referenced by ProcessEvent().

827 {
828 
829 
830  EQwDetectorPackage package = kPackageNull;
831  Int_t plane = 0;
832 
833  QwDetectorID local_id;
834  QwDetectorInfo * local_info;
835 
836  std::vector<QwHit>::iterator end=fTDCHits.end();
837  for(std::vector<QwHit>::iterator iter=fTDCHits.begin(); iter!=end; ++iter)
838  {
839 
840  local_id = iter->GetDetectorID();
841  package = local_id.fPackage;
842  plane = local_id.fPlane - 1;
843  // ahha, here is a hidden magic number 1.
844  local_info = fDetectorInfo.in(package).at(plane);
845  iter->SetDetectorInfo(local_info);
846  iter->ApplyTimeCalibration(fF1TDCResolutionNS); // Fill fTimeRes and fTimeNs in QwHit
847  }
848 
849 
850  return;
851 }
const QwGeometry in(const EQwRegionID &r) const
Get detectors in given region.
Definition: QwGeometry.h:92
EQwDetectorPackage
Definition: QwTypes.h:70
QwGeometry fDetectorInfo
Geometry information of this subsystem.
std::vector< QwHit > fTDCHits

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

Field Documentation

Int_t QwDriftChamberHDC::fR2Octant
protected

Definition at line 65 of file QwDriftChamberHDC.h.

Referenced by LoadChannelMap(), and ProcessOptions().

std::vector< std::vector< std::vector<Double_t> > > QwDriftChamberHDC::fTimeWireOffsets
protected

Definition at line 62 of file QwDriftChamberHDC.h.

Referenced by LoadTimeWireOffset(), and SubtractWireTimeOffset().

std::vector< Double_t> QwDriftChamberHDC::fTtoDNumbers
protected

Definition at line 63 of file QwDriftChamberHDC.h.

Referenced by CalculateDriftDistance(), and LoadTtoDParameters().

Double_t QwDriftChamberHDC::trig_h1
protected

Definition at line 61 of file QwDriftChamberHDC.h.


The documentation for this class was generated from the following files: