Hall A ROOT/C++ Analyzer (podd)
Loading...
Searching...
No Matches
THaTriggerTime.cxx
Go to the documentation of this file.
1
2// //
3// THaTriggerTime //
4// //
5// Simple processing: report the global offset of the common-start/stop //
6// for a given spectrometer, for use as a common offset to fix-up //
7// places where the absolute time (ie VDC) is used. //
8// //
9// Lookup in the database what the offset is for each trigger type, and //
10// report it for each appropriate event. //
11// //
13
14#include "THaTriggerTime.h"
15#include "THaDetMap.h"
16#include "THaEvData.h"
17#include "TMath.h"
18#include "VarDef.h"
19#include "VarType.h"
20#include "THaAnalyzer.h"
21#include <cstdio>
22#include <vector>
23#include <cassert>
24
25using namespace std;
27
28//____________________________________________________________________________
29THaTriggerTime::THaTriggerTime( const char* name, const char* desc )
30 : TimeCorrectionModule(name, desc, THaAnalyzer::kDecode),
31 fDetMap(new THaDetMap), fTDCRes(0.0), fCommonStop(0), fEvtType(-1)
32{
33 // basic do-nothing-else constructor
34}
35
36//____________________________________________________________________________
38{
39 // Normal destructor
41
42 delete fDetMap;
43}
44
45//____________________________________________________________________________
47{
48 // Reset event-by-event variables
49 TimeCorrectionModule::Clear(opt);
50 fEvtType = -1;
51 fTrgTimes.assign( fTrgTypes.size(), kBig );
52}
53
54//____________________________________________________________________________
56{
57 // Look through the channels to determine the trigger type.
58 // Be certain to decode only once per event.
59 if (fEvtType>0) return kOK;
60
61 // Each trigger is connected to a single TDC channel.
62 // All TDCs are started or stopped with the same signal.
63 // Whichever trigger occurs absolute first based on these TDCs is
64 // considered the event's "trigger type".
65 // Each TDC channel may have multiple hits, in which case the
66 // hit occurring first in time will be used and compared to the
67 // TDC values for the other triggers.
68 // Each trigger type is assigned a fixed time offset fToffsets,
69 // read from the database.
70 // The "event time" is reported as fToffsets[earliest] + fGlOffset.
71 // The actual TDC values of all the triggers are reported in fTrgTimes,
72 // primarily for debugging purposes.
73 UInt_t earliest = kMaxUInt;
74 Double_t reftime = kBig;
75 for( UInt_t imod = 0; imod < fDetMap->GetSize(); imod++ ) {
77
78 // Loop over our short-list of relevant channels.
79 // Continuing with the assumption that each channel has its own entry in
80 // fDetMap, so that fDetMap, fTrgTimes, fTrgTypes, and fToffsets are
81 // parallel.
82
83 // look through each channel and take the earliest hit
84 // In common stop mode, the earliest hit corresponds to the largest TDC data.
85 // This is taken care of by making fTDCRes and hence all fTrgTimes negative.
86 for( UInt_t ihit = 0; ihit < evdata.GetNumHits(d->crate, d->slot, d->lo); ihit++ ) {
87 Double_t t = fTDCRes * evdata.GetData(d->crate, d->slot, d->lo, ihit);
88 if( t < fTrgTimes[imod] )
89 fTrgTimes[imod] = t;
90 }
91 if( earliest == kMaxUInt || fTrgTimes[imod] < reftime ) {
92 earliest = imod;
93 reftime = fTrgTimes[imod];
94 }
95 }
96
98 if( earliest != kMaxUInt ) {
99 assert( earliest < fToffsets.size() );
100 assert( earliest < fTrgTypes.size() );
101 fEvtTime += fToffsets[earliest];
102 fEvtType = fTrgTypes[earliest];
103 }
104 return kOK;
105}
106
107//____________________________________________________________________________
109{
110 // Read this detector's parameters from the database.
111 // 'date' contains the date/time of the run being analyzed.
112
113 const char* const here = "ReadDatabase";
114
116 if( ret != kOK )
117 return ret;
118
119 FILE* file = OpenFile( date );
120 if( !file )
121 return kFileError;
122
123 fCommonStop = 0;
124 fTrgTypes.clear();
125 fToffsets.clear();
126 fDetMap->Clear();
127
128 // Read configuration parameters
129 vector<Double_t> trigdef;
130 DBRequest config_request[] = {
131 { "tdc_res", &fTDCRes },
132 { "trigdef", &trigdef, kDoubleV },
133 { "glob_off", &fGlOffset, kDouble, 0, true },
134 { "common_stop", &fCommonStop, kInt, 0, true },
135 { nullptr }
136 };
137 Int_t err = LoadDB( file, date, config_request, fPrefix );
138 fclose(file);
139 if( err )
140 return err;
141
142 // Sanity checks
143 if( trigdef.size() % 5 != 0 ) {
144 Error( Here(here), "Incorrect number of elements in \"trigdef\" "
145 "parameter = %u. Must be a multiple of 5. Fix database.",
146 static_cast<unsigned int>(trigdef.size()) );
147 return kInitError;
148 }
149
150 // Read in the time offsets, in the format below, to be subtracted from
151 // the times measured in other detectors.
152 //
153 // Trigger types NOT listed are assumed to have a zero offset.
154 //
155 // <TrgType> <time offset in seconds>
156 // eg:
157 // 1 0 crate slot chan
158 // 2 10.e-9 crate slot chan
159
160 UInt_t ndef = trigdef.size() / 5;
161 for( UInt_t i = 0; i<ndef; ++i ) {
162 UInt_t base = 5*i;
163 Int_t itrg = TMath::Nint( trigdef[base] );
164 Double_t toff = trigdef[base+1];
165 Int_t crate = TMath::Nint( trigdef[base+2] );
166 Int_t slot = TMath::Nint( trigdef[base+3] );
167 Int_t chan = TMath::Nint( trigdef[base+4] );
168 fTrgTypes.push_back( itrg );
169 fToffsets.push_back( toff );
170 fDetMap->AddModule( crate, slot, chan, chan, itrg );
171 }
172 assert( ndef == fDetMap->GetSize() );
173 assert( ndef == fTrgTypes.size() );
174 assert( ndef == fToffsets.size() );
175
176 // Use the TDC mode flag to set the sign of fTDCRes (see comments in Process)
177 fTDCRes = (fCommonStop ? 1.0 : -1.0) * TMath::Abs(fTDCRes);
178
179 // now construct the appropriate arrays
180 fTrgTimes.resize( fTrgTypes.size() );
181
182 return kOK;
183}
184
185//____________________________________________________________________________
187{
188 // Define/delete event-by-event global variables
189
190 RVarDef vars[] = {
191 { "evtype", "Earliest trg-bit for the event", "fEvtType" },
192 { "trgtimes","Times for each trg-type", "fTrgTimes" },
193 { nullptr }
194 };
195
196 return DefineVarsFromList( vars, mode );
197}
198
199//____________________________________________________________________________
int Int_t
unsigned int UInt_t
const Data_t kBig
Definition DataType.h:15
uint32_t chan
#define d(i)
const UInt_t kMaxUInt
double Double_t
const char Option_t
Option_t Option_t TPoint TPoint const char mode
char name[80]
static const char *const here
Definition THaVar.cxx:64
virtual Int_t ReadDatabase(const TDatime &date)
static Int_t LoadDB(FILE *file, const TDatime &date, const DBRequest *request, const char *prefix, Int_t search=0, const char *here="THaAnalysisObject::LoadDB")
static Int_t DefineVarsFromList(const void *list, EType type, EMode mode, const char *def_prefix, const TObject *obj, const char *prefix, const char *here, const char *comment_subst="")
virtual const char * Here(const char *) const
virtual FILE * OpenFile(const TDatime &date)
virtual Int_t AddModule(UInt_t crate, UInt_t slot, UInt_t chan_lo, UInt_t chan_hi, UInt_t first=0, Int_t model=0, Int_t refindex=-1, Int_t refchan=-1, UInt_t plane=0, UInt_t signal=0)
UInt_t GetSize() const
Definition THaDetMap.h:125
void Clear()
Definition THaDetMap.h:117
Module * GetModule(UInt_t i) const
Definition THaDetMap.h:248
UInt_t GetNumHits(UInt_t crate, UInt_t slot, UInt_t chan) const
Definition THaEvData.h:264
UInt_t GetData(UInt_t crate, UInt_t slot, UInt_t chan, UInt_t hit) const
Definition THaEvData.h:273
virtual Int_t Process(const THaEvData &)
std::vector< Double_t > fToffsets
virtual void Clear(Option_t *opt="")
virtual Int_t DefineVariables(EMode mode=kDefine)
virtual ~THaTriggerTime()
THaTriggerTime(const char *name="trg", const char *description="")
std::vector< Double_t > fTrgTimes
THaDetMap * fDetMap
virtual Int_t ReadDatabase(const TDatime &date)
std::vector< Int_t > fTrgTypes
virtual void Error(const char *method, const char *msgfmt,...) const
Int_t Nint(T x)
Double_t Abs(Double_t d)
STL namespace.
ClassImp(TPyArg)