Hall C ROOT/C++ Analyzer (hcana)
Loading...
Searching...
No Matches
THcRaster.cxx
Go to the documentation of this file.
1
13#include "TMath.h"
14
15#include "THcRaster.h"
16#include "THaEvData.h"
17#include "THaDetMap.h"
18#include "THcAnalyzer.h"
19
20#include "THcParmList.h"
21#include "THcGlobals.h"
22#include "THaGlobals.h"
23#include "THaApparatus.h"
24#include "THcRawAdcHit.h"
25#include "THcSignalHit.h"
26
27//#include "THcHitList.h"
28
29#include <cstring>
30#include <cstdio>
31#include <cstdlib>
32#include <iostream>
33#include <fstream>
34
35
36using namespace std;
37
38//_____________________________________________________________________________
39THcRaster::THcRaster( const char* name, const char* description,
40 THaApparatus* apparatus ) :
41 THaBeamDet(name,description,apparatus)
42{
43
46 FRXA_rawadc = 0;
47 FRXB_rawadc = 0;
48 FRYA_rawadc = 0;
49 FRYB_rawadc = 0;
50 fXA_ADC = 0;
51 fYA_ADC = 0;
52 fXB_ADC = 0;
53 fYB_ADC = 0;
54 fXA_pos = 0;
55 fYA_pos = 0;
56 fXB_pos = 0;
57 fYB_pos = 0;
58 for (Int_t nt=0;nt<4;nt++) {
59 fXbeam_prev[nt] = 0;
60 fYbeam_prev[nt] = 0;
61 }
62 fXpbeam_prev = 0;
63 fYpbeam_prev = 0;
64 BPMXA_raw = 0;
65 BPMYA_raw = 0;
66 BPMXB_raw = 0;
67 BPMYB_raw = 0;
68 BPMXC_raw = 0;
69 BPMYC_raw = 0;
70 BPMXA_pos = 0;
71 BPMYA_pos = 0;
72 BPMXB_pos = 0;
73 BPMYB_pos = 0;
74 BPMXC_pos = 0;
75 BPMYC_pos = 0;
76 fFrCalMom = 0;
77 fFrXA_ADCperCM = 1.0;
78 fFrYA_ADCperCM = 1.0;
79 fFrXB_ADCperCM = 1.0;
80 fFrYB_ADCperCM = 1.0;
85
89 fEbeamEpics=0.;
90 for(Int_t i=0;i<4;i++){
91
92 fPedADC[i] = 0;
93 }
94
95 InitArrays();
96}
97
98//_____________________________________________________________________________
100 THaBeamDet("THcRaster") // no default constructor available
101{
102
103 frPosAdcPulseIntRaw = NULL;
104
105 InitArrays();
106
107
108}
109//_____________________________________________________________________________
111{
112
113// Destructor
114
116
117 DeleteArrays();
118}
119
120//_____________________________________________________________________________
121
123{
124 //cout << "In THcRaster::Init()" << endl;
125
126 char EngineDID[] = "xRASTER";
127 EngineDID[0] = toupper(GetApparatus()->GetName()[0]);
128
129
130 // Fill detector map with RASTER type channels
131 if( gHcDetectorMap->FillMap(fDetMap, EngineDID) < 0 ) {
132 static const char* const here = "Init()";
133 Error( Here(here), "Error filling detectormap for %s.", EngineDID);
134
135 return kInitError;
136 }
137
138 InitHitList(fDetMap,"THcRasterRawHit",fDetMap->GetTotNumChan()+1);
139
140 EStatus status;
141 if( (status = THaBeamDet::Init( date )) )
142 return fStatus=status;
143
144 return fStatus = kOK;
145
146}
147
148//_____________________________________________________________________________
150
151{
152
153 //cout << "THcRaster::ReadDatabase()" << endl;
154
155 // Read parameters such as calibration factor, of this detector from the database.
156 //cout << "THcRaster::ReadDatabase()" << endl;
157
158 char prefix[2];
159
160 //cout << " THcRaster::ReadDatabase GetName() called " << GetName() << endl;
161 // prefix[0]=tolower(GetName()[0]);
162 // bpw- The prefix is hardcoded so that we don't have to change the gbeam.param file. o/w to get the following variables, we need to change to parameter names to rfr_cal_mom, etc where "r" comes from prefix[0]=tolower(GetName()[0]).
163 prefix[0]='g';
164 prefix[1]='\0';
165
166
167
168 // string names;
169 DBRequest list[]={
170 {"fr_cal_mom",&fFrCalMom, kDouble},
171 {"frxa_adcpercm",&fFrXA_ADCperCM, kDouble},
172 {"frya_adcpercm",&fFrYA_ADCperCM, kDouble},
173 {"frxb_adcpercm",&fFrXB_ADCperCM, kDouble},
174 {"fryb_adcpercm",&fFrYB_ADCperCM, kDouble},
175 {"frxa_adc_zero_offset",&fFrXA_ADC_zero_offset,kDouble},
176 {"frya_adc_zero_offset",&fFrYA_ADC_zero_offset,kDouble},
177 {"frxb_adc_zero_offset",&fFrXB_ADC_zero_offset,kDouble},
178 {"fryb_adc_zero_offset",&fFrYB_ADC_zero_offset,kDouble},
179 {"pbeam",&fgpbeam, kDouble,0,1},
180 {"frx_dist", &fgfrx_dist, kDouble},
181 {"fry_dist", &fgfry_dist, kDouble},
182 {"beam_x", &fgbeam_xoff, kDouble,0,1},
183 {"beam_xp", &fgbeam_xpoff, kDouble,0,1},
184 {"beam_y", &fgbeam_yoff, kDouble,0,1},
185 {"beam_yp", &fgbeam_ypoff, kDouble,0,1},
186 {"bpmxa_slope", &fgbpmxa_slope, kDouble,0,1},
187 {"bpmxa_off", &fgbpmxa_off, kDouble,0,1},
188 {"bpmxb_slope", &fgbpmxb_slope, kDouble,0,1},
189 {"bpmxb_off", &fgbpmxb_off, kDouble,0,1},
190 {"bpmxc_slope", &fgbpmxc_slope, kDouble,0,1},
191 {"bpmxc_off", &fgbpmxc_off, kDouble,0,1},
192 {"bpmya_slope", &fgbpmya_slope, kDouble,0,1},
193 {"bpmya_off", &fgbpmya_off, kDouble,0,1},
194 {"bpmyb_slope", &fgbpmyb_slope, kDouble,0,1},
195 {"bpmyb_off", &fgbpmyb_off, kDouble,0,1},
196 {"bpmyc_slope", &fgbpmyc_slope, kDouble,0,1},
197 {"bpmyc_off", &fgbpmyc_off, kDouble,0,1},
198 {"bpma_zpos", &fgbpma_zpos, kDouble,0,1},
199 {"bpmb_zpos", &fgbpmb_zpos, kDouble,0,1},
200 {"bpmc_zpos", &fgbpmc_zpos, kDouble,0,1},
201 {"usefr", &fgusefr, kInt,0,1},
202 {0}
203 };
204 fgpbeam = 0.001;
205 //
206 //positions of BPMs relative to target (from Fall 2023 survey)
207 fgbpma_zpos = 320.22;
208 fgbpmb_zpos = 224.62 ;// cm
209 fgbpmc_zpos = 124.36 ;// cm
210 // Default offsets to zero and slopes to +/- 1
211 fgbeam_xoff = -999.;
212 fgbeam_xpoff =-999.;
213 fgbeam_yoff = -999.;
214 fgbeam_ypoff =-999.;
215 //
216 fgbpmxa_off = 0.0;
217 fgbpmxa_slope = -1.0;
218 fgbpmxb_off = 0.0;
219 fgbpmxb_slope = -1.0;
220 fgbpmxc_off = 0.0;
221 fgbpmxc_slope = -1.0;
222 fgbpmya_off = 0.0;
223 fgbpmya_slope = 1.0;
224 fgbpmyb_off = 0.0;
225 fgbpmyb_slope = 1.0;
226 fgbpmyc_off = 0.0;
227 fgbpmyc_slope = 1.0;
228 fgusefr = 0;
229
230
231 // get the calibration factors from gbeam.param file
232 gHcParms->LoadParmValues((DBRequest*)&list,prefix);
233
234 frPosAdcPulseIntRaw = new TClonesArray("THcSignalHit", 4);
235
236 THcAnalyzer *analyzer = dynamic_cast<THcAnalyzer*>(THcAnalyzer::GetInstance());
237 fEpicsHandler = analyzer->GetEpicsEvtHandler();
238
239 //
241 if (fgbeam_xoff ==-999. && fgbeam_yoff ==-999.) {
243 cout << " THcRaster is using EPICS bpm for beam position" << endl;
244 } else {
245 if (fgbeam_xoff ==-999.) fgbeam_xoff = 0.;
246 if (fgbeam_yoff ==-999.) fgbeam_yoff = 0.;
247 if (fgbeam_xpoff ==-999.) fgbeam_xpoff = 0.;
248 if (fgbeam_ypoff ==-999.) fgbeam_ypoff = 0.;
249 cout << " THcRaster is using parameters for beam position" << " x = " << fgbeam_xoff<< " y = " << fgbeam_yoff<< endl;
250 cout << " THcRaster is using parameters for beam angle position" << " xp = " << fgbeam_xpoff<< " yp = " << fgbeam_ypoff<< endl;
251
252 }
253
254
255 return kOK;
256
257}
258
259
260//_____________________________________________________________________________
262{
263 // Initialize global variables for histogramming and tree
264
265 //cout << "THcRaster::DefineVariables called " << GetName() << endl;
266
267 if( mode == kDefine && fIsSetup ) return kOK;
268 fIsSetup = ( mode == kDefine );
269
270 // Register variables in global list
271
272 RVarDef vars[] = {
273 {"frxaRawAdc", "Raster XA raw ADC", "FRXA_rawadc"},
274 {"fryaRawAdc", "Raster YA raw ADC", "FRYA_rawadc"},
275 {"frxbRawAdc", "Raster XB raw ADC", "FRXB_rawadc"},
276 {"frybRawAdc", "Raster YB raw ADC", "FRYB_rawadc"},
277 {"frxa_adc", "Raster XA ADC", "fXA_ADC"},
278 {"frya_adc", "Raster YA ADC", "fYA_ADC"},
279 {"frxb_adc", "Raster XB ADC", "fXB_ADC"},
280 {"fryb_adc", "Raster YB ADC", "fYB_ADC"},
281 {"fr_xa", "Raster XA position", "fXA_pos"},
282 {"fr_ya", "Raster YA position", "fYA_pos"},
283 {"fr_xb", "Raster XB position", "fXB_pos"},
284 {"fr_yb", "Raster YB position", "fYB_pos"},
285 {"fr_xbpm_tar", "X BPM at target (+X is beam right)", "fXbpm_tar"},
286 {"fr_ybpm_tar", "Y BPM at target (+Y is up)", "fYbpm_tar"},
287 {"fr_xbpmA", "X BPM at BPMA (+X is beam right)", "fXbpm_A"},
288 {"fr_ybpmA", "Y BPM at BPMA (+Y is up)", "fYbpm_A"},
289 {"fr_xbpmB", "X BPM at BPMB (+X is beam right)", "fXbpm_B"},
290 {"fr_ybpmB", "Y BPM at BPMB (+Y is up)", "fYbpm_B"},
291 {"fr_xbpmC", "X BPM at BPMC (+X is beam right)", "fXbpm_C"},
292 {"fr_ybpmC", "Y BPM at BPMC (+Y is up)", "fYbpm_C"},
293 {"ebeam_epics", "Beam energy of epics variable HALLC:p", "fEbeamEpics"},
294 { 0 }
295 };
296
297 return DefineVarsFromList( vars, mode );
298}
299
300//_____________________________________________________________________________
301
302inline
304{
305
306 fNhits = 0;
307
309}
310//_____________________________________________________________________________
312{
313 /*
314 Extract data from the hit list, accumulating into arrays for
315 calculating pedestals.
316 From ENGINE/g_analyze_misc.f -
317
318 * JRA: Code to check FR pedestals. Since the raster is a fixed frequency
319 * and the pedestals come at a fixed rate, it is possible to keep getting
320 * the same value for each pedestal event, and get the wrong zero value.
321 * (see HCLOG #28325). So calculate pedestal from first 1000 REAL
322 * events and compare to value from pedestal events. Error on each
323 * measurement is RMS/sqrt(1000), error on diff is *sqrt(2), so 3 sigma
324 * check is 3*sqrt(2)*RMS/sqrt(1000) = .13*RMS
325 !
326 ! Can't use RMS, since taking sum of pedestal**2 for these signals
327 ! gives rollover for integer*4. Just assume signal is +/-2000
328 ! channels, gives sigma of 100 channels, so check for diff>130.
329 !
330 */
331 //cout << "THcRaster::AccumulatePedestels()" << endl;
332
333 Int_t nrawhits = rawhits->GetLast()+1;
334
335 Int_t ihit = 0;
336 UInt_t nrPosAdcHits=0;
337
338 while(ihit < nrawhits) {
340 THcRawAdcHit& rawPosAdcHit = hit->GetRawAdcHitPos();
341 Int_t nsig = hit->fCounter;
342
343 for (UInt_t thit=0; thit<rawPosAdcHit.GetNPulses(); ++thit) {
344 ((THcSignalHit*) frPosAdcPulseIntRaw->ConstructedAt(nrPosAdcHits))->Set(nsig, rawPosAdcHit.GetPulseIntRaw(thit));
345 ++nrPosAdcHits;
346 }
347 ihit++;
348 }
349
350
351 for(Int_t ielem = 0; ielem < frPosAdcPulseIntRaw->GetEntries(); ielem++) {
352 Int_t nraster = ((THcSignalHit*) frPosAdcPulseIntRaw->ConstructedAt(ielem))->GetPaddleNumber() - 1;
353 Double_t pulseIntRaw = ((THcSignalHit*) frPosAdcPulseIntRaw->ConstructedAt(ielem))->GetData();
354 if (nraster ==0) fPedADC[0] = pulseIntRaw;
355 if (nraster ==1) fPedADC[1] = pulseIntRaw;
356 if (nraster ==2) fPedADC[2] = pulseIntRaw;
357 if (nraster ==3) fPedADC[3] = pulseIntRaw;
358 }
359
360
361}
362
363
364//_____________________________________________________________________________
366{
367 /*
368 Use the accumulated pedestal data to calculate pedestals
369 From ENGINE/g_analyze_misc.f -
370
371 if (numfr.eq.1000) then
372 avefrx = sumfrx / float(numfr)
373 avefry = sumfry / float(numfr)
374 if (abs(avefrx-gfrx_adc_ped).gt.130.) then
375 write(6,*) 'FRPED: peds give <frx>=',gfrx_adc_ped,
376 $ ' realevents give <frx>=',avefrx
377 endif
378 if (abs(avefry-gfry_adc_ped).gt.130.) then
379 write(6,*) 'FRPED: peds give <fry>=',gfry_adc_ped,
380 $ ' realevents give <fry>=',avefry
381 endif
382 endif
383 */
384 //cout << "THcRaster::CalculatePedestels()" << endl;
385
390
391
392}
393
394
395//_____________________________________________________________________________
397{
398
399 //cout << "THcRaster::Decode()" << endl;
400 // Get the Hall C style hitlist (fRawHitList) for this event
401
402 fNhits = DecodeToHitList(evdata);
403
404
405 // Get the pedestals from the first 1000 events
406 //if(fNPedestalEvents < 10)
407 if(gHaCuts->Result("Pedestal_event") && (fNPedestalEvents < 1000)) {
409 fAnalyzePedestals = 1; // Analyze pedestals first normal events
411
412 return(0);
413 }
414
417 fAnalyzePedestals = 0; // Don't analyze pedestals next event
418 }
419
420
421 Int_t ihit = 0;
422 UInt_t nrPosAdcHits=0;
423
424 while(ihit < fNhits) {
426 THcRawAdcHit& rawPosAdcHit = hit->GetRawAdcHitPos();
427 Int_t nsig = hit->fCounter;
428
429 for (UInt_t thit=0; thit<rawPosAdcHit.GetNPulses(); ++thit) {
430 ((THcSignalHit*) frPosAdcPulseIntRaw->ConstructedAt(nrPosAdcHits))->Set(nsig, rawPosAdcHit.GetPulseIntRaw(thit));
431 ++nrPosAdcHits;
432 }
433 if (rawPosAdcHit.GetNPulses()==0 &&rawPosAdcHit.GetNSamples()>0 ) {
434 Int_t NSA= rawPosAdcHit.GetF250_NSA();
435 UInt_t LS = 0;
436 UInt_t HS = NSA;
437 Int_t rawdata = rawPosAdcHit.GetIntegral(LS,HS);
438 ((THcSignalHit*) frPosAdcPulseIntRaw->ConstructedAt(nrPosAdcHits))->Set(nsig,rawdata );
439 ++nrPosAdcHits;
440 }
441
442 ihit++;
443 }
444
445 for(Int_t ielem = 0; ielem < frPosAdcPulseIntRaw->GetEntries(); ielem++) {
446 Int_t nraster = ((THcSignalHit*) frPosAdcPulseIntRaw->ConstructedAt(ielem))->GetPaddleNumber() - 1;
447 Double_t pulseIntRaw = ((THcSignalHit*) frPosAdcPulseIntRaw->ConstructedAt(ielem))->GetData();
448 if (nraster ==0) FRYA_rawadc = pulseIntRaw;
449 if (nraster ==1) FRXA_rawadc = pulseIntRaw;
450 if (nraster ==2) FRYB_rawadc = pulseIntRaw;
451 if (nraster ==3) FRXB_rawadc = pulseIntRaw;
452 }
453
454 if (fEpicsHandler) {
455 if (fEpicsHandler->IsLoaded("IPM3H07A.XRAW")){
456 BPMXA_raw = atof(fEpicsHandler->GetString("IPM3H07A.XRAW"));
457 }
458 if (fEpicsHandler->IsLoaded("IPM3H07A.YRAW")){
459 BPMYA_raw = atof(fEpicsHandler->GetString("IPM3H07A.YRAW"));
460 }
461 if (fEpicsHandler->IsLoaded("IPM3H07B.XRAW")){
462 BPMXB_raw = atof(fEpicsHandler->GetString("IPM3H07B.XRAW"));
463 }
464 if (fEpicsHandler->IsLoaded("IPM3H07B.YRAW")){
465 BPMYB_raw = atof(fEpicsHandler->GetString("IPM3H07B.YRAW"));
466 }
467 if (fEpicsHandler->IsLoaded("IPM3H07C.XRAW")){
468 BPMXC_raw = atof(fEpicsHandler->GetString("IPM3H07C.XRAW"));
469 }
470 if (fEpicsHandler->IsLoaded("IPM3H07C.YRAW")){
471 BPMYC_raw = atof(fEpicsHandler->GetString("IPM3H07C.YRAW"));
472 }
473 if (fEpicsHandler->IsLoaded("HALLC:p")){
474 fEbeamEpics_read = atof(fEpicsHandler->GetString("HALLC:p"));
475 }
476 }
477
478 return 0;
479
480}
481
482//_____________________________________________________________________________
484
485 //cout << "In THcRaster::Process()" << endl;
486
487 /*
488 calculate raster position from ADC value.
489 From ENGINE/g_analyze_misc.f -
490
491 gfrx_adc = gfrx_raw_adc - gfrx_adc_ped
492 gfry_adc = gfry_raw_adc - gfry_adc_ped
493 */
494
495 //std::cout << "Beam energy = " << fgpbeam << std::endl;
496 //std::cout << "Raster Calibration Momentum = " << fFrCalMom << std::endl;
497 //std::cout << "Raster X calibration per cm = " << fFrXA_ADCperCM << std::endl;
498 //std::cout << "Raster Y calibration per cm = " << fFrYA_ADCperCM << std::endl;
499 //std::cout << "Raster Offsets: " << "XA = " << fFrXA_ADC_zero_offset << " YA = " << fFrYA_ADC_zero_offset<< std::endl;
500 //std::cout << "Raster Offsets: " << "XB = " << fFrXB_ADC_zero_offset << " YB = " << fFrYB_ADC_zero_offset<< std::endl;
501 // calculate the raster currents
506
507
508 /*
509 calculate the raster positions
510
511 gfrx = (gfrx_adc/gfrx_adcpercm)*(gfr_cal_mom/ebeam)
512 gfry = (gfry_adc/gfry_adcpercm)*(gfr_cal_mom/ebeam)
513 */
514
519
520
521 Bool_t checkBPM = (BPMXA_raw == 0) && (BPMYA_raw == 0) && (BPMXB_raw == 0) &&
522 (BPMYB_raw == 0) && (BPMXC_raw ==0) && (BPMYC_raw == 0);
523 //cout << "BPM's working? " << checkBPM << endl;
524
525 // Factor of 0.1 is to convert from mm to cm
526 // BPM's are typically calibrated in mm, whereas the standard distances are in cm
527 // in the analyzer.
528 //
535
536 //cout << "BPMXA_pos = " << BPMXA_pos << endl;
537 //cout << "BPMYA_pos = " << BPMYA_pos << endl;
538 //cout << "BPMXB_pos = " << BPMXB_pos << endl;
539 //cout << "BPMYB_pos = " << BPMYB_pos << endl;
540 //cout << "BPMXC_pos = " << BPMXC_pos << endl;
541 //cout << "BPMYC_pos = " << BPMYC_pos << endl;
542
543 // Calculate position and direction at target from BPM values
544 // Use the A and C BPM information, as these are located downstream of the raster
545 // magnets. If there is no BPM information available, assume zero offsets.
546 //
547 Double_t xbeam;
548 Double_t ybeam;
549 Double_t xpbeam;
550 Double_t ypbeam;
551 if (!checkBPM){
554 xpbeam = (xbeam-BPMXA_pos)/fgbpma_zpos;
555 ypbeam = (ybeam-BPMYA_pos)/fgbpma_zpos;
556 fXbeam_prev[0]=xbeam;
557 fYbeam_prev[0]=ybeam;
558 fXpbeam_prev=xpbeam;
559 fYpbeam_prev=ypbeam;
568 }else{
569 xbeam = fXbeam_prev[0];
570 ybeam = fYbeam_prev[0];
571 xpbeam = fXpbeam_prev;
572 ypbeam = fYpbeam_prev;
580 }
581
582
584 fXbpm_tar= -xbeam; // assumes the BPM calibration gives BPM with +X beam left (HARP coordinate system)
585 fYbpm_tar= ybeam;
586 fXpbpm_tar= -xpbeam; // assumes the BPM calibration gives BPM with +X beam left (HARP coordinate system)
587 fYpbpm_tar= ypbeam;
588
589 } else {
590 fXbpm_tar= fgbeam_xoff; // +X beam right (EPICS coordinate system)
592 fXpbpm_tar= fgbeam_xpoff; // +X beam right (EPICS coordinate system)
594 }
595
596
597 fXbpm_A= -fXbeam_prev[1];
599 fXbpm_B= -fXbeam_prev[2];
601 fXbpm_C= -fXbeam_prev[3];
603 //std::cout<<" XA = "<<fXA_pos<<" YA = "<<fYA_pos<<std::endl;
604 //std::cout<<" XB = "<<fXB_pos<<" YB = "<<fYB_pos<<std::endl;
605
606 //std::cout << "Use FR? " << fgusefr << std::endl;
607
608 //std::cout << "BPM X at Z=0 -> " << fgbeam_xoff << std::endl;
609 //std::cout << "BPM Y at Z=0 -> " << fgbeam_yoff << std::endl;
610 //std::cout << "BPM XP at Z=0 -> " << fgbeam_xpoff << std::endl;
611 //std::cout << "BPM YP at Z=0 -> " << fgbeam_ypoff << std::endl;
612
613 //std::cout << "Raster X distance = " << fgfrx_dist << std::endl;
614 //std::cout << "Raster Y distance = " << fgfry_dist << std::endl;
615
617 Double_t tp = fYpbpm_tar;
618
619 if(fgusefr != 0) {
625 fDirection.SetXYZ(tt, tp ,1.0); // Set arbitrarily to avoid run time warnings
626 fDirection *= 1.0/TMath::Sqrt(1.0+tt*tt+tp*tp);
627 } else { // Just use fixed beam position and angle
631 fDirection.SetXYZ(tt, tp ,1.0); // Set arbitrarily to avoid run time warnings
632 fDirection *= 1.0/TMath::Sqrt(1.0+tt*tt+tp*tp);
633 }
634 /*
635 fXA_pos = fPosition[1](0);
636 fYA_pos = fPosition[1](1);
637 fXB_pos = fPosition[2](0);
638 fYB_pos = fPosition[2](1);
639 */
640 //std::cout<<" Setting fPosition and fDirection ... " << std::endl;
641 //std::cout<< fPosition[0](0) << " " << fPosition[0](1) << " " << fPosition[0](2) << std::endl;
642 //std::cout<< fPosition[1](0) << " " << fPosition[1](1) << " " << fPosition[1](2) << std::endl;
643 //std::cout<< fPosition[2](0) << " " << fPosition[2](1) << " " << fPosition[2](2) << std::endl;
644 //std::cout<< fDirection(0) << " " << fDirection(1) << " " << fDirection(2) << std::endl;
645
646 //std::cout<<" FRXA = "<<fXA_pos<<" FRYA = "<<fYA_pos<<std::endl;
647 //std::cout<<" FRXB = "<<fXB_pos<<" FRYB = "<<fYB_pos<<std::endl;
648
649 return 0;
650}
651
652
653
int Int_t
unsigned int UInt_t
uint32_t NSA
bool Bool_t
const Bool_t kFALSE
double Double_t
const Bool_t kTRUE
const char Option_t
Option_t Option_t TPoint TPoint const char mode
R__EXTERN class THaCutList * gHaCuts
ClassImp(VDC::AnalyticTTDConv) using namespace std
R__EXTERN class THcParmList * gHcParms
Definition THcGlobals.h:11
R__EXTERN class THcDetectorMap * gHcDetectorMap
Definition THcGlobals.h:12
void Clear(Option_t *option="") override
TObject * ConstructedAt(Int_t idx)
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
static THaAnalyzer * GetInstance()
THaEpicsEvtHandler * GetEpicsEvtHandler() const
virtual Int_t Result(const char *cutname="", EWarnMode mode=kWarn)
UInt_t GetTotNumChan() const
THaDetMap * fDetMap
THaApparatus * GetApparatus() const
TString GetString(const char *tag, UInt_t event=0) const
Bool_t IsLoaded(const char *tag) const
Hall C analyzer class.
Definition THcAnalyzer.h:12
virtual Int_t FillMap(THaDetMap *detmap, const char *detectorname)
virtual Int_t DecodeToHitList(const THaEvData &evdata, Bool_t suppress=kFALSE)
Populate the hitlist from the raw event data.
TClonesArray * fRawHitList
Definition THcHitList.h:51
void InitHitList(THaDetMap *detmap, const char *hitclass, Int_t maxhits, Int_t tdcref_cut=0, Int_t adcref_cut=0)
Save the electronics module to detector mapping and initialize a hit array of hits of class hitclass.
Int_t LoadParmValues(const DBRequest *list, const char *prefix="")
Retrieve parameter values from the parameter cache.
Class representing a single raw hit for the raster.
Detector class for fast raster.
Definition THcRaster.h:23
Double_t BPMYB_raw
Definition THcRaster.h:93
Double_t BPMXC_raw
Definition THcRaster.h:94
Double_t BPMYC_pos
Definition THcRaster.h:101
Double_t fYA_ADC
Definition THcRaster.h:103
Double_t fFrYB_ADC_zero_offset
Definition THcRaster.h:134
Double_t fYbpm_B
Definition THcRaster.h:117
Double_t fXpbpm_tar
Definition THcRaster.h:112
Double_t fYpbpm_tar
Definition THcRaster.h:113
void DeleteArrays()
Definition THcRaster.h:42
Double_t fgbpmyb_off
Definition THcRaster.h:78
Double_t fYbpm_C
Definition THcRaster.h:119
Double_t fgbpmyb_slope
Definition THcRaster.h:77
Double_t fPedADC[4]
Definition THcRaster.h:137
Double_t fFrXB_ADC_zero_offset
Definition THcRaster.h:133
Double_t FRXB_rawadc
Definition THcRaster.h:88
Double_t BPMYA_raw
Definition THcRaster.h:91
Double_t BPMXB_raw
Definition THcRaster.h:92
Double_t fFrYA_ADCperCM
Definition THcRaster.h:152
Double_t fgbeam_xpoff
Definition THcRaster.h:66
Double_t fEbeamEpics_prev
Definition THcRaster.h:128
Double_t fgbeam_ypoff
Definition THcRaster.h:68
Double_t fYbpm_A
Definition THcRaster.h:115
Double_t fgfry_dist
Definition THcRaster.h:64
Double_t fXB_pos
Definition THcRaster.h:108
Double_t BPMXC_pos
Definition THcRaster.h:100
Double_t fgbpmya_off
Definition THcRaster.h:76
Double_t fgbpmxa_slope
Definition THcRaster.h:69
Double_t fFrYA_ADC_zero_offset
Definition THcRaster.h:132
THaEpicsEvtHandler * fEpicsHandler
Definition THcRaster.h:145
Double_t BPMXB_pos
Definition THcRaster.h:98
Double_t fEbeamEpics
Definition THcRaster.h:126
Double_t fgbpmyc_slope
Definition THcRaster.h:79
TVector3 fDirection
Definition THcRaster.h:142
Int_t Decode(const THaEvData &)
Double_t fXbpm_tar
Definition THcRaster.h:110
Double_t fgbeam_xoff
Definition THcRaster.h:65
void AccumulatePedestals(TClonesArray *rawhits)
Double_t fXA_pos
Definition THcRaster.h:106
TClonesArray * frPosAdcPulseIntRaw
Definition THcRaster.h:144
Double_t fXbpm_A
Definition THcRaster.h:114
Double_t fgbpmyc_off
Definition THcRaster.h:80
Int_t fNPedestalEvents
Definition THcRaster.h:149
void Clear(Option_t *opt="")
Int_t ReadDatabase(const TDatime &date)
Double_t fXbeam_prev[4]
Definition THcRaster.h:120
Double_t fEbeamEpics_read
Definition THcRaster.h:127
Double_t fgbpmxa_off
Definition THcRaster.h:70
Int_t Process()
Double_t fFrXA_ADC_zero_offset
Definition THcRaster.h:131
Double_t FRYB_rawadc
Definition THcRaster.h:89
void CalculatePedestals()
Bool_t fFlag_use_EPICS_bpm
Definition THcRaster.h:124
Double_t FRXA_rawadc
Definition THcRaster.h:86
Bool_t fAnalyzePedestals
Definition THcRaster.h:148
Double_t fYpbeam_prev
Definition THcRaster.h:123
Double_t fgfrx_dist
Definition THcRaster.h:63
Double_t fgbeam_yoff
Definition THcRaster.h:67
void InitArrays()
Definition THcRaster.h:41
Double_t fYA_pos
Definition THcRaster.h:107
Double_t fgbpmxb_off
Definition THcRaster.h:72
Double_t FRYA_rawadc
Definition THcRaster.h:87
Double_t fXA_ADC
Definition THcRaster.h:102
Double_t fgbpmya_slope
Definition THcRaster.h:75
Double_t fFrYB_ADCperCM
Definition THcRaster.h:154
Double_t fYbeam_prev[4]
Definition THcRaster.h:121
Double_t BPMYB_pos
Definition THcRaster.h:99
Double_t BPMYA_pos
Definition THcRaster.h:97
Double_t fgpbeam
Definition THcRaster.h:62
Double_t fYB_pos
Definition THcRaster.h:109
Double_t fXbpm_B
Definition THcRaster.h:116
Double_t fXbpm_C
Definition THcRaster.h:118
Double_t fgbpmxb_slope
Definition THcRaster.h:71
Int_t DefineVariables(EMode mode)
Double_t fXB_ADC
Definition THcRaster.h:104
Double_t fFrCalMom
Definition THcRaster.h:150
Double_t fYbpm_tar
Definition THcRaster.h:111
Double_t BPMXA_raw
Definition THcRaster.h:90
Double_t fgbpmc_zpos
Definition THcRaster.h:83
Double_t fgbpmxc_off
Definition THcRaster.h:74
Double_t BPMXA_pos
Definition THcRaster.h:96
Double_t fgbpmxc_slope
Definition THcRaster.h:73
Double_t fFrXA_ADCperCM
Definition THcRaster.h:151
Double_t fXpbeam_prev
Definition THcRaster.h:122
Int_t fNhits
Definition THcRaster.h:58
TVector3 fPosition[3]
Definition THcRaster.h:141
Double_t BPMYC_raw
Definition THcRaster.h:95
Double_t fgbpmb_zpos
Definition THcRaster.h:82
Int_t fgusefr
Definition THcRaster.h:84
Double_t fFrXB_ADCperCM
Definition THcRaster.h:153
Double_t fgbpma_zpos
Definition THcRaster.h:81
Double_t fYB_ADC
Definition THcRaster.h:105
Class representing a single raw ADC hit.
Definition THcRawAdcHit.h:7
UInt_t GetNPulses() const
Gets number of set pulses.
UInt_t GetNSamples() const
Gets number of set samples.
Int_t GetPulseIntRaw(UInt_t iPulse=0) const
Gets raw pulse integral. In channels.
Int_t GetF250_NSA() const
Int_t GetIntegral(UInt_t iSampleLow, UInt_t iSampleHigh) const
Gets integral of raw samples. In channels.
Int_t fCounter
Definition THcRawHit.h:53
THcRawAdcHit & GetRawAdcHitPos()
const char * GetName() const override
Int_t GetEntries() const override
TObject * At(Int_t idx) const override
Int_t GetLast() const override
virtual void Error(const char *method, const char *msgfmt,...) const
void SetXYZ(Double_t x, Double_t y, Double_t z)
Double_t Sqrt(Double_t x)
STL namespace.
auto * tt