SND@LHC Software
Loading...
Searching...
No Matches
ScalerUnpack.cxx
Go to the documentation of this file.
1#include <cassert>
2
3// ROOT headers
4#include "TTree.h"
5#include "TFile.h"
6#include "ROOT/TSeq.hxx"
7
8// Fair headers
9#include "FairLogger.h"
10#include "FairRootManager.h"
11
12// SHiP headers
13#include "ScalerUnpack.h"
15
16// ScalerUnpack: Constructor
18
19// Virtual ScalerUnpack: Public method
21
22// Init: Public method
24{
25 Register();
26 fMan = FairRootManager::Instance();
27 return kTRUE;
28}
29
30// Register: Protected method
32{
33 LOG(DEBUG) << "ScalerUnpack : Registering..." ;
34}
35
38 uint16_t PSW;
39 uint16_t SPW;
40 uint32_t PinStatus;
41 uint32_t scalars[16];
42 uint32_t slices[0];
43 int getSliceCount() { return (header.size - 88) / sizeof(uint32_t); }
44};
45
46// DoUnpack: Public method
47Bool_t ScalerUnpack::DoUnpack(Int_t *data, Int_t size)
48{
49 static_assert(sizeof(ScalarFrame) == 88, "Scaler frame size incorrect!");
50 auto df = reinterpret_cast<ScalarFrame *>(data);
51 LOG(DEBUG) << "ScalerUnpack : Unpacking frame... size/bytes = " << size ;
52 LOG(DEBUG) << "PSW = " << df->PSW ;
53 LOG(DEBUG) << "SPW = " << df->SPW ;
54 LOG(DEBUG) << "POT from SPS = " << df->scalars[0] ;
55 LOG(DEBUG) << "S1raw = " << df->scalars[1] ;
56 LOG(DEBUG) << "S1strobed = " << df->scalars[2] ;
57 LOG(DEBUG) << "S1*S2 TrgRaw = " << df->scalars[3] ;
58 LOG(DEBUG) << "S1*S2 TrgStrobed = " << df->scalars[4] ;
59 for (auto i : ROOT::MakeSeq(df->getSliceCount())) {
60 LOG(DEBUG) << "Slice " << i << "= " << df->slices[i] ;
61 }
62 auto f = fMan->GetOutFile();
63 tree = dynamic_cast<TTree *>(f->Get("scalers"));
64 if (tree == nullptr) {
65 tree = new TTree("scalers", "scalers");
66 }
67 tree->Branch("PSW", &(df->PSW));
68 tree->Branch("SPW", &(df->SPW));
69 for (auto i : ROOT::MakeSeq(16)) {
70 switch (i) {
71 case 11: {
72 fGoliath = int(int(df->scalars[i]) / 0x10000);
73 fDavid = int(int(df->scalars[i]) % 0x10000);
74 tree->Branch("Goliath", &fGoliath);
75 tree->Branch("David", &fDavid);
76 LOG(INFO) << "David: " << fDavid ;
77 LOG(INFO) << "Goliath: " << fGoliath ;
78 break;
79 }
80 case 12:
81 LOG(INFO) << "Spill type: " << df->scalars[i] ;
82 tree->Branch("spill_type", &(df->scalars[i]));
83 break;
84 default: tree->Branch(TString::Format("SC%.2d", i), &(df->scalars[i]));
85 }
86 }
87 auto slices = df->getSliceCount() > 0 ? std::vector<uint32_t>(df->slices, df->slices + df->getSliceCount())
88 : std::vector<uint32_t>();
89 tree->Branch("slices", "vector<uint32_t>", &slices);
90 tree->Fill();
91 tree->Write("", TObject::kOverwrite);
92 return kTRUE;
93}
94
95// Reset: Public method
97{
98 LOG(DEBUG) << "ScalerUnpack : Clearing Data Structure" ;
99}
100
virtual void Register() override
virtual Bool_t DoUnpack(Int_t *data, Int_t size) override
FairRootManager * fMan
virtual ~ScalerUnpack()
virtual Bool_t Init() override
virtual void Reset() override
ClassImp(ecalContFact) ecalContFact
int getSliceCount()
uint32_t scalars[16]
uint32_t PinStatus
DataFrameHeader header
uint32_t slices[0]