Hall C ROOT/C++ Analyzer (hcana)
Loading...
Searching...
No Matches
THcHelicityReader.cxx
Go to the documentation of this file.
1//*-- Author : Ole Hansen August 2006
2// Extracted from Bob Michaels' THaHelicity CVS 1.19
4//
5// THcHelicityReader
6//
8
9#include "THcHelicityReader.h"
10#include "THaEvData.h"
11#include "THcGlobals.h"
12#include "THcRun.h"
13#include "THcParmList.h"
14#include "TMath.h"
15#include "TError.h"
16#include "VarDef.h"
17#include "THaAnalysisObject.h" // For LoadDB
18#include <iostream>
19#include <vector>
20#include "TH1F.h"
21
22#define DEFAULT_FADCTHRESHOLD 8000
23
24using namespace std;
25
26//____________________________________________________________________
28 : fTITime(0), fTITime_last(0), fTITime_rollovers(0), fFADCModule(0),
29 fTTimeDiff(0), fHaveROCs(kFALSE)
30{
31 // Default constructor
32
33}
34//____________________________________________________________________
36{
37 // Destructor
38
39 // Histograms will be deleted by ROOT
40 // for( Int_t i = 0; i < NHISTR; ++i ) {
41 // delete fHistoR[i];
42 // }
43}
44
45//____________________________________________________________________
52
53//_____________________________________________________________________________
54Int_t THcHelicityReader::ReadDatabase( const char* /*dbfilename*/,
55 const char* /*prefix*/,
56 const TDatime& /*date*/,
57 int /*debug_flag*/ )
58{
59
60 // Eventually get these from the parameter file
61
62 // SHMS settings see https://logbooks.jlab.org/entry/3614445
63 cout << "THcHelicityReader: Helicity information from ROC 2 (SHMS)" << endl;
64 SetROCinfo(kHel,2,14,9);
65 SetROCinfo(kHelm,2,14,8);
66 SetROCinfo(kMPS,2,14,10);
67 SetROCinfo(kQrt,2,14,7); // Starting about run 5818
68 SetROCinfo(kTime,2,21,2);
69
72
73 DBRequest list[] = {
74 {"helicity_adcthreshold",&fADCThreshold, kInt, 0, 1},
75 {"helicity_adcrawsamples",&fADCRawSamples, kInt, 0, 1},
76 {0}
77 };
78
79 gHcParms->LoadParmValues(list, "");
80
82}
83
84//_____________________________________________________________________________
86{
87 // static const char* const here = "THcHelicityReader::Begin";
88 // cout<<here<<endl;
89
90 fTITime_last = 0;
91 fTITime = 0;
93
94 return;
95}
96//____________________________________________________________________
98{
99 // static const char* const here = "THcHelicityReader::End";
100 // cout<<here<<endl;
101
102 return;
103}
104//____________________________________________________________________
106{
107 // Obtain the present data from the event for QWEAK helicity mode.
108
109 static const char* here = "THcHelicityReader::ReadData";
110
111 // std::cout<<" kHel, kTime, kRing="<< kHel<<" "<<kTime<<" "<<kRing<<endl;
112 // for (int jk=0; jk<3; jk++)
113 // {
114 // std::cout<<" which="<<jk
115 // <<" roc="<<fROCinfo[jk].roc
116 // <<" header="<<fROCinfo[jk].header
117 // <<" index="<<fROCinfo[jk].index
118 // <<endl;
119 // }
120
121 // std::cout<<" fHaveROCs="<<fHaveROCs<<endl;
122 if( !fHaveROCs ) {
123 ::Error( here, "ROC data (detector map) not properly set up." );
124 return -1;
125 }
126
127 // Check if ROC info is correct
128
129 if(!evdata.GetModule(fROCinfo[kTime].roc, fROCinfo[kTime].slot)) {
130 // cout << "ROC, Slot: " << fROCinfo[kTime].roc << " " << fROCinfo[kTime].slot << endl;
131 // cout << "Evtype = " << evdata.GetEvType() << endl;
132 // cout << "THcHelicityReader: ROC 2 not found" << endl;
133 //cout << "Changing to ROC 1 (HMS)" << endl;
134 SetROCinfo(kHel,1,18,9);
135 SetROCinfo(kHelm,1,18,8);
136 SetROCinfo(kMPS,1,18,10);
137 SetROCinfo(kQrt,1,18,7);
138 SetROCinfo(kTime,1,21,2);
139 }
140
141 auto coda_version=(gHaRun) ? gHaRun->GetDataVersion() :3;
142
143 // Get the TI Data
144 ULong64_t titime;
145
146 if(coda_version<3){
147 titime = (ULong64_t) evdata.GetData(fROCinfo[kTime].roc,
149 fROCinfo[kTime].index, 0);
150 }else{
151 titime = (ULong64_t) evdata.GetEvTime();
152 }
153
154 if(!fFADCModule) {
155 fFADCModule = dynamic_cast<Decoder::Fadc250Module*>
156 (evdata.GetModule(fROCinfo[kHel].roc, fROCinfo[kHel].slot));
158 }
159
160 if(titime == 0) {
161 cout << "Event " << evdata.GetEvNum() << " TI Trigger time missing." << endl;
162 cout << "Event " << evdata.GetEvNum() << " TI Trigger time missing. Using FADC trigger time" << endl;
163 if(coda_version<3){
164 cout << " kTime = " << kTime << " TI Roc slot : " << fROCinfo[kTime].roc << " " << fROCinfo[kTime].slot << " " << fROCinfo[kTime].index << " " << evdata.GetData(fROCinfo[kTime].roc,fROCinfo[kTime].slot, 0, 0)
165 << " " << evdata.GetData(fROCinfo[kTime].roc,fROCinfo[kTime].slot, 1, 0) << " " << evdata.GetData(fROCinfo[kTime].roc,fROCinfo[kTime].slot, 2, 0)<< endl;
166 }
168 }
169
170 // cout << "Time " << titime << " " << fTITime_last << " " <<
171 // evdata.GetData(1,21,2,0) << " " <<
172 // evdata.GetData(2,21,2,0) << " " << fFADCModule->GetTriggerTime() << endl;
173
174 // Check again if ROC info is correct
175 // cout << "Evtype = " << evdata.GetEvType() << endl;
176#if 0
177 if(titime == 0 && fTITime_last==0) {
178 cout << "B ROC, Slot: " << fROCinfo[kTime].roc << " " << fROCinfo[kTime].slot << endl;
179 cout << "THcHelicityReader: ROC 2 not found" << endl;
180 cout << "Changing to ROC 1 (HMS)" << endl;
181 SetROCinfo(kHel,1,18,9);
182 SetROCinfo(kHelm,1,18,8);
183 SetROCinfo(kMPS,1,18,10);
184 SetROCinfo(kQrt,1,18,7);
185 SetROCinfo(kTime,1,21,2);
186 titime = (UInt_t) evdata.GetData(fROCinfo[kTime].roc,
188 fROCinfo[kTime].index, 0);
189 }
190#endif
191
192 // cout << fTITime_last << " " << titime << endl;
193 if(titime < fTITime_last && coda_version<3) {
195 }
196
197 fTITime = titime + fTITime_rollovers*TMath::Power(2,32);
198 fTITime_last = titime;
199
200 if(coda_version<3)const_cast<THaEvData&>(evdata).SetEvTime(fTITime);
201
202 // Get the helicity control signals. These are from the pedestals
203 // acquired by FADC channels. If helpraw and helmraw both
204 // come back zero, assume that we have raw sample data instead
205 // and switch to using raw sample data for the rest of the run.
206 // If the FADC threshold has not be set by the user, set the threshold
207 // to the default threshold divided by 4.
208
209 Int_t helpraw, helmraw, mpsraw, qrtraw;
210 if(!fADCRawSamples) {
211 helpraw = evdata.GetData(Decoder::kPulsePedestal,
212 fROCinfo[kHel].roc,
213 fROCinfo[kHel].slot,
214 fROCinfo[kHel].index, 0);
215 helmraw = evdata.GetData(Decoder::kPulsePedestal,
216 fROCinfo[kHelm].roc,
217 fROCinfo[kHelm].slot,
218 fROCinfo[kHelm].index, 0);
219 mpsraw = evdata.GetData(Decoder::kPulsePedestal,
220 fROCinfo[kMPS].roc,
221 fROCinfo[kMPS].slot,
222 fROCinfo[kMPS].index, 0);
223 qrtraw = evdata.GetData(Decoder::kPulsePedestal,
224 fROCinfo[kQrt].roc,
225 fROCinfo[kQrt].slot,
226 fROCinfo[kQrt].index, 0);
227 if(helpraw <= 0 && helmraw <=0) {
228 fADCRawSamples = 1;
231 }
232 }
233 }
234 if(fADCRawSamples) {
235 helpraw = evdata.GetData(Decoder::kSampleADC,
236 fROCinfo[kHel].roc,
237 fROCinfo[kHel].slot,
238 fROCinfo[kHel].index, 0);
239 helmraw = evdata.GetData(Decoder::kSampleADC,
240 fROCinfo[kHelm].roc,
241 fROCinfo[kHelm].slot,
242 fROCinfo[kHelm].index, 0);
243 mpsraw = evdata.GetData(Decoder::kSampleADC,
244 fROCinfo[kMPS].roc,
245 fROCinfo[kMPS].slot,
246 fROCinfo[kMPS].index, 0);
247 qrtraw = evdata.GetData(Decoder::kSampleADC,
248 fROCinfo[kQrt].roc,
249 fROCinfo[kQrt].slot,
250 fROCinfo[kQrt].index, 0);
251 }
252
253 fIsQrt = qrtraw > fADCThreshold;
254 fIsMPS = mpsraw > fADCThreshold;
255 fIsHelp = helpraw > fADCThreshold;
256 fIsHelm = helmraw > fADCThreshold;
257
258 // cout << helpraw << " " << helmraw << " " << mpsraw << " " << qrtraw << endl;
259 // cout << fADCThreshold << " " << fADCRawSamples << " " << fIsQrt << " " << fIsMPS << " " << fIsHelp << " " << fIsHelm << endl;
260
261 return 0;
262}
263
264//TODO: this should not be needed once LoadDB can fill fROCinfo directly
265//____________________________________________________________________
267 Int_t slot, Int_t index )
268{
269
270 // Define source and offset of data. Normally called by ReadDatabase
271 // of the detector that is a THcHelicityReader.
272 //
273 // "which" is one of { kHel, kKelm, kQrt, kTime }.
274 // You must define at least the kHel and kTime ROCs.
275 // Returns <0 if parameter error, 0 if success
276
277 if( which<kHel || which>=kCount )
278 return -1;
279 if( roc <= 0 || roc > 255 )
280 return -2;
281
282 fROCinfo[which].roc = roc;
283 fROCinfo[which].slot = slot;
284 fROCinfo[which].index = index;
285
286 //cout << "SetROCInfo: " << which << " " << fROCinfo[which].roc << " " << fROCinfo[which].slot <<
287 // " " << fROCinfo[which].index << endl;
289
290 return 0;
291}
292
293//____________________________________________________________________
295
int Int_t
unsigned int UInt_t
const Bool_t kFALSE
const char Option_t
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t index
R__EXTERN class THaRunBase * gHaRun
ClassImp(VDC::AnalyticTTDConv) using namespace std
R__EXTERN class THcParmList * gHcParms
Definition THcGlobals.h:11
#define DEFAULT_FADCTHRESHOLD
virtual UInt_t GetTriggerTime() const
virtual ULong64_t GetEvTime() const
UInt_t GetData(Decoder::EModuleType type, UInt_t crate, UInt_t slot, UInt_t chan, UInt_t hit) const
UInt_t GetEvNum() const
virtual Decoder::Module * GetModule(UInt_t roc, UInt_t slot) const
virtual Int_t GetDataVersion()
ROCinfo fROCinfo[kCount]
virtual void Clear(Option_t *opt="")
Int_t SetROCinfo(EROC which, Int_t roc, Int_t slot, Int_t index)
Decoder::Fadc250Module * fFADCModule
virtual Int_t ReadData(const THaEvData &evdata)
Int_t ReadDatabase(const char *dbfilename, const char *prefix, const TDatime &date, int debug_flag=0)
Int_t LoadParmValues(const DBRequest *list, const char *prefix="")
Retrieve parameter values from the parameter cache.
unsigned long long ULong64_t
Double_t Power(Double_t x, Double_t y)
STL namespace.