KaliVeda
Toolkit for HIC analysis
Loading...
Searching...
No Matches
GTGanilData.cpp
1/***************************************************************************
2* GTGanilData.cpp - Main Header to ROOTGAnilTape
3* -------------------
4* begin : Thu Jun 14 2001
5* copyright : (C) 2001 by Garp
6* email : patois@ganil.fr
7***************************************************************************
8 * *
9 * This program is free software; you can redistribute it and/or modify *
10 * it under the terms of the GNU General Public License as published by *
11 * the Free Software Foundation; either version 2 of the License, or *
12 * (at your option) any later version. *
13 * *
14 ***************************************************************************/
15
16#include "Riostream.h"
17#include "KVConfig.h"
18#ifdef CCIN2P3_XRD
19#include "TSystem.h"
20#endif
21// ganil_tape lib headers
22
24
25#include <ERR_GAN.H>
26#include <gan_tape_erreur.h>
27#include <gan_tape_general.h>
28#include <gan_tape_param.h>
29#include <gan_tape_test.h>
30#include <gan_tape_alloc.h>
31#include <gan_tape_mount.h>
32#include <gan_tape_file.h>
33#include <acq_mt_fct_ganil.h>
34#include <acq_ebyedat_get_next_event.h>
35#include <gan_tape_get_parametres.h>
36extern "C" {
37 typedef struct SCALE {
38 UNSINT32 Length; /* Nb bytes utils suivant ce mot */
39 UNSINT32 Nb_channel;
40 UNSINT32 Acq_status;
41 UNSINT32 Reserve[2];
42 union JBUS_SCALE {
43 scale_struct UnScale[NB_MAX_CHANNEL];
44 UNSINT16 Info_jbus [NB_MAX_JBUS];
46 } scale;
47}
48
49
50#include "GTGanilData.h"
51
52bool AutoswBuf, Swbufon;
53Int_t BufSize;
54
55using namespace std;
56
58
59
60
61
63
65{
66 // Default constructor
67
68 InitDefault();
69}
70
71
72
75
77{
78 // Create the reading class over the given file with all defaults values
79
82}
83
84
85
89
91{
92 //Set the name of the file to read (use if default ctor is used to create object).
93 //We check if 'filename' begins with "rfio:"; if so, we remove it
94
96 if (filename.BeginsWith("rfio:")) fFileName.Remove(0, 5);
97 strncpy(fDevice->DevName, fFileName, MAX_CARACTERES);
98}
99
100
101
103
105{
106 if (fScaler)delete fScaler;
108 if (fDevice)delete fDevice;
109 if (fBuffer)delete fBuffer;
110 if (fStructEvent)delete[] fStructEvent;
111 if (fDataArray) delete[] fDataArray;
112 if (fFired) delete[] fFired;
113 if (fEventBrut) delete[] fEventBrut;
114 if (fEventCtrl) delete[] fEventCtrl;
115}
116
117
118
119
127
128void GTGanilData::InitDefault(const Int_t argc, char** argv)
129{
130 // PRIVATE
131 // Called by every constructors.
132 // Init parameters to default values
133 // Data are taken first from command line,
134 // Second from environment variables
135 // if none works, takes defaults values.
136#ifdef CCIN2P3_XRD
137 //if we are using version of GanTape compiled with Xrootd support
138 //we need to load the ROOT XrdPosix library now
139 gSystem->Load("libXrdPosix");
140#endif
141 fStatus = ACQ_OK;
143 fEventNumber = -1;
144 fEventCount = 0;
145 fIsCtrl = false;
146 fIsScalerBuffer = false;
147 fDataArray = 0;
148 fFired = nullptr;
149
151 fBuffer = new in2p3_buffer_struct;
152 fStructEvent = new char[STRUCTEVENTSIZE];
153
154 fScalerTree = NULL; // By default, no scaler tree
155 fScaler = new GTScalers; // Scaler class
156
157 val_ret DeRetour; // Union to communicate with ganil_tape lib routines
158 // Getting Logical Unit ID number
159 fStatus = acq_get_param_env(lun_id, &DeRetour, argc, argv);
160 if (DeRetour.Val_INT >= 0) fLun = DeRetour.Val_INT;
161 else fLun = 10; // Random
162 fDevice->Lun = fLun;
163 // Getting tape/file name
164 fStatus = acq_get_param_env(device_id, &DeRetour, argc, argv);
165 strncpy(fDevice->DevName, DeRetour.Val_CAR, MAX_CARACTERES);
166 fFileName = DeRetour.Val_CAR;
167
168 // GLOBAL VARIABLES ASSIGNED HERE
169 fStatus = acq_get_param_env(autoswbuf_id, &DeRetour, argc, argv);
170 AutoswBuf = DeRetour.Val_BOL;
171
172 fStatus = acq_get_param_env(swbufon_id, &DeRetour, argc, argv);
173 Swbufon = DeRetour.Val_BOL;
174
175 // Variable/fix length events
176 fStatus = acq_get_param_env(ctrlform_id, &DeRetour, argc, argv);
177 fCtrlForm = DeRetour.Val_INT;
178
179 fStatus = acq_get_param_env(density_id, &DeRetour, argc, argv);
180 fDensity = DeRetour.Val_INT;
181
182 fStatus = acq_get_param_env(bufsize_id, &DeRetour, argc, argv);
183 fBufSize = DeRetour.Val_INT;
184 // GLOBAL VARIABLES ASSIGNED HERE
185 BufSize = fBufSize;
186
187 fStatus = acq_get_param_env(evbsize_id, &DeRetour, argc, argv);
188 fEvbsize = DeRetour.Val_INT;
189
190 fStatus = acq_get_param_env(evcsize_id, &DeRetour, argc, argv);
191 fEvcsize = DeRetour.Val_INT;
192
195
196 // calculate size of CTRL_EVNT header to be subtracted from size of event
197 // returned by acq_ebyedat_get_next_event in order to calculate number of
198 // label-value parameter pairs
199 fCTRLEVNT_HD = sizeof(CTRL_EVENT) / 2 - 1;
200}
201
202
203
206
208{
209 // Print every class parameters, for now a simple dump.
211}
212
213
214
218
220{
221 // Not implemented.
222 // ...
223}
224
225
226
233
235{
236 // Open the data file (could be on a tape) with a few checks.
237 // Use IsOpen() to check whether the file is opened successfully.
238 // After successfully opening the file, we call ReadParameters()
239 // to fill the parameter list from the first data buffer
240 // which *has to* be a parameter buffer.
241
242 Int_t structEvent_size = STRUCTEVENTSIZE; // bidon
243
244 fStatus = acq_dev_is_alloc_c(*fDevice); // Allocation test
245 if (fStatus == ACQ_OK) fAllocated = true;
246 else {
247 gan_tape_erreur(fStatus, "test d'allocation");
248 return;
249 }
250
251 fStatus = acq_mt_alloc_c(*fDevice); // Allocation
252 if (fStatus != ACQ_OK) {
253 if (fStatus == ACQ_NOTALLOC)
254 cout << "This tape drive is already in use !" << endl;
255 gan_tape_erreur(fStatus, "allocation");
256 return;
257 }
258
259 fStatus = acq_mt_mount_c(*fDevice, fDensity, fBufSize);
260 if (fStatus != ACQ_OK) { // There might be some imprecisions here
261 if (fStatus == ACQ_ALREADYMOUNT)
262 cout << "The tape you want to use is already mounted" << endl;
263 gan_tape_erreur(fStatus, "mount");
264 return;
265 }
266
267 fStatus = acq_mt_open_c(fDevice, o_read, &fLun);
268 if (fStatus != ACQ_OK) { // There might be some imprecisions here
269 gan_tape_erreur(fStatus, "open");
270 return;
271 }
272 fRunNumber = 0;
273 Char_t Date[21];
274 fStatus = acq_mt_ini_run_c(fLun, fBuffer, fBufSize, &fRunNumber,
275 fStructEvent, structEvent_size, Date);
276 fDateStart = Date;
277
278 if (fStatus != ACQ_OK) {
279 gan_tape_erreur(fStatus, "Init Run");
280 return;
281 }
282 cout << "Run " << fRunNumber << " opened" << endl;
284 return;
285}
286
287
288
289
294
296{
297 // PRIVATE
298 // Read the data parameters from the current buffer, put it in the
299 // parameter list.
300
301 do {
302 ReadBuffer();
303 if (strcmp(fHeader, PARAM_Id) == 0)
304 fDataArraySize = fDataParameters->Fill(fBuffer->les_donnees.cas.Buf_param);
305 }
306 while (strcmp(fHeader, PARAM_Id) == 0);
307
308 fDataArray = new UShort_t[fDataArraySize + 1]; // Data buffer is allocated
309 fFired = new Bool_t[fDataArraySize + 1]; // Data buffer is allocated
310 for (Int_t i = 1; i <= fDataArraySize; i++) {
311 fFired[i] = false;
312 }
313}
314
315
316
319
320void GTGanilData::Connect(const Int_t index, UShort_t** p) const
321{
322 // Connect a pointer to a data to the defined index in the Data Array
323
324 if ((index < 1) || (index > fDataArraySize)) {
325 cout << "Invalid connexion:" << index << ". Valid only in 1<=index<="
326 << fDataArraySize << endl;
327 return;
328 }
329 *p = &(fDataArray[index]);
330}
331
332
333
336
337void GTGanilData::ConnectFired(const Int_t index, Bool_t** p) const
338{
339 // Connect a pointer to a data to the defined index in the Data Array
340
341 if ((index < 1) || (index > fDataArraySize)) {
342 cout << "Invalid connexion:" << index << ". Valid only in 1<=index<="
343 << fDataArraySize << endl;
344 return;
345 }
346 *p = &(fFired[index]);
347}
348
349
350
386
388{
389 // Read an event on tape/file and put in into the event array
390 // return true until read fails
391 //
392 // If scaler buffer management (fWhatScaler) has been set to
393 // kSkipScaler, kDumpScaler or kAutoWriteScaler, then every time this method
394 // returns kTRUE a new event has been read (and perhaps 1 or more scaler buffers
395 // were read and dealt with internally).
396 // In this case, a loop over all events will look like this:
397 //
398 // while( my_gtganildata->Next() ){
399 // // new event read from file
400 // my_gtganilData->GetFiredDataParameters()->ls(); // or whatever
401 // }
402 //
403 // If fWhatScaler = kReportScaler then this method also returns kTRUE after
404 // reading a scaler buffer (no new event read), so that the user can do
405 // something with the scalers.
406 // In this case, a loop over all events including treatment of scaler buffers
407 // will look like this:
408 //
409 // while( my_gtganildata->Next() ){
410 // if( my_gtganildata->IsScalerBuffer() ){
411 // // scaler buffer read from file
412 // GTScalers* scalers = my_gtganildata->GetScalers();
413 // // N.B. GTGanilData::GetScalers() also resets the IsScalerBuffer()
414 // // flag ready for next event/buffer, so even if you don't
415 // // do anything with the scalers, you should call it
416 // }
417 // else
418 // {
419 // // new event read from file
420 // my_gtganilData->GetFiredDataParameters()->ls(); // or whatever
421 // }
422 // }
423
424 if (fEventNumber == -1) ReadBuffer();
425 if (fStatus) return (false); // Maybe more to say here
426 TString _heads(fHeader);
427 _heads.ToUpper();
428 if (_heads.Contains("ENDRUN")) {
429 //we have reached the end of the file
430 return (false);
431 }
432 if (strcmp(fHeader, SCALER_Id) == 0) {
433 switch (fWhatScaler) {
434 case kSkipScaler:
435 return (Next()); // Skip
436 case kDumpScaler:
437 fScaler->Fill((scale*) & (fBuffer->les_donnees.cas.scale));
439 return (Next());
440
441 case kReportScaler:
442 fScaler->Fill((scale*) & (fBuffer->les_donnees.cas.scale));
443 fIsScalerBuffer = true;
444 return kTRUE;
445
446 case kAutoWriteScaler: {
447 fScaler->Fill((scale*) & (fBuffer->les_donnees.cas.scale));
449 fScalerTree->Fill();
450 return (Next());
451 }
452 break;
453 default:
454 cout << "GTGanilData::Next: unexpected value for fWhatScaler" << endl;
455 }
456 }
457 if (strcmp(fHeader, EVENTDB_Id) == 0) {
458 fIsScalerBuffer = false;
459 fIsCtrl = false;
460 fEventCount++;
461 return (ReadNextEvent());
462 }
463 else if (strcmp(fHeader, EVENTCT_Id) == 0) {
464 fIsScalerBuffer = false;
465 fIsCtrl = true;
466 fEventCount++;
467 return (ReadNextEvent());
468 }
469 else if (strcmp(fHeader, EBYEDAT_Id) == 0) {
470 fIsScalerBuffer = false;
471 fIsCtrl = false;
472 fEventCount++;
473 return (ReadNextEvent_EBYEDAT());
474 }
475 cout << "Unknown header in GTGanilData::Next:" << fHeader << endl;
476 return (true);
477}
478
479
480
485
486void GTGanilData::MakeTree(const TString filename, UInt_t nEvents)
487{
488 // Automatically create and fill a tree created from ganil data.
489 // Can be use to convert a ganil tape or file to a ROOT file.
490 // The number of actual converted events is set with the nEvents parameter.
491
492 TString theFilename;
493 if (filename == "")
494 theFilename = fFileName + ".root";
495 else theFilename = filename;
496 TFile theTreeFile(theFilename, "recreate");
497 TTree theTree("AutoTree", "Automatic filled Tree");
498 for (Int_t i = 1; i <= fDataArraySize; i++) {
499 TString parName = fDataParameters->GetParName(i);
500 TString parType = parName + "/s";
501 theTree.Branch(parName, &(fDataArray[i]), parType);
502 }
503
504 while (Next()) {
505 theTree.Fill();
506 nEvents--;
507 if (nEvents <= 0) break;
508 }
509 theTree.Write();
510}
511
512
513
514
518
520{
521 // PRIVATE
522 // Utility routine to read events from buffers.
523
524 if (!fIsCtrl) {
525 fStatus = get_next_event(fBuffer, fBufSize,
526 (short*)fEventBrut, fEvbsize, &fEventNumber);
527 if (fStatus == ACQ_OK) {
528 fStatus = s_evctrl((short*)fEventBrut, (short*)fEventCtrl,
530 if (fStatus == GR_OK)
531 return (EventUnravelling((CTRL_EVENT*)fEventCtrl));
532 else {
533 puts("\n>>> Erreur de reconstruction");
534 return (false);
535 }
536 }
537 else if (fStatus == ACQ_ENDOFBUFFER) {
538 fEventNumber = -1;
539 return (Next());
540 }
541 else {
542 gan_tape_erreur(fStatus, "obtention d'evenement");
543 return (false);
544 }
545 }
546 else {
547 fStatus = get_next_event(fBuffer, fBufSize, (short*)fEventCtrl,
549 if (fStatus == ACQ_OK) {
550 return (EventUnravelling((CTRL_EVENT*)fEventCtrl));
551 }
552 else if (fStatus == ACQ_ENDOFBUFFER) {
553 fEventNumber = -1;
554 return (Next());
555 }
556 else {
557 gan_tape_erreur(fStatus, "obtention d'evenement");
558 return (false);
559 }
560 }
561}
562
563
564
569
571{
572 // PRIVATE
573 // Utility routine to read EBYEDAT events from buffers.
574 //Shamelessly copied from Luc Legeard's GEvent
575 if (!fIsCtrl) {
576 UNSINT16* ebyeEventCtrl;
577 fStatus = acq_ebyedat_get_next_event((UNSINT16*)fBuffer, &ebyeEventCtrl, &fEventNumber, EVCT_VAR);
578 if (fStatus == ACQ_OK) return (EventUnravelling((CTRL_EVENT*)ebyeEventCtrl));
579 else if (fStatus == ACQ_ENDOFBUFFER) {
580 fEventNumber = -1;
581 return (Next());
582 }
583 else {
584 cout << " in ReadNextEvent__EBEYEDAT error status: " << fStatus << " n";
585 return (false);
586 }
587 }
588 else {
589 fStatus = ACQ_BADEVENTFORM;
590 }
591 return (false);
592}
593
594
598
600{
601 // PRIVATE
602 // Read a single buffer from file
603
604 fStatus = acq_mt_read_c(fLun, fBuffer->Buffer, &fBufSize);
605 //if ( fStatus != ACQ_OK )
606 //cout << "Not OK"<<endl; // Maybe stuff to do to handle the error
607
608 strncpy(fHeader, fBuffer->les_donnees.Ident, 8);
609 fHeader[8] = '\0';
610}
611
612
613
621
622bool GTGanilData::EventUnravelling(CTRL_EVENT* pCtrlEvent)
623{
624 // PRIVATE
625 // If mode is variable length event, we have to reconstruct the Data buffer
626 // from the given event.
627 // WARNING: temporary the default: we dont check that it's really the case
628 // Before reading event, all parameters have their value set to -1 (65535 - fDataArray is UShort_t)
629 // Parameters which are not fired in the event will have value -1 (65535 - cast back to Short_t for real value)
630
631 Short_t* brutData = &(pCtrlEvent->ct_par);
632
633 Int_t eventLength = pCtrlEvent->ct_len - fCTRLEVNT_HD;
634 //printf("EventUnravelling : eventLength=%d sizeof(CTRL_EVENT)=%d\n",eventLength,sizeof(CTRL_EVENT));
635
636 for (Int_t i = 1; i <= fDataArraySize; i++) {
637 fFired[i] = false;
638 }
639
640 for (Int_t i = 0; i < eventLength; i += 2) {
641 if (brutData[i] <= fDataArraySize && brutData[i] >= 1) {
642 fDataArray[brutData[i]] = brutData[i + 1];
643 fFired[brutData[i]] = true;
644 }
645 }
646 return (true);
647}
648
649
650
651
652
653
656
658{
659 // Dump parameter index, name and value for the current event.
660
661 cout << "--------- DUMPING EVENT ----------------------------" << endl;
662 cout << "------- number:" << fEventCount << endl;
663 for (Int_t i = 1; i <= fDataArraySize; i++) {
664 if (fFired[i]) cout << "index :" << i << " " << fDataParameters->GetParName(i) << " : " << fDataArray[i] << endl;
665 }
666 cout << "--------- END ----------------------------" << endl;
667}
668
669
670
673
675{
676 // Dump parameter index and name
677
678 cout << "--------- DUMPING PARAMETERS ----------------------------" << endl;
679 cout << "------- number:" << fEventCount << endl;
680 for (Int_t i = 1; i <= fDataArraySize; i++) {
681 // WARNING : index start at one, the first value is boggus
682 cout << "index :" << i << " " << fDataParameters->GetParName(i) << endl;
683 }
684 cout << "--------- END ----------------------------" << endl;
685}
686
687
688
702
704{
705 // Set scaler buffers management. It can be:
706 // GTGanilData::kSkipScaler : Skip scaler buffers
707 // GTGanilData::kDumpScaler : Dump all scaler buffers on stdout
708 // GTGanilData::kAutoWriteScaler : Automatic scaler buffer management, all scalers written in a TTree
709 // To use this, the current TFile (i.e. gFile) must be writable.
710 // i.e. you should do TFile file("somefile.root", "create")
711 // and then toto.SetScalerBuffersManagement(GTGanilData::kReportScaler)
712 // GTGanilData::kReportScaler : when Next() encounters a scaler buffer, the IsScalerBuffer()
713 // flag is set to kTRUE and the data can be retrieved by
714 // calling GetScalers() (returns a pointer to a GScalers object).
715 // WARNING: this option changes the logic of a loop over all events
716 // in the file (see GTGanilData::Next()).
717
718 fWhatScaler = sc;
719
721 if (gFile && gFile->IsWritable()) {
722 fScalerTree = new TTree("Scalers", "Automatic filled scalers");
723 fScalerTree->Branch("scalers", fScaler, 8000, 99);
724 }
725 else {
726 cout << "Error in <GTGanilData::SetScalerBuffersManagement> : ";
727 cout << "You must open a writable TFile before calling SetScalerBuffersManagement(kAutoWriteScaler)" << endl;
729 }
730 }
731}
732
733
734
735
738
740{
741 //Returns current run number
742 return fRunNumber;
743}
744
745
746
747
749
751{
752 return (fStatus == ACQ_OK && fLun);
753}
754
755
756
757
760
762{
763 //Not used
764}
765
766
int Int_t
unsigned int UInt_t
bool Bool_t
unsigned short UShort_t
char Char_t
short Short_t
constexpr Bool_t kTRUE
#define gFile
winID h TVirtualViewer3D TVirtualGLPainter p
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char filename
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 TSystem * gSystem
const char * GetParName(const int index) const
Return the text label name corresponding to a parameter number.
int Fill(const char *buffParam)
Read GANIL formatted tapes or files.
Definition GTGanilData.h:64
Int_t fLun
Logical number of the device.
bool ReadNextEvent(void)
Int_t fRunNumber
current file run number
void ReadBuffer(void)
char fHeader[9]
Buffer header.
void ConnectFired(const Int_t index, Bool_t **p) const
Connect a pointer to a data to the defined index in the Data Array.
UShort_t * fEventBrut
Brut data buffer.
void InitDefault(const Int_t argc=0, char **argv=NULL)
bool fAllocated
True if the Tape is allocated (useless?)
Int_t fEventCount
Our event counter.
Int_t fEventNumber
Local event number in current buffer (should be renamed)
GTScalers * fScaler
Scaler array.
Int_t fCtrlForm
Fix/variable length events.
void DumpEvent(void) const
Dump parameter index, name and value for the current event.
gan_tape_desc * fDevice
Device structure (ganil_tape interface)
void Connect(const Int_t index, UShort_t **p) const
Connect a pointer to a data to the defined index in the Data Array.
Int_t fDensity
Tape density.
char * fStructEvent
??? (ganil_tape interface)
Int_t GetRunNumber(void) const
Returns current run number.
Int_t fStatus
Status, 0 is OK, any other value suspect.
void SetFileName(const TString filename)
TString fDateStart
Date/time of start of run read from file.
UShort_t * fEventCtrl
Control data buffer.
bool ReadNextEvent_EBYEDAT(void)
void PrintRunParameters(void) const
Print every class parameters, for now a simple dump.
virtual void SetUserTree(TTree *)
Not used.
bool fIsCtrl
We are currently in a CTRL buffer.
bool fIsScalerBuffer
The current buffer is a scaler buffer.
void DumpParameterName(void) const
Dump parameter index and name.
in2p3_buffer_struct * fBuffer
Brut data buffer (ganil_tape interface)
Bool_t IsOpen(void) const
Bool_t * fFired
fired parameters in each event
virtual void ReadParameters(void)
GTGanilData()
Default constructor.
GTDataParameters * fDataParameters
Data parameters names class.
TString fFileName
Filename, can be a tape drive.
Int_t fBufSize
Buffer size.
Int_t fDataArraySize
Data array size.
void SetScalerBuffersManagement(const ScalerWhat_t sc)
UShort_t * fDataArray
Physical data array.
bool Next(void)
Int_t fEvbsize
Size of the brut data buffer.
ScalerWhat_t fWhatScaler
What do we do with scalers buffers.
void MakeTree(const TString filename="", UInt_t nEvents=kMaxUInt)
virtual ~GTGanilData(void)
virtual bool EventUnravelling(CTRL_EVENT *)
Int_t fEvcsize
Size of the ctrl data buffer.
Int_t fCTRLEVNT_HD
size (in 16-bit words) of CTRL_EVNT header
void PrintDataParameters(void) const
TTree * fScalerTree
Scaler tree for automatic filling.
void Open(void)
ScalerWhat_t
What to do with scaler buffer.
Definition GTGanilData.h:67
@ kAutoWriteScaler
have to take care of it.
Definition GTGanilData.h:72
Handle scaler buffers in GANIL DAQ data.
Definition GTScalers.h:34
void Fill(void *)
Definition GTScalers.cpp:61
void DumpScalers(void)
Dump scalers value on cout.
Definition GTScalers.cpp:83
void ToUpper()
TString & Remove(EStripType s, char c)
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
virtual int Load(const char *module, const char *entry="", Bool_t system=kFALSE)
virtual Int_t Fill()
Int_t Write(const char *name=nullptr, Int_t option=0, Int_t bufsize=0) const override
virtual Int_t Branch(const char *folder, Int_t bufsize=32000, Int_t splitlevel=99)
UNSINT32 Reserve[2]
UNSINT32 Length
UNSINT32 Acq_status
union SCALE::JBUS_SCALE jbus_scale
UNSINT32 Nb_channel
ROOT headers.
Definition GTGanilData.h:42
ClassImp(TPyArg)
UNSINT16 Info_jbus[NB_MAX_JBUS]
scale_struct UnScale[NB_MAX_CHANNEL]