Mantid
Loading...
Searching...
No Matches
LoadDNSEvent.cpp
Go to the documentation of this file.
1// Mantid Repository : https://github.com/mantidproject/mantid
2//
3// Copyright © 2022 ISIS Rutherford Appleton Laboratory UKRI,
4// NScD Oak Ridge National Laboratory, European Spallation Source,
5// Institut Laue - Langevin & CSNS, Institute of High Energy Physics, CAS
6// SPDX - License - Identifier: GPL - 3.0 +
7
12#include "MantidAPI/Run.h"
23#include <stdexcept>
24#include <vector>
25
26#include <chrono>
27#include <iostream>
28
29using separator_t = std::array<uint8_t, 8>;
30static constexpr separator_t header_sep{0x00, 0x00, 0x55, 0x55, 0xAA, 0xAA, 0xFF, 0xFF};
31static constexpr separator_t block_sep = {0x00, 0x00, 0xFF, 0xFF, 0x55, 0x55, 0xAA, 0xAA}; // 0xAAAA5555FFFF0000; //
32static constexpr separator_t closing_sig = {0xFF, 0xFF, 0xAA, 0xAA, 0x55, 0x55, 0x00, 0x00}; // 0x00005555AAAAFFFF; //
33
34using namespace Mantid::DataObjects;
35using namespace Mantid::Kernel;
36using namespace Mantid::API;
37
38namespace {
39
40template <typename V1, typename V2> bool startsWith(const V1 &sequence, const V2 &subSequence) {
41 return std::equal(std::begin(subSequence), std::end(subSequence), std::begin(sequence));
42}
43
44} // namespace
45
46namespace Mantid::DataHandling {
47
49
50
56int LoadDNSEvent::confidence(Kernel::FileDescriptor &descriptor) const {
57 const std::string &extn = descriptor.extension();
58 if (extn != ".mdat")
59 return 0;
60 auto &file = descriptor.data();
61
62 std::string fileline;
63 std::getline(file, fileline);
64 if (fileline.find("mesytec psd listmode data") != std::string::npos) {
65 return 80;
66 } else
67 return 0;
68}
69
70const unsigned int MAX_BUFFER_BYTES_SIZE = 1500; // maximum buffer size in data file
71const unsigned int PIXEL_PER_TUBE = 1024; // maximum buffer size in data file
72
75
76 const std::vector<std::string> exts{".mdat"};
77 declareProperty(std::make_unique<FileProperty>("InputFile", "", FileProperty::Load, exts),
78 "A DNS mesydaq listmode event datafile.");
79
80 declareProperty<uint32_t>("ChopperChannel", 2u, std::make_shared<BoundedValidator<uint32_t>>(1, 4),
81 "The Chopper Channel (1 to 4)", Kernel::Direction::Input);
82
83 declareProperty<uint32_t>("NumberOfTubes", 128, std::make_shared<BoundedValidator<uint32_t>>(1, 128),
84 "The number of tubes, each tube has 1024 pixels (1 to 128)", Kernel::Direction::Input);
85
86 declareProperty<bool>("DiscardPreChopperEvents", true, std::make_shared<BoundedValidator<bool>>(0, 1),
87 "Discards events before first chopper trigger (turn off for elastic)",
89
90 declareProperty<bool>("SetBinBoundary", true, std::make_shared<BoundedValidator<bool>>(0, 1),
91 "Sets all bin boundaries to include all events (can be turned off to save time).",
93
95 std::make_unique<WorkspaceProperty<DataObjects::EventWorkspace>>("OutputWorkspace", "", Direction::Output),
96 "The name of the output workspace.");
97}
98
101 // loadProperties:
102 const std::string fileName = getPropertyValue("InputFile");
103 m_chopperChannel = static_cast<uint32_t>(getProperty("ChopperChannel"));
104 m_discardPreChopperEvents = getProperty("DiscardPreChopperEvents");
105 m_setBinBoundary = getProperty("SetBinBoundary");
106 m_detectorPixelCount = static_cast<uint32_t>(getProperty("NumberOfTubes")) * PIXEL_PER_TUBE;
107
108 // create workspace
109 EventWorkspace_sptr outputWS = std::dynamic_pointer_cast<EventWorkspace>(
110 WorkspaceFactory::Instance().create("EventWorkspace", m_detectorPixelCount, 1, 1));
111 outputWS->switchEventType(Mantid::API::EventType::TOF);
112 outputWS->getAxis(0)->setUnit("TOF");
113 outputWS->setYUnit("Counts");
114
115 // g_log.notice() << "ChopperChannel: " << m_chopperChannel << std::endl;
116
117 FileByteStream file(static_cast<std::string>(fileName));
118
119 auto finalEventAccumulator = parse_File(file, fileName);
120 populate_EventWorkspace(outputWS, finalEventAccumulator);
121 if (m_setBinBoundary && (outputWS->getNumberEvents() > 0)) {
122 outputWS->setAllX({0, outputWS->getEventXMax()});
123 }
124 setProperty("OutputWorkspace", outputWS);
125}
126
128 static const unsigned EVENTS_PER_PROGRESS = 100;
129 // The number of steps depends on the type of input file
130 Progress progress(this, 0.0, 1.0, finalEventAccumulator.neutronEvents.size() / EVENTS_PER_PROGRESS);
131
132 // Sort reversed (latest event first, most early event last):
133 std::sort(finalEventAccumulator.triggerEvents.begin(), finalEventAccumulator.triggerEvents.end(),
134 [](auto l, auto r) { return l.timestamp > r.timestamp; });
135
136 std::atomic<unsigned int> oversizedChanelIndexCounterA(0);
137 std::atomic<unsigned int> oversizedPosCounterA(0);
138
140 for (int eventIndex = 0; eventIndex < static_cast<int>(finalEventAccumulator.neutronEvents.size()); ++eventIndex) {
141 // uint64_t chopperTimestamp = 0;
142 unsigned int oversizedChanelIndexCounter = 0;
143 unsigned int oversizedPosCounter = 0;
144 const auto wsIndex = eventIndex;
145 auto &eventList = finalEventAccumulator.neutronEvents[eventIndex];
146 if (eventList.size() != 0) {
147 std::sort(eventList.begin(), eventList.end(), [](auto l, auto r) { return l.timestamp < r.timestamp; });
148 }
149
150 auto chopperIt = finalEventAccumulator.triggerEvents.cbegin();
151
152 auto &spectrum = eventWS->getSpectrum(wsIndex);
154
155 uint64_t numProcessed = 0;
156 for (const auto &event : eventList) {
157 numProcessed++;
158 if (numProcessed % EVENTS_PER_PROGRESS == 0) {
159 progress.report();
160 if (this->getCancel()) {
161 throw CancelException();
162 }
163 }
164
165 chopperIt =
166 std::lower_bound(finalEventAccumulator.triggerEvents.cbegin(), finalEventAccumulator.triggerEvents.cend(),
167 event.timestamp, [](auto l, auto r) { return l.timestamp > r; });
168 const uint64_t chopperTimestamp = chopperIt != finalEventAccumulator.triggerEvents.cend()
169 ? chopperIt->timestamp
170 : 0; // before first chopper trigger
171 if ((chopperTimestamp == 0) && m_discardPreChopperEvents) {
172 // throw away events before first chopper trigger
173 continue;
174 }
175 spectrum.addEventQuickly(Types::Event::TofEvent(double(event.timestamp - chopperTimestamp) / 10.0));
176 }
177
179 oversizedChanelIndexCounterA += oversizedChanelIndexCounter;
180 oversizedPosCounterA += oversizedPosCounter;
181 }
183
184 if (oversizedChanelIndexCounterA > 0) {
185 g_log.warning() << "Bad chanel indices: " << oversizedChanelIndexCounterA << std::endl;
186 }
187 if (oversizedPosCounterA > 0) {
188 g_log.warning() << "Bad position values: " << oversizedPosCounterA << std::endl;
189 }
190}
191
192std::vector<uint8_t> LoadDNSEvent::parse_Header(FileByteStream &file) {
193 // using Boyer-Moore String Search:
194
195 // search for header_sep and store actual header:
196 std::vector<uint8_t> header;
197
198 std::array<uint8_t, header_sep.size()> current_window;
199 file.readRaw(current_window);
200 while (!file.eof()) {
201 if (current_window == header_sep) {
202 return header;
203 } else {
204 const auto orig_header_size = header.size();
205 header.resize(header.size() + 1);
206 const auto win_data = current_window.data();
207 std::copy(win_data, win_data + 1, header.data() + orig_header_size);
208 const std::array<uint8_t, header_sep.size()> orig_window = current_window;
209 file.readRaw(current_window, 1);
210 std::copy(orig_window.data() + 1, orig_window.data() + header_sep.size(), win_data);
211 }
212 }
213 return header;
214}
215
216std::vector<std::vector<uint8_t>> LoadDNSEvent::split_File(FileByteStream &file, const unsigned maxChunckCount) {
217
218 const uint64_t minChunckSize = MAX_BUFFER_BYTES_SIZE;
219 const uint64_t chunckSize = std::max(minChunckSize, file.fileSize() / maxChunckCount);
220
221 std::vector<std::vector<uint8_t>> result;
222
223 while (!file.eof()) {
224 result.push_back(std::vector<uint8_t>());
225 // read a big chunck of file:
226 auto &data = result.back();
227 data.resize(chunckSize);
228 try {
229 file.readRaw(*data.begin(), chunckSize);
230 } catch (std::ifstream::failure &) {
231 data.resize(file.gcount());
232 return result;
233 }
234
235 // search for a block_separator, and append everything up to it :
236 static const auto windowSize = block_sep.size();
237 uint8_t *current_window = nullptr;
238 std::array<uint8_t, windowSize> *windowAsArray =
239 reinterpret_cast<std::array<uint8_t, windowSize> *>(current_window);
240
241 try {
242 data.resize(data.size() + windowSize);
243 // accomodate for possible relocation of vector...:
244 current_window = &*(data.end() - windowSize);
245 windowAsArray = reinterpret_cast<std::array<uint8_t, windowSize> *>(current_window);
246 file.readRaw(current_window[0], windowSize);
247
248 while ((*windowAsArray != block_sep) && (!file.eof())) {
249 const auto orig_data_size = data.size();
250 data.resize(orig_data_size + 1);
251 // accomodate for possible relocation of vector...:
252 current_window = (&data.back() - windowSize + 1);
253 windowAsArray = reinterpret_cast<std::array<uint8_t, windowSize> *>(current_window);
254 file.readRaw(current_window[windowSize - 1], 1);
255 }
256 } catch (std::ifstream::failure &) {
257 return result;
258 }
259 } // while
260 return result;
261}
262
264 // File := Header Body
265 std::vector<uint8_t> header = parse_Header(file);
266
267 // check it is actually a mesytec psd listmode file:
268 if (!startsWith(header, std::string("mesytec psd listmode data"))) {
269 g_log.error() << "This seems not to be a mesytec psd listmode file: " << fileName;
270 throw Exception::FileError("This seems not to be a mesytec psd listmode file: ", fileName);
271 }
272
273 // Parse actual data:
274 const int threadCount = PARALLEL_GET_MAX_THREADS;
275 // Split File:
276 std::vector<std::vector<uint8_t>> filechuncks = split_File(file, threadCount);
277 g_log.debug() << "filechuncks count = " << filechuncks.size() << std::endl;
278
279 std::vector<EventAccumulator> eventAccumulators(filechuncks.size());
280 for (auto &evtAcc : eventAccumulators) {
281 evtAcc.neutronEvents.resize(m_detectorPixelCount);
282 }
283
284 // parse individual file chuncks:
286 for (int i = 0; i < static_cast<int>(filechuncks.size()); ++i) {
287 auto filechunck = filechuncks[static_cast<size_t>(i)];
288 g_log.debug() << "filechunck.size() = " << filechunck.size() << std::endl;
289 auto vbs = VectorByteStream(filechunck);
290 parse_BlockList(vbs, eventAccumulators[static_cast<size_t>(i)]);
291 }
292
293 EventAccumulator finalEventAccumulator;
294 finalEventAccumulator.neutronEvents.resize(m_detectorPixelCount);
295
296 // combine eventAccumulators:
297
298 // combine triggerEvents:
299 for (const auto &v : eventAccumulators) {
300 finalEventAccumulator.triggerEvents.insert(finalEventAccumulator.triggerEvents.end(), v.triggerEvents.begin(),
301 v.triggerEvents.end());
302 }
303
304 // combine neutronEvents:
306 for (int i = 0; i < static_cast<int>(finalEventAccumulator.neutronEvents.size()); ++i) {
307 auto &allNeutronEvents = finalEventAccumulator.neutronEvents[static_cast<size_t>(i)];
308 for (const auto &v : eventAccumulators) {
309 allNeutronEvents.insert(allNeutronEvents.end(), v.neutronEvents[static_cast<size_t>(i)].begin(),
310 v.neutronEvents[static_cast<size_t>(i)].end());
311 }
312 }
313
314 return finalEventAccumulator;
315}
316
318 // BlockList := DataBuffer BlockListTrail
319 while (!file.eof() && file.peek() != 0xFF) {
320 parse_Block(file, eventAccumulator);
321 }
322}
323
325 // Block := DataBufferHeader DataBuffer
326 parse_DataBuffer(file, eventAccumulator);
328}
329
331 auto separator = file.readRaw(separator_t());
332 if (separator != block_sep) {
333 throw std::runtime_error(std::string("File Integrety LOST. 0x"));
334 }
335}
336
338 const auto bufferHeader = parse_DataBufferHeader(file);
339
340 const uint16_t dataLength = uint16_t(bufferHeader.bufferLength - 21);
341 const auto event_count = dataLength / 3;
342
343 for (uint16_t i = 0; i < event_count; i++) {
344 parse_andAddEvent(file, bufferHeader, eventAccumulator);
345 }
346}
347
349 uint16_t ts1 = 0;
350 uint16_t ts2 = 0;
351 uint16_t ts3 = 0;
352 BufferHeader header = {};
353 file.read<2>(header.bufferLength);
354 file.read<2>(header.bufferVersion);
355 file.read<2>(header.headerLength);
356 file.read<2>(header.bufferNumber);
357 file.read<2>(header.runId);
358 file.read<1>(header.mcpdId);
359 file.read<1>(header.deviceStatus);
360 file.read<2>(ts1);
361 file.read<2>(ts2);
362 file.read<2>(ts3);
363 // 48 bit timestamp is 3 2-byte MSB words ordered LSB
364 header.timestamp = (uint64_t)ts3 << 32 | (uint64_t)ts2 << 16 | ts1;
365 file.skip<24>();
366 return header;
367}
368
370 const LoadDNSEvent::BufferHeader &bufferHeader) {
371 TriggerEvent triggerEvent = {};
372 uint8_t trigId = 0;
373 uint8_t dataId = 0;
374 trigId = (data >> 44) & 0b111; // 3 bit
375 dataId = (data >> 40) & 0b1111; // 4 bit
376 triggerEvent.event.timestamp = data & 0x7ffff; // 19bit
377 triggerEvent.event.timestamp += bufferHeader.timestamp;
378 triggerEvent.isChopperTrigger = ((dataId == m_chopperChannel - 1) && (trigId == 7));
379 return triggerEvent;
380}
381
383 const LoadDNSEvent::BufferHeader &bufferHeader) {
384 NeutronEvent neutronEvent = {};
385 uint16_t position = 0;
386 uint8_t modid = 0;
387 uint8_t slotid = 0;
388 modid = (data >> 44) & 0b111; // 3bit
389 slotid = (data >> 39) & 0b11111; // 5bit
390 position = (data >> 19) & 0x3ff; // 10bit
391 neutronEvent.event.timestamp = (data & 0x7ffff); // 19bit
392 neutronEvent.event.timestamp += bufferHeader.timestamp;
393 neutronEvent.wsIndex = (bufferHeader.mcpdId << 6u | modid << 3u | slotid) << 10u | position;
394 return neutronEvent;
395}
396
398 LoadDNSEvent::EventAccumulator &eventAccumulator) {
399 uint16_t data1 = 0;
400 uint16_t data2 = 0;
401 uint16_t data3 = 0;
402 uint64_t data = 0;
403 event_id_e eventId;
404 file.read<2>(data1);
405 file.read<2>(data2);
406 file.read<2>(data3);
407 // 48 bit event is 3 2-byte MSB words ordered LSB
408 data = (uint64_t)data3 << 32 | (uint64_t)data2 << 16 | data1;
409 eventId = static_cast<event_id_e>((data >> 47) & 0b1);
410 switch (eventId) {
411 case event_id_e::TRIGGER: {
412 const auto triggerEvent = processTrigger(data, bufferHeader);
413 if (triggerEvent.isChopperTrigger) {
414 eventAccumulator.triggerEvents.push_back(triggerEvent.event);
415 }
416 } break;
417 case event_id_e::NEUTRON: {
418 const auto neutronEvent = processNeutron(data, bufferHeader);
419 eventAccumulator.neutronEvents[neutronEvent.wsIndex].push_back(neutronEvent.event);
420 } break;
421 default:
422 g_log.error() << "unknown event id " << eventId << "\n";
423 break;
424 }
425}
426
428 auto separator = file.readRaw(separator_t());
429 if (separator != closing_sig) {
430 throw std::runtime_error(std::string("File Integrety LOST. 0x"));
431 }
432}
433
434} // namespace Mantid::DataHandling
double position
Definition GetAllEi.cpp:154
static constexpr separator_t closing_sig
static constexpr separator_t block_sep
std::array< uint8_t, 8 > separator_t
static constexpr separator_t header_sep
#define PARALLEL_START_INTERRUPT_REGION
Begins a block to skip processing is the algorithm has been interupted Note the end of the block if n...
#define PARALLEL_FOR_NO_WSP_CHECK()
#define PARALLEL_END_INTERRUPT_REGION
Ends a block to skip processing is the algorithm has been interupted Note the start of the block if n...
#define PARALLEL_FOR_IF(condition)
Empty definitions - to enable set your complier to enable openMP.
#define PARALLEL_GET_MAX_THREADS
#define PARALLEL_CHECK_INTERRUPT_REGION
Adds a check after a Parallel region to see if it was interupted.
#define DECLARE_FILELOADER_ALGORITHM(classname)
DECLARE_FILELOADER_ALGORITHM should be used in place of the standard DECLARE_ALGORITHM macro when wri...
bool eof() const
Definition BitStream.h:61
FileByteStream & readRaw(T &result, const std::size_t &bytecount)
Definition BitStream.h:41
uint64_t fileSize() const
Definition BitStream.h:71
std::streamsize gcount() const
Definition BitStream.h:68
CancelException is thrown to cancel execution of the algorithm.
Definition Algorithm.h:136
void declareProperty(std::unique_ptr< Kernel::Property > p, const std::string &doc="") override
Add a property to the list of managed properties.
std::string getPropertyValue(const std::string &name) const override
Get the value of a property as a string.
TypedValue getProperty(const std::string &name) const override
Get the value of a property.
Kernel::Logger & g_log
Definition Algorithm.h:422
void progress(double p, const std::string &msg="", double estimatedTime=0.0, int progressPrecision=0)
Sends ProgressNotification.
bool getCancel() const
Returns the cancellation state.
@ Load
allowed here which will be passed to the algorithm
Helper class for reporting progress from algorithms.
Definition Progress.h:25
A property class for workspaces.
TriggerEvent processTrigger(const uint64_t &data, const LoadDNSEvent::BufferHeader &bufferHeader)
std::vector< uint8_t > parse_Header(FileByteStream &file)
void parse_BlockList(VectorByteStream &file, EventAccumulator &eventAccumulator)
void populate_EventWorkspace(Mantid::DataObjects::EventWorkspace_sptr &eventWS, EventAccumulator &finalEventAccumulator)
void parse_andAddEvent(VectorByteStream &file, const BufferHeader &bufferHeader, EventAccumulator &eventAccumulator)
NeutronEvent processNeutron(const uint64_t &data, const LoadDNSEvent::BufferHeader &bufferHeader)
void exec() override
Run the algorithm.
EventAccumulator parse_File(FileByteStream &file, const std::string &fileName)
void parse_BlockSeparator(VectorByteStream &file)
BufferHeader parse_DataBufferHeader(VectorByteStream &file)
void parse_EndSignature(FileByteStream &file)
void parse_DataBuffer(VectorByteStream &file, EventAccumulator &eventAccumulator)
std::vector< std::vector< uint8_t > > split_File(FileByteStream &file, const unsigned maxChunckCount)
void init() override
Virtual method - must be overridden by concrete algorithm.
void parse_Block(VectorByteStream &file, EventAccumulator &eventAccumulator)
BoundedValidator is a validator that requires the values to be between upper or lower bounds,...
Records the filename and the description of failure.
Definition Exception.h:98
Defines a wrapper around an open file.
IPropertyManager * setProperty(const std::string &name, const T &value)
Templated method to set the value of a PropertyWithValue.
void debug(const std::string &msg)
Logs at debug level.
Definition Logger.cpp:145
void error(const std::string &msg)
Logs at error level.
Definition Logger.cpp:108
void warning(const std::string &msg)
Logs at warning level.
Definition Logger.cpp:117
static T & Instance()
Return a reference to the Singleton instance, creating it if it does not already exist Creation is do...
VectorByteStream & read(T &result)
Definition BitStream.h:117
uint8_t peek()
Definition BitStream.h:123
VectorByteStream & skip()
Definition BitStream.h:127
VectorByteStream & readRaw(T &result, const std::size_t &bytecount)
Definition BitStream.h:98
bool eof() const
Definition BitStream.h:125
bool startsWith(const string &str, const string &prefix)
Returns true if str starts with prefix.
const unsigned int PIXEL_PER_TUBE
const unsigned int MAX_BUFFER_BYTES_SIZE
std::shared_ptr< EventWorkspace > EventWorkspace_sptr
shared pointer to the EventWorkspace class
std::unique_ptr< T > create(const P &parent, const IndexArg &indexArg, const HistArg &histArg)
This is the create() method that all the other create() methods call.
std::enable_if< std::is_pointer< Arg >::value, bool >::type threadSafe(Arg workspace)
Thread-safety check Checks the workspace to ensure it is suitable for multithreaded access.
std::vector< std::vector< CompactEvent > > neutronEvents
@ Input
An input workspace.
Definition Property.h:53
@ Output
An output workspace.
Definition Property.h:54