Mantid
Loading...
Searching...
No Matches
TimeSeriesProperty.cpp
Go to the documentation of this file.
1// Mantid Repository : https://github.com/mantidproject/mantid
2//
3// Copyright © 2018 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 +
10#include "MantidKernel/Logger.h"
13#include "MantidNexus/NexusFile.h"
14
15#include <json/value.h>
16
17#include <boost/regex.hpp>
18#include <numeric>
19
20namespace Mantid {
21using namespace Types::Core;
22namespace Kernel {
23namespace {
25Logger g_log("TimeSeriesProperty");
26
33template <typename TYPE> bool allValuesAreSame(const std::vector<TimeValueUnit<TYPE>> &values) {
34 assert(values.size() > 1);
35 const auto &first_value = values.front().value();
36 return std::all_of(std::next(values.cbegin()), values.cend(),
37 [&first_value](const auto &v) { return v.value() == first_value; });
38}
39} // namespace
40
45template <typename TYPE>
47 : Property(name, typeid(std::vector<TimeValueUnit<TYPE>>)), m_values(), m_size(), m_propSortedFlag() {}
48
55template <typename TYPE>
57 const std::vector<Types::Core::DateAndTime> &times,
58 const std::vector<TYPE> &values)
60 addValues(times, values);
61}
62
64template <typename TYPE> TimeSeriesProperty<TYPE>::~TimeSeriesProperty() = default;
65
70 return new TimeSeriesProperty<TYPE>(*this);
71}
72
77template <typename TYPE>
79 : Property(*p), m_values(), m_size(), m_propSortedFlag() {}
80
88template <typename TYPE> Property *TimeSeriesProperty<TYPE>::cloneInTimeROI(const TimeROI &timeROI) const {
89 auto filteredTS = new TimeSeriesProperty<TYPE>(this);
90
91 createFilteredData(timeROI, filteredTS->m_values);
92
93 filteredTS->m_size = static_cast<int>(filteredTS->m_values.size());
94
95 return filteredTS;
96}
97
102template <typename TYPE> Property *TimeSeriesProperty<TYPE>::cloneWithTimeShift(const double timeShift) const {
103 auto timeSeriesProperty = this->clone();
104 auto values = timeSeriesProperty->valuesAsVector();
105 auto times = timeSeriesProperty->timesAsVector();
106 // Shift the time
107 for (auto it = times.begin(); it != times.end(); ++it) {
108 // There is a known issue which can cause cloneWithTimeShift to be called
109 // with a large (~9e+9 s) shift. Actual shifting is capped to be ~4.6e+19
110 // seconds in DateAndTime::operator+=
111 (*it) += timeShift;
112 }
113 timeSeriesProperty->clear();
114 timeSeriesProperty->addValues(times, values);
115 return timeSeriesProperty;
116}
117
125template <typename TYPE> std::unique_ptr<TimeSeriesProperty<double>> TimeSeriesProperty<TYPE>::getDerivative() const {
126
127 if (this->m_values.size() < 2) {
128 throw std::runtime_error("Derivative is not defined for a time-series "
129 "property with less then two values");
130 }
131
132 this->sortIfNecessary();
133 auto it = this->m_values.begin();
134 int64_t t0 = it->time().totalNanoseconds();
135 TYPE v0 = it->value();
136
137 it++;
138 auto timeSeriesDeriv = std::make_unique<TimeSeriesProperty<double>>(this->name() + "_derivative");
139 timeSeriesDeriv->reserve(this->m_values.size() - 1);
140 for (; it != m_values.end(); it++) {
141 TYPE v1 = it->value();
142 int64_t t1 = it->time().totalNanoseconds();
143 if (t1 != t0) {
144 double deriv = 1.e+9 * (double(v1 - v0) / double(t1 - t0));
145 auto tm = static_cast<int64_t>((t1 + t0) / 2);
146 timeSeriesDeriv->addValue(Types::Core::DateAndTime(tm), deriv);
147 }
148 t0 = t1;
149 v0 = v1;
150 }
151 return timeSeriesDeriv;
152}
154template <> std::unique_ptr<TimeSeriesProperty<double>> TimeSeriesProperty<std::string>::getDerivative() const {
155 throw std::runtime_error("Time series property derivative is not defined for strings");
157
161template <typename TYPE> size_t TimeSeriesProperty<TYPE>::getMemorySize() const {
162 // Rough estimate
163 return m_values.size() * (sizeof(TYPE) + sizeof(DateAndTime));
164}
165
181 auto const *rhs = dynamic_cast<TimeSeriesProperty<TYPE> const *>(right);
183 if (rhs) {
184 if (this->operator!=(*rhs)) {
185 m_values.insert(m_values.end(), rhs->m_values.begin(), rhs->m_values.end());
186 m_propSortedFlag = TimeSeriesSortStatus::TSUNKNOWN;
187 } else {
188 // Do nothing if appending yourself to yourself. The net result would be
189 // the same anyway
190 ;
192
193 // Count the REAL size.
194 m_size = static_cast<int>(m_values.size());
195
196 } else
197 g_log.warning() << "TimeSeriesProperty " << this->name()
198 << " could not be added to another property of the same "
199 "name but incompatible type.\n";
200
201 return *this;
202}
209template <typename TYPE> bool TimeSeriesProperty<TYPE>::operator==(const TimeSeriesProperty<TYPE> &right) const {
210 sortIfNecessary();
211
212 if (this->name() != right.name()) // should this be done?
213 {
214 return false;
215 }
216
217 if (this->m_size != right.m_size) {
218 return false;
220
221 if (this->realSize() != right.realSize()) {
222 return false;
223 } else {
224 const std::vector<DateAndTime> lhsTimes = this->timesAsVector();
225 const std::vector<DateAndTime> rhsTimes = right.timesAsVector();
226 if (!std::equal(lhsTimes.begin(), lhsTimes.end(), rhsTimes.begin())) {
227 return false;
228 }
230 const std::vector<TYPE> lhsValues = this->valuesAsVector();
231 const std::vector<TYPE> rhsValues = right.valuesAsVector();
232 if (!std::equal(lhsValues.begin(), lhsValues.end(), rhsValues.begin())) {
233 return false;
235 }
237 return true;
239
245template <typename TYPE> bool TimeSeriesProperty<TYPE>::operator==(const Property &right) const {
246 auto rhs_tsp = dynamic_cast<const TimeSeriesProperty<TYPE> *>(&right);
247 if (!rhs_tsp)
248 return false;
249 return this->operator==(*rhs_tsp);
251
257template <typename TYPE> bool TimeSeriesProperty<TYPE>::operator!=(const TimeSeriesProperty<TYPE> &right) const {
258 return !(*this == right);
259}
266template <typename TYPE> bool TimeSeriesProperty<TYPE>::operator!=(const Property &right) const {
267 return !(*this == right);
268}
269
273template <typename TYPE> void TimeSeriesProperty<TYPE>::setName(const std::string &name) { m_name = name; }
274
283template <typename TYPE>
285 std::vector<TimeValueUnit<TYPE>> &filteredData) const {
286 filteredData.clear();
288 // Expediently treat a few special cases
290 // Nothing to copy
291 if (m_values.empty()) {
292 return;
293 }
295 // Copy the only value
296 if (m_values.size() == 1) {
297 filteredData.push_back(m_values.front());
298 return;
300
301 // Copy the first value only, if all values are the same
302 // Exclude "proton_charge" logs from consideration, because in a real measurement those values can't be the same,
303 // Removing some of them just because they are equal will cause wrong total proton charge results.
304 if (allValuesAreSame(m_values) && this->name() != "proton_charge") {
305 filteredData.push_back(m_values.front());
306 return;
308
309 // Copy everything
310 if (timeROI.useAll()) {
311 std::copy(m_values.cbegin(), m_values.cend(), std::back_inserter(filteredData));
312 return;
313 }
314
315 // Copy the first value only
316 if (timeROI.useNone()) {
317 filteredData.push_back(m_values.front());
318 return;
319 }
321 // Now treat the general case
322
323 // Get all ROI time boundaries. Every other value is start/stop of an ROI "use" region.
324 const std::vector<Types::Core::DateAndTime> &roiTimes = timeROI.getAllTimes();
325 auto itROI = roiTimes.cbegin();
326 const auto itROIEnd = roiTimes.cend();
328 auto itValue = m_values.cbegin();
329 const auto itValueEnd = m_values.cend();
330 auto itLastValueUsed = itValue; // last value used up to the moment
331
332 while (itROI != itROIEnd && itValue != itValueEnd) {
333 // Try fast-forwarding the current ROI "use" region towards the current time value. Note, the current value might
334 // be in an ROI "ignore" region together with one or more following values.
335 while (std::distance(itROI, itROIEnd) > 2 && *(std::next(itROI, 2)) <= itValue->time())
336 std::advance(itROI, 2);
337 // Try finding the first value equal or past the beginning of the current ROI "use" region.
338 itValue = std::lower_bound(itValue, itValueEnd, *itROI,
339 [](const auto &value, const auto &roi_time) { return value.time() < roi_time; });
340 // Calculate a [begin,end) range for the values to use
341 auto itBeginUseValue = itValue;
342 auto itEndUseValue = itValue;
343 // If there are no values past the current ROI "use" region, get the previous value
344 if (itValue == itValueEnd) {
345 itBeginUseValue =
346 std::prev(itValue); // std::prev is safe here, because "m_values is empty" case has already been treated above
347 itEndUseValue = itValueEnd;
348 }
349 // If the value is inside the current ROI "use" region, look for other values in the same ROI "use" region
350 else if (itValue->time() <= *(std::next(itROI))) {
351 // First, try including a value immediately preceding the first value in the ROI "use" region.
352 itBeginUseValue = itValue == m_values.begin() ? itValue : std::prev(itValue);
353 // Now try finding the first value past the end of the current ROI "use" region.
354 while (itValue != itValueEnd && itValue->time() <= *(std::next(itROI)))
355 itValue++;
356 // Include the current value, therefore, advance itEndUseValue, because std::copy works as [begin,end).
357 itEndUseValue = itValue == itValueEnd ? itValue : std::next(itValue);
359 // If we are at the last ROI "use" region or the value is not past the beginning of the next ROI "use" region, keep
360 // it for the current ROI "use" region.
361 else if (std::distance(itROI, itROIEnd) == 2 ||
362 (std::distance(itROI, itROIEnd) > 2 && itValue->time() < *(std::next(itROI, 2)))) {
363 // Try including the value immediately preceding the current value
364 itBeginUseValue = itValue == m_values.begin() ? itValue : std::prev(itValue);
365 itEndUseValue = std::next(itValue);
367 // Do not use a value already copied for the previous ROI
368 if (!filteredData.empty()) {
369 itBeginUseValue = std::max(itBeginUseValue, std::next(itLastValueUsed));
371
372 // Copy all [begin,end) values and mark the last value copied
373 if (itBeginUseValue < itEndUseValue) {
374 std::copy(itBeginUseValue, itEndUseValue, std::back_inserter(filteredData));
375 itLastValueUsed = std::prev(itEndUseValue);
376 }
377
378 // Move to the next ROI "use" region
379 std::advance(itROI, 2);
380 }
381}
382
388template <typename TYPE> void TimeSeriesProperty<TYPE>::removeDataOutsideTimeROI(const TimeROI &timeROI) {
389 std::vector<TimeValueUnit<TYPE>> mp_copy;
390 createFilteredData(timeROI, mp_copy);
391
392 m_values.clear();
393 m_values = mp_copy;
394 mp_copy.clear();
395
396 m_size = static_cast<int>(m_values.size());
397}
398
399// The makeFilterByValue & expandFilterToRange methods generate a bunch of
400// warnings when the template type is the wider integer types
401// (when it's being assigned back to a double such as in a call to minValue or
402// firstValue)
403// However, in reality these methods are only used for TYPE=int or double (they
404// are only called from FilterByLogValue) so suppress the warnings
405#ifdef _WIN32
406#pragma warning(push)
407#pragma warning(disable : 4244)
408#pragma warning(disable : 4804) // This one comes about for TYPE=bool - again
409 // the method is never called for this type
410#endif
411#if defined(__GNUC__) && !(defined(__INTEL_COMPILER))
412#pragma GCC diagnostic ignored "-Wconversion"
413#endif
414
428template <typename TYPE>
429void TimeSeriesProperty<TYPE>::makeFilterByValue(std::vector<SplittingInterval> &split, double min, double max,
430 double TimeTolerance, bool centre) const {
431 const bool emptyMin = (min == EMPTY_DBL());
432 const bool emptyMax = (max == EMPTY_DBL());
433
434 if (!emptyMin && !emptyMax && max < min) {
435 std::stringstream ss;
436 ss << "TimeSeriesProperty::makeFilterByValue: 'max' argument must be "
437 "greater than 'min' "
438 << "(got min=" << min << " max=" << max << ")";
439 throw std::invalid_argument(ss.str());
440 }
441
442 // If min or max were unset ("empty") in the algorithm, set to the min or max
443 // value of the log
444 if (emptyMin)
445 min = static_cast<double>(minValue());
446 if (emptyMax)
447 max = static_cast<double>(maxValue());
448
449 // Make sure the splitter starts out empty
450 split.clear();
451
452 // Do nothing if the log is empty.
453 if (m_values.empty())
454 return;
455
456 // 1. Sort
457 sortIfNecessary();
458
459 // 2. Do the rest
460 bool lastGood(false);
461 time_duration tol = DateAndTime::durationFromSeconds(TimeTolerance);
462 int numgood = 0;
463 DateAndTime t;
464 DateAndTime start, stop;
465
466 for (size_t i = 0; i < m_values.size(); ++i) {
467 const DateAndTime lastGoodTime = t;
468 // The new entry
469 t = m_values[i].time();
470 TYPE val = m_values[i].value();
471
472 // A good value?
473 const bool isGood = ((val >= min) && (val <= max));
474 if (isGood)
475 numgood++;
476
477 if (isGood != lastGood) {
478 // We switched from bad to good or good to bad
479
480 if (isGood) {
481 // Start of a good section. Subtract tolerance from the time if
482 // boundaries are centred.
483 start = centre ? t - tol : t;
484 } else {
485 // End of the good section. Add tolerance to the LAST GOOD time if
486 // boundaries are centred.
487 // Otherwise, use the first 'bad' time.
488 stop = centre ? lastGoodTime + tol : t;
489 split.emplace_back(start, stop, 0);
490 // Reset the number of good ones, for next time
491 numgood = 0;
492 }
493 lastGood = isGood;
494 }
495 }
496
497 if (numgood > 0) {
498 // The log ended on "good" so we need to close it using the last time we
499 // found
500 stop = t + tol;
501 split.emplace_back(start, stop, 0);
502 }
503}
504
508template <>
509void TimeSeriesProperty<std::string>::makeFilterByValue(std::vector<SplittingInterval> & /*split*/, double /*min*/,
510 double /*max*/, double /*TimeTolerance*/,
511 bool /*centre*/) const {
512 throw Exception::NotImplementedError("TimeSeriesProperty::makeFilterByValue "
513 "is not implemented for string "
514 "properties");
515}
516
533template <typename TYPE>
534TimeROI TimeSeriesProperty<TYPE>::makeFilterByValue(double min, double max, bool expand,
535 const TimeInterval &expandRange, double TimeTolerance, bool centre,
536 const TimeROI *existingROI) const {
537 const bool emptyMin = (min == EMPTY_DBL());
538 const bool emptyMax = (max == EMPTY_DBL());
539
540 if (!emptyMin && !emptyMax && max < min) {
541 std::stringstream ss;
542 ss << "TimeSeriesProperty::makeFilterByValue: 'max' argument must be "
543 "greater than 'min' "
544 << "(got min=" << min << " max=" << max << ")";
545 throw std::invalid_argument(ss.str());
546 }
547
548 // If min or max were unset ("empty") in the algorithm, set to the min or max
549 // value of the log
550 if (emptyMin)
551 min = static_cast<double>(minValue());
552 if (emptyMax)
553 max = static_cast<double>(maxValue());
554
555 TimeROI newROI;
556
557 // Do nothing if the log is empty.
558 if (m_values.empty())
559 return newROI;
560
561 // 1. Sort
562 sortIfNecessary();
563
564 // 2. Do the rest
565 const time_duration tol = DateAndTime::durationFromSeconds(TimeTolerance);
566 DateAndTime stop_t;
567 DateAndTime start, stop;
568
569 bool isGood = false;
570 for (size_t i = 0; i < m_values.size(); ++i) {
571 TYPE val = m_values[i].value();
572
573 if ((val >= min) && (val <= max)) {
574 if (isGood) {
575 stop_t = m_values[i].time();
576 } else {
577 isGood = true;
578 stop_t = m_values[i].time();
579 start = centre ? m_values[i].time() - tol : m_values[i].time();
580 }
581 } else if (isGood) {
582 stop = centre ? stop_t + tol : m_values[i].time();
583 if (start < stop)
584 newROI.addROI(start, stop);
585 isGood = false;
586 }
587 }
588 if (isGood) {
589 stop = centre ? stop_t + tol : stop_t;
590 if (start < stop)
591 newROI.addROI(start, stop);
592 }
593
594 if (expand) {
595 if (expandRange.start() < firstTime()) {
596 auto val = static_cast<double>(firstValue());
597 if ((val >= min) && (val <= max)) {
598 newROI.addROI(expandRange.start(), firstTime());
599 }
600 }
601 if (expandRange.stop() > lastTime()) {
602 auto val = static_cast<double>(lastValue());
603 if ((val >= min) && (val <= max)) {
604 newROI.addROI(lastTime(), expandRange.stop());
605 }
606 }
607 }
608
609 // If the TimeROI is empty there are no values inside the filter
610 // so we should return USE_NONE
611 if (newROI.useAll()) {
612 return TimeROI::USE_NONE;
613 }
614
615 if (existingROI != nullptr && !existingROI->useAll()) {
616 newROI.update_intersection(*existingROI);
617 }
618 return newROI;
619}
620
624template <>
625TimeROI TimeSeriesProperty<std::string>::makeFilterByValue(double /*min*/, double /*max*/, bool /*expand*/,
626 const TimeInterval & /*expandRange*/,
627 double /*TimeTolerance*/, bool /*centre*/,
628 const TimeROI * /*existingROI*/) const {
629 throw Exception::NotImplementedError("TimeSeriesProperty::makeFilterByValue "
630 "is not implemented for string "
631 "properties");
632}
633
644template <typename TYPE>
645void TimeSeriesProperty<TYPE>::expandFilterToRange(std::vector<SplittingInterval> &split, double min, double max,
646 const TimeInterval &range) const {
647 const bool emptyMin = (min == EMPTY_DBL());
648 const bool emptyMax = (max == EMPTY_DBL());
649
650 if (!emptyMin && !emptyMax && max < min) {
651 std::stringstream ss;
652 ss << "TimeSeriesProperty::expandFilterToRange: 'max' argument must be "
653 "greater than 'min' "
654 << "(got min=" << min << " max=" << max << ")";
655 throw std::invalid_argument(ss.str());
656 }
657
658 // If min or max were unset ("empty") in the algorithm, set to the min or max
659 // value of the log
660 if (emptyMin)
661 min = static_cast<double>(minValue());
662 if (emptyMax)
663 max = static_cast<double>(maxValue());
664
665 if (range.start() < firstTime()) {
666 // Assume everything before the 1st value is constant
667 double val = static_cast<double>(firstValue());
668 if ((val >= min) && (val <= max)) {
669 SplittingIntervalVec extraFilter;
670 extraFilter.emplace_back(range.start(), firstTime(), 0);
671 // Include everything from the start of the run to the first time measured
672 // (which may be a null time interval; this'll be ignored)
673 split = split | extraFilter;
674 }
675 }
676
677 if (lastTime() < range.stop()) {
678 // Assume everything after the LAST value is constant
679 double val = static_cast<double>(lastValue());
680 if ((val >= min) && (val <= max)) {
681 SplittingIntervalVec extraFilter;
682 extraFilter.emplace_back(lastTime(), range.stop(), 0);
683 // Include everything from the start of the run to the first time measured
684 // (which may be a null time interval; this'll be ignored)
685 split = split | extraFilter;
686 }
687 }
688}
689
693template <>
694void TimeSeriesProperty<std::string>::expandFilterToRange(std::vector<SplittingInterval> & /*split*/, double /*min*/,
695 double /*max*/, const TimeInterval & /*range*/) const {
696 throw Exception::NotImplementedError("TimeSeriesProperty::makeFilterByValue "
697 "is not implemented for string "
698 "properties");
699}
700
705template <typename TYPE> double TimeSeriesProperty<TYPE>::timeAverageValue(const TimeROI *timeRoi) const {
706 double retVal = 0.0;
707 try {
708 if ((timeRoi == nullptr) || (timeRoi->useAll())) {
709 const auto &intervals = getTimeIntervals();
710 retVal = this->averageValueInFilter(intervals);
711 } else if (timeRoi->useNone()) {
712 // if TimeROI bans everything, use the simple mean
713 const auto stats =
715 return stats.mean;
716 } else {
717 const auto &filter = timeRoi->toTimeIntervals();
718 retVal = this->averageValueInFilter(filter);
719 g_log.warning("Calls to TimeSeriesProperty::timeAverageValue should be replaced with "
720 "Run::getTimeAveragedValue");
721 }
722 } catch (std::exception &) {
723 // just return nan
724 retVal = std::numeric_limits<double>::quiet_NaN();
725 }
726 return retVal;
727}
728
729template <> double TimeSeriesProperty<std::string>::timeAverageValue(const TimeROI * /*timeRoi*/) const {
730 throw Exception::NotImplementedError("TimeSeriesProperty::timeAverageValue is not implemented for string properties");
731}
732
739template <typename TYPE>
740double TimeSeriesProperty<TYPE>::averageValueInFilter(const std::vector<TimeInterval> &filter) const {
741 // TODO: Consider logs that aren't giving starting values.
742
743 // If there's just a single value in the log, return that.
744 if (size() == 1) {
745 return static_cast<double>(this->firstValue());
746 }
747
748 // First of all, if the log or the filter is empty, return NaN
749 if (realSize() == 0 || filter.empty()) {
750 return std::numeric_limits<double>::quiet_NaN();
751 }
752
753 sortIfNecessary();
754
755 double numerator(0.0), totalTime(0.0);
756 // Loop through the filter ranges
757 for (const auto &time : filter) {
758 // Calculate the total time duration (in seconds) within by the filter
759 totalTime += time.duration();
760
761 // Get the log value and index at the start time of the filter
762 int index;
763 double currentValue = static_cast<double>(getSingleValue(time.start(), index));
764 DateAndTime startTime = time.start();
765
766 while (index < realSize() - 1 && m_values[index + 1].time() < time.stop()) {
767 ++index;
768 numerator += DateAndTime::secondsFromDuration(m_values[index].time() - startTime) * currentValue;
769 startTime = m_values[index].time();
770 currentValue = static_cast<double>(m_values[index].value());
771 }
772
773 // Now close off with the end of the current filter range
774 numerator += DateAndTime::secondsFromDuration(time.stop() - startTime) * currentValue;
775 }
776
777 if (totalTime > 0) {
778 // 'Normalise' by the total time
779 return numerator / totalTime;
780 } else {
781 // give simple mean
782 const auto stats = Mantid::Kernel::getStatistics(this->valuesAsVector(), Mantid::Kernel::Math::StatisticType::Mean);
783 return stats.mean;
784 }
785}
786
790template <>
791double TimeSeriesProperty<std::string>::averageValueInFilter(const std::vector<TimeInterval> & /*filter*/) const {
792 throw Exception::NotImplementedError("TimeSeriesProperty::"
793 "averageValueInFilter is not "
794 "implemented for string properties");
795}
796
797template <typename TYPE>
798std::pair<double, double>
799TimeSeriesProperty<TYPE>::averageAndStdDevInFilter(const std::vector<TimeInterval> &intervals) const {
800 double mean_prev, mean_current(0.0), s(0.0), variance, duration, weighted_sum(0.0);
801
802 // First of all, if the log or the intervals are empty or is a single value,
803 // return NaN for the uncertainty
804 if (realSize() <= 1 || intervals.empty()) {
805 return std::pair<double, double>{this->averageValueInFilter(intervals), std::numeric_limits<double>::quiet_NaN()};
806 }
807 auto real_size = realSize();
808 for (const auto &time : intervals) {
809 int index;
810 auto currentValue = static_cast<double>(getSingleValue(time.start(), index));
811 DateAndTime startTime = time.start();
812 while (index < realSize() - 1 && m_values[index + 1].time() < time.stop()) {
813 index++;
814 if (index == real_size) {
815 duration = DateAndTime::secondsFromDuration(time.stop() - startTime);
816 } else {
817 duration = DateAndTime::secondsFromDuration(m_values[index].time() - startTime);
818 startTime = m_values[index].time();
819 }
820 mean_prev = mean_current;
821 if (duration > 0.) {
822 weighted_sum += duration;
823
824 mean_current = mean_prev + (duration / weighted_sum) * (currentValue - mean_prev);
825 s += duration * (currentValue - mean_prev) * (currentValue - mean_current);
826 }
827 currentValue = static_cast<double>(m_values[index].value());
828 }
829
830 // Now close off with the end of the current filter range
831 duration = DateAndTime::secondsFromDuration(time.stop() - startTime);
832 if (duration > 0.) {
833 weighted_sum += duration;
834 mean_prev = mean_current;
835
836 mean_current = mean_prev + (duration / weighted_sum) * (currentValue - mean_prev);
837 s += duration * (currentValue - mean_prev) * (currentValue - mean_current);
838 }
839 }
840 variance = s / weighted_sum;
841 // Normalise by the total time
842 return std::pair<double, double>{mean_current, std::sqrt(variance)};
843}
844
848template <>
849std::pair<double, double>
850TimeSeriesProperty<std::string>::averageAndStdDevInFilter(const std::vector<TimeInterval> & /*filter*/) const {
851 throw Exception::NotImplementedError("TimeSeriesProperty::"
852 "averageAndStdDevInFilter is not "
853 "implemented for string properties");
854}
855
856template <typename TYPE>
857std::pair<double, double> TimeSeriesProperty<TYPE>::timeAverageValueAndStdDev(const Kernel::TimeROI *timeRoi) const {
858 // time series with less than two entries are conner cases
859 if (this->realSize() == 0)
860 return std::pair<double, double>{std::numeric_limits<double>::quiet_NaN(),
861 std::numeric_limits<double>::quiet_NaN()};
862 else if (this->realSize() == 1)
863 return std::pair<double, double>(static_cast<double>(this->firstValue()), 0.0);
864
865 // Derive splitting intervals from either the roi or from the first/last entries in the time series
866 std::vector<TimeInterval> intervals;
867 if (timeRoi && !timeRoi->useAll()) {
868 intervals = timeRoi->toTimeIntervals(this->firstTime());
869 } else {
870 intervals = this->getTimeIntervals();
871 }
872
873 return this->averageAndStdDevInFilter(intervals);
874}
875
879template <>
880std::pair<double, double>
883 "TimeSeriesProperty::timeAverageValueAndStdDev is not implemented for string properties");
884}
885
886// Re-enable the warnings disabled before makeFilterByValue
887#ifdef _WIN32
888#pragma warning(pop)
889#endif
890#if defined(__GNUC__) && !(defined(__INTEL_COMPILER))
891#pragma GCC diagnostic warning "-Wconversion"
892#endif
893
900template <typename TYPE> std::map<DateAndTime, TYPE> TimeSeriesProperty<TYPE>::valueAsCorrectMap() const {
901 // 1. Sort if necessary
902 sortIfNecessary();
903
904 // 2. Data Strcture
905 std::map<DateAndTime, TYPE> asMap;
906
907 if (!m_values.empty()) {
908 for (size_t i = 0; i < m_values.size(); i++)
909 asMap[m_values[i].time()] = m_values[i].value();
910 }
911
912 return asMap;
913}
914
919template <typename TYPE> std::vector<TYPE> TimeSeriesProperty<TYPE>::valuesAsVector() const {
920 sortIfNecessary();
921
922 std::vector<TYPE> out;
923 out.reserve(m_values.size());
924
925 for (size_t i = 0; i < m_values.size(); i++)
926 out.emplace_back(m_values[i].value());
927
928 return out;
929}
930
937template <typename TYPE> std::multimap<DateAndTime, TYPE> TimeSeriesProperty<TYPE>::valueAsMultiMap() const {
938 std::multimap<DateAndTime, TYPE> asMultiMap;
939
940 if (!m_values.empty()) {
941 for (size_t i = 0; i < m_values.size(); i++)
942 asMultiMap.insert(std::make_pair(m_values[i].time(), m_values[i].value()));
943 }
944
945 return asMultiMap;
946}
947
952template <typename TYPE> std::vector<DateAndTime> TimeSeriesProperty<TYPE>::timesAsVector() const {
953 sortIfNecessary();
954
955 std::vector<DateAndTime> out;
956 out.reserve(m_values.size());
957
958 for (size_t i = 0; i < m_values.size(); i++) {
959 out.emplace_back(m_values[i].time());
960 }
961
962 return out;
963}
964
974template <typename TYPE>
975std::vector<DateAndTime> TimeSeriesProperty<TYPE>::filteredTimesAsVector(const Kernel::TimeROI *roi) const {
976 if (roi && !roi->useAll()) {
977 this->sortIfNecessary();
978 std::vector<DateAndTime> filteredTimes;
979 if (roi->firstTime() > this->m_values.back().time()) {
980 // Since the ROI starts after everything, just return the last time in the log
981 filteredTimes.emplace_back(roi->firstTime());
982 } else { // only use the times in the filter - this is very similar to FilteredTimeSeriesProperty::applyFilter
983 // the index into the m_values array of the time, or -1 (before) or m_values.size() (after)
984 std::size_t index_current_log{0};
985
986 for (const auto &splitter : roi->toTimeIntervals()) {
987 const auto endTime = splitter.stop();
988
989 // check if the splitter starts too early
990 if (endTime < this->m_values[index_current_log].time()) {
991 continue; // skip to the next splitter
992 }
993
994 // cache values to reduce number of method calls
995 const auto beginTime = splitter.start();
996
997 // find the first log that should be added
998 if (this->m_values.back().time() < beginTime) {
999 // skip directly to the end if the filter starts after the last log
1000 index_current_log = this->m_values.size() - 1;
1001 } else {
1002 // search for the right starting point
1003 while ((this->m_values[index_current_log].time() <= beginTime)) {
1004 if (index_current_log + 1 > this->m_values.size())
1005 break;
1006 index_current_log++;
1007 }
1008 // need to back up by one
1009 if (index_current_log > 0)
1010 index_current_log--;
1011 // go backwards more while times are equal to the one being started at
1012 while (index_current_log > 0 &&
1013 this->m_values[index_current_log].time() == this->m_values[index_current_log - 1].time()) {
1014 index_current_log--;
1015 }
1016 }
1017
1018 // add everything up to the end time
1019 for (; index_current_log < this->m_values.size(); ++index_current_log) {
1020 if (this->m_values[index_current_log].time() >= endTime)
1021 break;
1022
1023 // start time is when this value was created or when the filter started
1024 filteredTimes.emplace_back(std::max(beginTime, this->m_values[index_current_log].time()));
1025 }
1026 // go back one so the next splitter can add a value
1027 if (index_current_log > 0)
1028 index_current_log--;
1029 }
1030 }
1031
1032 return filteredTimes;
1033 } else {
1034 return this->timesAsVector();
1035 }
1036}
1037
1038template <typename TYPE> std::vector<DateAndTime> TimeSeriesProperty<TYPE>::filteredTimesAsVector() const {
1039 return this->timesAsVector();
1040}
1041
1046template <typename TYPE> std::vector<double> TimeSeriesProperty<TYPE>::timesAsVectorSeconds() const {
1047 // 1. Sort if necessary
1048 sortIfNecessary();
1049
1050 // 2. Output data structure
1051 std::vector<double> out;
1052 if (!m_values.empty()) {
1053 out.reserve(m_values.size());
1054
1055 Types::Core::DateAndTime start = m_values[0].time();
1056 for (size_t i = 0; i < m_values.size(); i++) {
1057 out.emplace_back(DateAndTime::secondsFromDuration(m_values[i].time() - start));
1058 }
1059 }
1060
1061 return out;
1062}
1063
1069template <typename TYPE>
1070std::vector<double> TimeSeriesProperty<TYPE>::timesAsVectorSeconds(Types::Core::DateAndTime start) const {
1071 // 1. Sort if necessary
1072 sortIfNecessary();
1073
1074 // 2. Output data structure
1075 std::vector<double> out;
1076 if (!m_values.empty()) {
1077 out.reserve(m_values.size());
1078 for (size_t i = 0; i < m_values.size(); i++) {
1079 out.emplace_back(DateAndTime::secondsFromDuration(m_values[i].time() - start));
1080 }
1081 }
1082
1083 return out;
1084}
1085
1091template <typename TYPE>
1092void TimeSeriesProperty<TYPE>::addValue(const Types::Core::DateAndTime &time, const TYPE &value) {
1093 TimeValueUnit<TYPE> newvalue(time, value);
1094 // Add the value to the back of the vector
1095 m_values.emplace_back(newvalue);
1096 // Increment the separate record of the property's size
1097 m_size++;
1098
1099 // Toggle the sorted flag if necessary
1100 // (i.e. if the flag says we're sorted and the added time is before the prior
1101 // last time)
1102 if (m_size == 1) {
1103 // First item, must be sorted.
1104 m_propSortedFlag = TimeSeriesSortStatus::TSSORTED;
1105 } else if (m_propSortedFlag == TimeSeriesSortStatus::TSUNKNOWN && m_values.back() < *(m_values.rbegin() + 1)) {
1106 // Previously unknown and still unknown
1107 m_propSortedFlag = TimeSeriesSortStatus::TSUNSORTED;
1108 } else if (m_propSortedFlag == TimeSeriesSortStatus::TSSORTED && m_values.back() < *(m_values.rbegin() + 1)) {
1109 // Previously sorted but last added is not in order
1110 m_propSortedFlag = TimeSeriesSortStatus::TSUNSORTED;
1111 }
1112
1113 // m_filterApplied = false;
1114}
1115
1121template <typename TYPE> void TimeSeriesProperty<TYPE>::addValue(const std::string &time, const TYPE &value) {
1122 return addValue(Types::Core::DateAndTime(time), value);
1123}
1124
1130template <typename TYPE> void TimeSeriesProperty<TYPE>::addValue(const std::time_t &time, const TYPE &value) {
1131 Types::Core::DateAndTime dt;
1132 dt.set_from_time_t(time);
1133 return addValue(dt, value);
1134}
1135
1141template <typename TYPE>
1142void TimeSeriesProperty<TYPE>::addValues(const std::vector<Types::Core::DateAndTime> &times,
1143 const std::vector<TYPE> &values) {
1144 size_t length = std::min(times.size(), values.size());
1145 m_size += static_cast<int>(length);
1146 for (size_t i = 0; i < length; ++i) {
1147 m_values.emplace_back(times[i], values[i]);
1148 }
1149
1150 if (!values.empty())
1151 m_propSortedFlag = TimeSeriesSortStatus::TSUNKNOWN;
1152}
1153
1159template <typename TYPE>
1160void TimeSeriesProperty<TYPE>::replaceValues(const std::vector<Types::Core::DateAndTime> &times,
1161 const std::vector<TYPE> &values) {
1162 clear();
1163 addValues(times, values);
1164}
1165
1170template <typename TYPE> DateAndTime TimeSeriesProperty<TYPE>::lastTime() const {
1171 if (m_values.empty()) {
1172 const std::string error("lastTime(): TimeSeriesProperty '" + name() + "' is empty");
1173 g_log.debug(error);
1174 throw std::runtime_error(error);
1175 }
1176
1177 sortIfNecessary();
1178
1179 return m_values.rbegin()->time();
1180}
1181
1185template <typename TYPE> TYPE TimeSeriesProperty<TYPE>::firstValue() const {
1186 if (m_values.empty()) {
1187 const std::string error("firstValue(): TimeSeriesProperty '" + name() + "' is empty");
1188 g_log.debug(error);
1189 throw std::runtime_error(error);
1190 }
1191
1192 sortIfNecessary();
1193
1194 return m_values[0].value();
1195}
1196
1197template <typename TYPE> TYPE TimeSeriesProperty<TYPE>::firstValue(const Kernel::TimeROI &roi) const {
1198 const auto startTime = roi.firstTime();
1199 if (startTime <= this->firstTime()) {
1200 return this->firstValue();
1201 } else if (startTime >= this->lastTime()) {
1202 return this->lastValue();
1203 } else {
1204 const auto times = this->timesAsVector();
1205 auto iter = std::lower_bound(times.cbegin(), times.cend(), startTime);
1206 if (*iter > startTime)
1207 iter--;
1208 const auto index = std::size_t(std::distance(times.cbegin(), iter));
1209
1210 const auto values = this->valuesAsVector();
1211 const TYPE ret = values[index];
1212 return ret;
1213 }
1214}
1215
1219template <typename TYPE> DateAndTime TimeSeriesProperty<TYPE>::firstTime() const {
1220 if (m_values.empty()) {
1221 const std::string error("firstTime(): TimeSeriesProperty '" + name() + "' is empty");
1222 g_log.debug(error);
1223 throw std::runtime_error(error);
1224 }
1225
1226 sortIfNecessary();
1227
1228 return m_values[0].time();
1229}
1230
1235template <typename TYPE> TYPE TimeSeriesProperty<TYPE>::lastValue() const {
1236 if (m_values.empty()) {
1237 const std::string error("lastValue(): TimeSeriesProperty '" + name() + "' is empty");
1238 g_log.debug(error);
1239 throw std::runtime_error(error);
1240 }
1241
1242 sortIfNecessary();
1243
1244 return m_values.rbegin()->value();
1245}
1246
1247template <typename TYPE> TYPE TimeSeriesProperty<TYPE>::lastValue(const Kernel::TimeROI &roi) const {
1248 const auto stopTime = roi.lastTime();
1249 const auto times = this->timesAsVector();
1250 if (stopTime <= times.front()) {
1251 return this->firstValue();
1252 } else if (stopTime >= times.back()) {
1253 return this->lastValue();
1254 } else {
1255 auto iter = std::lower_bound(times.cbegin(), times.cend(), stopTime);
1256 if ((iter != times.cbegin()) && (*iter > stopTime))
1257 --iter;
1258 const auto index = std::size_t(std::distance(times.cbegin(), iter));
1259
1260 const auto values = this->valuesAsVector();
1261 const TYPE ret = values[index];
1262 return ret;
1263 }
1264}
1265
1275template <typename TYPE> double TimeSeriesProperty<TYPE>::durationInSeconds(const Kernel::TimeROI *roi) const {
1276 if (this->size() == 0)
1277 return std::numeric_limits<double>::quiet_NaN();
1278 if (roi && !roi->useAll()) {
1279 Kernel::TimeROI seriesSpan(*roi);
1280 const auto thisFirstTime = this->firstTime();
1281 // remove everything before the start time
1282 if (thisFirstTime > DateAndTime::GPS_EPOCH) {
1283 seriesSpan.addMask(DateAndTime::GPS_EPOCH, thisFirstTime);
1284 }
1285 return seriesSpan.durationInSeconds();
1286 } else {
1287 const auto &intervals = this->getTimeIntervals();
1288 const double duration_sec =
1289 std::accumulate(intervals.cbegin(), intervals.cend(), 0.,
1290 [](double sum, const auto &interval) { return sum + interval.duration(); });
1291 return duration_sec;
1292 }
1293}
1294
1295template <typename TYPE> TYPE TimeSeriesProperty<TYPE>::minValue() const {
1296 return std::min_element(m_values.begin(), m_values.end(), TimeValueUnit<TYPE>::valueCmp)->value();
1297}
1298
1299template <typename TYPE> TYPE TimeSeriesProperty<TYPE>::maxValue() const {
1300 return std::max_element(m_values.begin(), m_values.end(), TimeValueUnit<TYPE>::valueCmp)->value();
1301}
1302
1303template <typename TYPE> double TimeSeriesProperty<TYPE>::mean() const {
1304 Mantid::Kernel::Statistics raw_stats =
1305 Mantid::Kernel::getStatistics(this->filteredValuesAsVector(), StatOptions::Mean);
1306 return raw_stats.mean;
1307}
1308
1311template <typename TYPE> int TimeSeriesProperty<TYPE>::size() const { return m_size; }
1312
1317template <typename TYPE> int TimeSeriesProperty<TYPE>::realSize() const { return static_cast<int>(m_values.size()); }
1318
1319/*
1320 * Get the time series property as a string of 'time value'
1321 * @return time series property as a string
1322 */
1323template <typename TYPE> std::string TimeSeriesProperty<TYPE>::value() const {
1324 sortIfNecessary();
1325
1326 std::stringstream ins;
1327 for (size_t i = 0; i < m_values.size(); i++) {
1328 try {
1329 ins << m_values[i].time().toSimpleString();
1330 ins << " " << m_values[i].value() << "\n";
1331 } catch (...) {
1332 // Some kind of error; for example, invalid year, can occur when
1333 // converting boost time.
1334 ins << "Error Error"
1335 << "\n";
1336 }
1337 }
1338
1339 return ins.str();
1340}
1341
1346template <typename TYPE> std::vector<std::string> TimeSeriesProperty<TYPE>::time_tValue() const {
1347 sortIfNecessary();
1348
1349 std::vector<std::string> values;
1350 values.reserve(m_values.size());
1351
1352 for (size_t i = 0; i < m_values.size(); i++) {
1353 std::stringstream line;
1354 line << m_values[i].time().toSimpleString() << " " << m_values[i].value();
1355 values.emplace_back(line.str());
1356 }
1357
1358 return values;
1359}
1360
1368template <typename TYPE> std::map<DateAndTime, TYPE> TimeSeriesProperty<TYPE>::valueAsMap() const {
1369 // 1. Sort if necessary
1370 sortIfNecessary();
1371
1372 // 2. Build map
1373
1374 std::map<DateAndTime, TYPE> asMap;
1375 if (m_values.empty())
1376 return asMap;
1377
1378 TYPE d = m_values[0].value();
1379 asMap[m_values[0].time()] = d;
1380
1381 for (size_t i = 1; i < m_values.size(); i++) {
1382 if (m_values[i].value() != d) {
1383 // Only put entry with different value from last entry to map
1384 asMap[m_values[i].time()] = m_values[i].value();
1385 d = m_values[i].value();
1386 }
1387 }
1388 return asMap;
1389}
1390
1396template <typename TYPE> std::string TimeSeriesProperty<TYPE>::setValue(const std::string & /*unused*/) {
1397 throw Exception::NotImplementedError("TimeSeriesProperty<TYPE>::setValue - "
1398 "Cannot extract TimeSeries from a "
1399 "std::string");
1400}
1401
1407template <typename TYPE> std::string TimeSeriesProperty<TYPE>::setValueFromJson(const Json::Value & /*unused*/) {
1408 throw Exception::NotImplementedError("TimeSeriesProperty<TYPE>::setValue - "
1409 "Cannot extract TimeSeries from a "
1410 "Json::Value");
1411}
1412
1417template <typename TYPE>
1418std::string TimeSeriesProperty<TYPE>::setDataItem(const std::shared_ptr<DataItem> & /*unused*/) {
1419 throw Exception::NotImplementedError("TimeSeriesProperty<TYPE>::setValue - "
1420 "Cannot extract TimeSeries from "
1421 "DataItem");
1422}
1423
1426template <typename TYPE> void TimeSeriesProperty<TYPE>::clear() {
1427 m_size = 0;
1428 m_values.clear();
1429
1430 m_propSortedFlag = TimeSeriesSortStatus::TSSORTED;
1431 // m_filterApplied = false;
1432}
1433
1441template <typename TYPE> void TimeSeriesProperty<TYPE>::clearOutdated() {
1442 if (realSize() > 1) {
1443 auto lastValueInVec = m_values.back();
1444 clear();
1445 m_values.emplace_back(lastValueInVec);
1446 m_size = 1;
1447 }
1448}
1449
1450//--------------------------------------------------------------------------------------------
1460template <typename TYPE>
1461void TimeSeriesProperty<TYPE>::create(const Types::Core::DateAndTime &start_time, const std::vector<double> &time_sec,
1462 const std::vector<TYPE> &new_values) {
1463 if (time_sec.size() != new_values.size()) {
1464 std::stringstream msg;
1465 msg << "TimeSeriesProperty \"" << name() << "\" create: mismatched size "
1466 << "for the time and values vectors.";
1467 throw std::invalid_argument(msg.str());
1468 }
1469
1470 clear();
1471 const std::size_t num = new_values.size();
1472 m_values.reserve(num);
1473
1474 // set the sorted flag
1475 if (std::is_sorted(time_sec.cbegin(), time_sec.cend()))
1476 m_propSortedFlag = TimeSeriesSortStatus::TSSORTED;
1477 else
1478 m_propSortedFlag = TimeSeriesSortStatus::TSUNSORTED;
1479 // set the values
1480 constexpr double SEC_TO_NANO{1000000000.0};
1481 const int64_t start_time_ns = start_time.totalNanoseconds();
1482 for (std::size_t i = 0; i < num; i++) {
1483 m_values.emplace_back(start_time_ns + static_cast<int64_t>(time_sec[i] * SEC_TO_NANO), new_values[i]);
1484 }
1485
1486 // reset the size
1487 m_size = static_cast<int>(m_values.size());
1488}
1489
1490//--------------------------------------------------------------------------------------------
1499template <typename TYPE>
1500void TimeSeriesProperty<TYPE>::create(const std::vector<Types::Core::DateAndTime> &new_times,
1501 const std::vector<TYPE> &new_values) {
1502 if (new_times.size() != new_values.size())
1503 throw std::invalid_argument("TimeSeriesProperty::create: mismatched size "
1504 "for the time and values vectors.");
1505
1506 clear();
1507 // nothing to do without values
1508 if (new_times.empty()) {
1509 m_propSortedFlag = TimeSeriesSortStatus::TSSORTED;
1510 return;
1511 }
1512
1513 const std::size_t num = new_values.size();
1514 m_values.reserve(num);
1515
1516 // set the sorted flag
1517 if (std::is_sorted(new_times.cbegin(), new_times.cend()))
1518 m_propSortedFlag = TimeSeriesSortStatus::TSSORTED;
1519 else
1520 m_propSortedFlag = TimeSeriesSortStatus::TSUNSORTED;
1521 // add the values
1522 for (std::size_t i = 0; i < num; i++) {
1523 m_values.emplace_back(new_times[i], new_values[i]);
1524 }
1525
1526 // reset the size
1527 m_size = static_cast<int>(m_values.size());
1528}
1529
1534template <typename TYPE> TYPE TimeSeriesProperty<TYPE>::getSingleValue(const Types::Core::DateAndTime &t) const {
1535 if (m_values.empty()) {
1536 const std::string error("getSingleValue(): TimeSeriesProperty '" + name() + "' is empty");
1537 g_log.debug(error);
1538 throw std::runtime_error(error);
1539 }
1540
1541 // 1. Get sorted
1542 sortIfNecessary();
1543
1544 // 2.
1545 TYPE valueAtTime;
1546 if (t < m_values[0].time()) {
1547 // 1. Out side of lower bound
1548 valueAtTime = m_values[0].value();
1549 } else if (t >= m_values.back().time()) {
1550 // 2. Out side of upper bound
1551 valueAtTime = m_values.back().value();
1552 } else {
1553 // 3. Within boundary
1554 int index = this->findIndex(t);
1555
1556 if (index < 0) {
1557 // If query time "t" is earlier than the begin time of the series
1558 index = 0;
1559 } else if (index == int(m_values.size())) {
1560 // If query time "t" is later than the end time of the series
1561 index = static_cast<int>(m_values.size()) - 1;
1562 } else if (index > int(m_values.size())) {
1563 std::stringstream errss;
1564 errss << "TimeSeriesProperty.findIndex() returns index (" << index << " ) > maximum defined value "
1565 << m_values.size();
1566 throw std::logic_error(errss.str());
1567 }
1568
1569 valueAtTime = m_values[static_cast<size_t>(index)].value();
1570 }
1571
1572 return valueAtTime;
1573} // END-DEF getSinglevalue()
1574
1580template <typename TYPE>
1581TYPE TimeSeriesProperty<TYPE>::getSingleValue(const Types::Core::DateAndTime &t, int &index) const {
1582 if (m_values.empty()) {
1583 const std::string error("getSingleValue(): TimeSeriesProperty '" + name() + "' is empty");
1584 g_log.debug(error);
1585 throw std::runtime_error(error);
1586 }
1587
1588 // 1. Get sorted
1589 sortIfNecessary();
1590
1591 // 2.
1592 TYPE valueAtTime;
1593 if (t < m_values[0].time()) {
1594 // 1. Out side of lower bound
1595 valueAtTime = m_values[0].value();
1596 index = 0;
1597 } else if (t >= m_values.back().time()) {
1598 // 2. Out side of upper bound
1599 valueAtTime = m_values.back().value();
1600 index = int(m_values.size()) - 1;
1601 } else {
1602 // 3. Within boundary
1603 index = this->findIndex(t);
1604
1605 if (index < 0) {
1606 // If query time "t" is earlier than the begin time of the series
1607 index = 0;
1608 } else if (index == int(m_values.size())) {
1609 // If query time "t" is later than the end time of the series
1610 index = static_cast<int>(m_values.size()) - 1;
1611 } else if (index > int(m_values.size())) {
1612 std::stringstream errss;
1613 errss << "TimeSeriesProperty.findIndex() returns index (" << index << " ) > maximum defined value "
1614 << m_values.size();
1615 throw std::logic_error(errss.str());
1616 }
1617
1618 valueAtTime = m_values[static_cast<size_t>(index)].value();
1619 }
1620
1621 return valueAtTime;
1622} // END-DEF getSinglevalue()
1623
1634template <typename TYPE> TimeInterval TimeSeriesProperty<TYPE>::nthInterval(int n) const {
1635 // 0. Throw exception
1636 if (m_values.empty()) {
1637 const std::string error("nthInterval(): TimeSeriesProperty '" + name() + "' is empty");
1638 g_log.debug(error);
1639 throw std::runtime_error(error);
1640 }
1641
1642 // 1. Sort
1643 sortIfNecessary();
1644
1645 // 2. Calculate time interval
1646
1647 Kernel::TimeInterval deltaT;
1648
1649 // No filter
1650 if (n >= static_cast<int>(m_values.size()) || (n == static_cast<int>(m_values.size()) - 1 && m_values.size() == 1)) {
1651 // Out of bound
1652 ;
1653 } else if (n == static_cast<int>(m_values.size()) - 1) {
1654 // Last one by making up an end time.
1655 DateAndTime endTime = getFakeEndTime();
1656
1657 deltaT = Kernel::TimeInterval(m_values.rbegin()->time(), endTime);
1658 } else {
1659 // Regular
1660 DateAndTime startT = m_values[static_cast<std::size_t>(n)].time();
1661 DateAndTime endT = m_values[static_cast<std::size_t>(n) + 1].time();
1662 TimeInterval dt(startT, endT);
1663 deltaT = dt;
1664 }
1665
1666 return deltaT;
1667}
1668
1669template <typename TYPE> Types::Core::DateAndTime TimeSeriesProperty<TYPE>::getFakeEndTime() const {
1670 sortIfNecessary();
1671
1672 // the last time is the last thing known
1673 const auto ultimate = m_values.rbegin()->time();
1674
1675 // go backwards from the time before it that is different
1676 int counter = 0;
1677 while (DateAndTime::secondsFromDuration(ultimate - (m_values.rbegin() + counter)->time()) == 0.) {
1678 counter += 1;
1679 }
1680
1681 // get the last time that is different
1682 time_duration lastDuration = m_values.rbegin()->time() - (m_values.rbegin() + counter)->time();
1683
1684 // the last duration is equal to the previous, non-zero, duration
1685 return m_values.rbegin()->time() + lastDuration;
1686}
1687
1688//-----------------------------------------------------------------------------------------------
1694template <typename TYPE> TYPE TimeSeriesProperty<TYPE>::nthValue(int n) const {
1695
1696 // 1. Throw error if property is empty
1697 if (m_values.empty()) {
1698 const std::string error("nthValue(): TimeSeriesProperty '" + name() + "' is empty");
1699 g_log.debug(error);
1700 throw std::runtime_error(error);
1701 }
1702
1703 // 2. Sort and apply filter
1704 sortIfNecessary();
1705
1706 TYPE nthValue;
1707
1708 // 3. Situation 1: No filter
1709 if (static_cast<size_t>(n) < m_values.size()) {
1710 const auto entry = m_values[static_cast<std::size_t>(n)];
1711 nthValue = entry.value();
1712 } else {
1713 const auto entry = m_values[static_cast<std::size_t>(m_size) - 1];
1714 nthValue = entry.value();
1715 }
1716
1717 return nthValue;
1718}
1719
1725template <typename TYPE> Types::Core::DateAndTime TimeSeriesProperty<TYPE>::nthTime(int n) const {
1726 sortIfNecessary();
1727
1728 if (m_values.empty()) {
1729 const std::string error("nthTime(): TimeSeriesProperty '" + name() + "' is empty");
1730 g_log.debug(error);
1731 throw std::runtime_error(error);
1732 }
1733
1734 if (n < 0 || n >= static_cast<int>(m_values.size()))
1735 n = static_cast<int>(m_values.size()) - 1;
1736
1737 return m_values[static_cast<size_t>(n)].time();
1738}
1739
1743template <typename TYPE> void TimeSeriesProperty<TYPE>::countSize() const {
1744 // 1. Not filter
1745 m_size = int(m_values.size());
1746}
1747
1752template <typename TYPE> bool TimeSeriesProperty<TYPE>::isTimeString(const std::string &str) {
1753 static const boost::regex re("^[0-9]{4}.[0-9]{2}.[0-9]{2}.[0-9]{2}.[0-9]{2}.[0-9]{2}");
1754 return boost::regex_search(str.begin(), str.end(), re);
1755}
1756
1761template <typename TYPE> std::string TimeSeriesProperty<TYPE>::isValid() const { return ""; }
1762
1763/*
1764 * A TimeSeriesProperty never has a default, so return empty string
1765 * @returns Empty string as no defaults can be provided
1766 */
1767template <typename TYPE> std::string TimeSeriesProperty<TYPE>::getDefault() const {
1768 return ""; // No defaults can be provided=empty string
1769}
1770
1774template <typename TYPE> bool TimeSeriesProperty<TYPE>::isDefault() const { return false; }
1775
1783template <typename TYPE>
1785 // Start with statistics that are not time-weighted
1786 TimeSeriesPropertyStatistics out(Mantid::Kernel::getStatistics(this->filteredValuesAsVector(roi)));
1787 out.duration = this->durationInSeconds(roi);
1788 if (out.standard_deviation == 0.) {
1789 // if the simple std-dev is zero, just copy the simple mean
1790 out.time_mean = out.mean;
1791 out.time_standard_deviation = 0.;
1792 } else {
1793 // follow with time-weighted statistics
1794 auto avAndDev = this->timeAverageValueAndStdDev(roi);
1795 out.time_mean = avAndDev.first;
1796 out.time_standard_deviation = avAndDev.second;
1797 }
1798 return out;
1799}
1800
1801template <>
1803 // statistics of a string property doesn't make sense
1805 out.setAllToNan();
1806
1807 return out;
1808}
1809
1815template <typename TYPE>
1817 using namespace Kernel::Math;
1818 double singleValue = 0;
1819 switch (selection) {
1820 case FirstValue:
1821 if (roi && !roi->useAll())
1822 singleValue = double(this->firstValue(*roi));
1823 else
1824 singleValue = double(this->nthValue(0));
1825 break;
1826 case LastValue:
1827 if (roi && !roi->useAll())
1828 singleValue = double(this->lastValue(*roi));
1829 else
1830 singleValue = double(this->nthValue(this->size() - 1));
1831 break;
1832 case Minimum:
1833 singleValue = static_cast<double>(this->getStatistics(roi).minimum);
1834 break;
1835 case Maximum:
1836 singleValue = static_cast<double>(this->getStatistics(roi).maximum);
1837 break;
1838 case Mean:
1839 singleValue = this->getStatistics(roi).mean;
1840 break;
1841 case Median:
1842 singleValue = this->getStatistics(roi).median;
1843 break;
1844 case TimeAveragedMean:
1845 singleValue = this->getStatistics(roi).time_mean;
1846 break;
1847 case StdDev:
1848 singleValue = this->getStatistics(roi).standard_deviation;
1849 break;
1850 case TimeAverageStdDev:
1851 singleValue = this->getStatistics(roi).time_standard_deviation;
1852 break;
1853 default:
1854 throw std::invalid_argument("extractStatistic - Unknown statistic type: " + boost::lexical_cast<std::string>(this));
1855 };
1856 return singleValue;
1857}
1858
1862template <>
1864 UNUSED_ARG(selection);
1865 UNUSED_ARG(roi);
1866 throw Exception::NotImplementedError("TimeSeriesProperty::"
1867 "extractStatistic is not "
1868 "implemented for string properties");
1869}
1870
1871/*
1872 * Detects whether there are duplicated entries (of time) in property
1873 * If there is any, keep one of them
1874 */
1875template <typename TYPE> void TimeSeriesProperty<TYPE>::eliminateDuplicates() {
1876 // ensure that the values are sorted
1877 sortIfNecessary();
1878
1879 // cache the original size so the number removed can be reported
1880 const auto origSize{m_size};
1881
1882 // remove the first n-repeats
1883 // taken from
1884 // https://stackoverflow.com/questions/21060636/using-stdunique-and-vector-erase-to-remove-all-but-last-occurrence-of-duplicat
1885 auto it = std::unique(m_values.rbegin(), m_values.rend(),
1886 [](const auto &a, const auto &b) { return a.time() == b.time(); });
1887 m_values.erase(m_values.begin(), it.base());
1888
1889 // update m_size
1890 countSize();
1891
1892 // log how many values were removed
1893 const auto numremoved = origSize - m_size;
1894 if (numremoved > 0)
1895 g_log.notice() << "Log \"" << this->name() << "\" has " << numremoved << " entries removed due to duplicated time. "
1896 << "\n";
1897}
1898
1899/*
1900 * Print the content to string
1901 */
1902template <typename TYPE> std::string TimeSeriesProperty<TYPE>::toString() const {
1903 std::stringstream ss;
1904 for (size_t i = 0; i < m_values.size(); ++i)
1905 ss << m_values[i].time() << "\t\t" << m_values[i].value() << "\n";
1906
1907 return ss.str();
1908}
1909
1910//-------------------------------------------------------------------------
1911// Private methods
1912//-------------------------------------------------------------------------
1913
1914//----------------------------------------------------------------------------------
1915/*
1916 * Sort vector mP and set the flag. Only sorts if the values are not already
1917 * sorted.
1918 */
1919template <typename TYPE> void TimeSeriesProperty<TYPE>::sortIfNecessary() const {
1920 if (m_propSortedFlag == TimeSeriesSortStatus::TSUNKNOWN) {
1921 bool sorted = is_sorted(m_values.begin(), m_values.end());
1922 if (sorted)
1923 m_propSortedFlag = TimeSeriesSortStatus::TSSORTED;
1924 else
1925 m_propSortedFlag = TimeSeriesSortStatus::TSUNSORTED;
1926 }
1927
1928 if (m_propSortedFlag == TimeSeriesSortStatus::TSUNSORTED) {
1929 g_log.information() << "TimeSeriesProperty \"" << this->name()
1930 << "\" is not sorted. Sorting is operated on it. \n";
1931 std::stable_sort(m_values.begin(), m_values.end());
1932 m_propSortedFlag = TimeSeriesSortStatus::TSSORTED;
1933 }
1934}
1935
1942template <typename TYPE> int TimeSeriesProperty<TYPE>::findIndex(Types::Core::DateAndTime t) const {
1943 // 0. Return with an empty container
1944 if (m_values.empty())
1945 return 0;
1946
1947 // 1. Sort
1948 sortIfNecessary();
1949
1950 // 2. Extreme value
1951 if (t <= m_values[0].time()) {
1952 return -1;
1953 } else if (t >= m_values.back().time()) {
1954 return (int(m_values.size()));
1955 }
1956
1957 // 3. Find by lower_bound()
1958 typename std::vector<TimeValueUnit<TYPE>>::const_iterator fid;
1959 TimeValueUnit<TYPE> temp(t, m_values[0].value());
1960 fid = std::lower_bound(m_values.begin(), m_values.end(), temp);
1961
1962 int newindex = int(fid - m_values.begin());
1963 if (fid->time() > t)
1964 newindex--;
1965
1966 return newindex;
1967}
1968
1975template <typename TYPE>
1976int TimeSeriesProperty<TYPE>::upperBound(Types::Core::DateAndTime t, int istart, int iend) const {
1977 // 0. Check validity
1978 if (istart < 0) {
1979 throw std::invalid_argument("Start Index cannot be less than 0");
1980 }
1981 if (iend >= static_cast<int>(m_values.size())) {
1982 throw std::invalid_argument("End Index cannot exceed the boundary");
1983 }
1984 if (istart > iend) {
1985 throw std::invalid_argument("Start index cannot be greater than end index");
1986 }
1987
1988 // 1. Return instantly if it is out of boundary
1989 if (t < (m_values.begin() + istart)->time()) {
1990 return -1;
1991 }
1992 if (t > (m_values.begin() + iend)->time()) {
1993 return static_cast<int>(m_values.size());
1994 }
1995
1996 // 2. Sort
1997 sortIfNecessary();
1998
1999 // 3. Construct the pair for comparison and do lower_bound()
2000 const TimeValueUnit<TYPE> temppair(t, m_values[0].value());
2001 const auto first = m_values.cbegin() + istart;
2002 const auto last = m_values.cbegin() + iend + 1;
2003 const auto iter = std::lower_bound(first, last, temppair);
2004
2005 // 4. Calculate return value
2006 if (iter == last)
2007 throw std::runtime_error("Cannot find data");
2008 return static_cast<int>(std::distance(m_values.cbegin(), iter));
2009}
2010
2018template <typename TYPE> std::string TimeSeriesProperty<TYPE>::setValueFromProperty(const Property &right) {
2019 auto prop = dynamic_cast<const TimeSeriesProperty<TYPE> *>(&right);
2020 if (!prop) {
2021 return "Could not set value: properties have different type.";
2022 }
2023 m_values = prop->m_values;
2024 m_size = prop->m_size;
2025 m_propSortedFlag = prop->m_propSortedFlag;
2026 // m_filter = prop->m_filter;
2027 // m_filterQuickRef = prop->m_filterQuickRef;
2028 // m_filterApplied = prop->m_filterApplied;
2029 return "";
2030}
2031
2032//----------------------------------------------------------------------------------------------
2034template <typename TYPE> void TimeSeriesProperty<TYPE>::saveTimeVector(Nexus::File *file) {
2035 std::vector<DateAndTime> times = this->timesAsVector();
2036 const DateAndTime &start = times.front();
2037 std::vector<double> timeSec(times.size());
2038 for (size_t i = 0; i < times.size(); i++)
2039 timeSec[i] = static_cast<double>(times[i].totalNanoseconds() - start.totalNanoseconds()) * 1e-9;
2040 file->writeData("time", timeSec);
2041 file->openData("time");
2042 file->putAttr("start", start.toISO8601String());
2043 file->closeData();
2044}
2045
2046//----------------------------------------------------------------------------------------------
2048template <> void TimeSeriesProperty<std::string>::saveProperty(Nexus::File *file) {
2049 std::vector<std::string> values = this->valuesAsVector();
2050 if (values.empty())
2051 return;
2052 file->makeGroup(this->name(), "NXlog", true);
2053
2054 // Find the max length of any string
2055 auto max_it = std::max_element(values.begin(), values.end(),
2056 [](const std::string &a, const std::string &b) { return a.size() < b.size(); });
2057 // Increment by 1 to have the 0 terminator
2058 size_t maxlen = max_it->size() + 1;
2059 // Copy into one array
2060 std::vector<char> strs(values.size() * maxlen);
2061 size_t index = 0;
2062 for (const auto &prop : values) {
2063 std::copy(prop.begin(), prop.end(), &strs[index]);
2064 index += maxlen;
2065 }
2066
2067 const Mantid::Nexus::DimVector dims{values.size(), maxlen};
2068 file->makeData("value", NXnumtype::CHAR, dims, true);
2069 file->putData(strs.data());
2070 file->closeData();
2071 saveTimeVector(file);
2072 file->closeGroup();
2073}
2074
2081template <> void TimeSeriesProperty<bool>::saveProperty(Nexus::File *file) {
2082 std::vector<bool> value = this->valuesAsVector();
2083 if (value.empty())
2084 return;
2085 std::vector<uint8_t> asUint(value.begin(), value.end());
2086 file->makeGroup(this->name(), "NXlog", true);
2087 file->writeData("value", asUint);
2088 file->putAttr("boolean", "1");
2089 saveTimeVector(file);
2090 file->closeGroup();
2091}
2092
2093template <typename TYPE> void TimeSeriesProperty<TYPE>::saveProperty(Nexus::File *file) {
2094 auto values = this->valuesAsVector();
2095 if (values.empty())
2096 return;
2097 file->makeGroup(this->name(), "NXlog", 1);
2098 file->writeData("value", values);
2099 file->openData("value");
2100 file->putAttr("units", this->units());
2101 file->closeData();
2102 saveTimeVector(file);
2103 file->closeGroup();
2104}
2105
2110template <typename TYPE> Json::Value TimeSeriesProperty<TYPE>::valueAsJson() const { return Json::Value(value()); }
2111
2122template <typename TYPE>
2123void TimeSeriesProperty<TYPE>::histogramData(const Types::Core::DateAndTime &tMin, const Types::Core::DateAndTime &tMax,
2124 std::vector<double> &counts) const {
2125
2126 size_t nPoints = counts.size();
2127 if (nPoints == 0)
2128 return; // nothing to do
2129
2130 auto t0 = static_cast<double>(tMin.totalNanoseconds());
2131 auto t1 = static_cast<double>(tMax.totalNanoseconds());
2132 if (t0 > t1)
2133 throw std::invalid_argument("invalid arguments for histogramData; tMax<tMin");
2134
2135 double dt = (t1 - t0) / static_cast<double>(nPoints);
2136
2137 for (auto const &ev : m_values) {
2138 auto time = static_cast<double>(ev.time().totalNanoseconds());
2139 if (time < t0 || time >= t1)
2140 continue;
2141 auto ind = static_cast<size_t>((time - t0) / dt);
2142 counts[ind] += static_cast<double>(ev.value());
2143 }
2144}
2145
2146template <>
2147void TimeSeriesProperty<std::string>::histogramData(const Types::Core::DateAndTime &tMin,
2148 const Types::Core::DateAndTime &tMax,
2149 std::vector<double> &counts) const {
2150 UNUSED_ARG(tMin);
2151 UNUSED_ARG(tMax);
2152 UNUSED_ARG(counts);
2153 throw std::runtime_error("histogramData is not implememnted for time series "
2154 "properties containing strings");
2155}
2156
2166template <typename TYPE> std::vector<TYPE> TimeSeriesProperty<TYPE>::filteredValuesAsVector(const TimeROI *roi) const {
2167 if (roi && !roi->useAll()) {
2168 this->sortIfNecessary();
2169 std::vector<TYPE> filteredValues;
2170 if (roi->firstTime() > this->m_values.back().time()) {
2171 // Since the ROI starts after everything, just return the last value in the log
2172 filteredValues.emplace_back(this->m_values.back().value());
2173 } else { // only use the values in the filter - this is very similar to FilteredTimeSeriesProperty::applyFilter
2174 // the index into the m_values array of the time, or -1 (before) or m_values.size() (after)
2175 std::size_t index_current_log{0};
2176 for (const auto &splitter : roi->toTimeIntervals()) {
2177 const auto endTime = splitter.stop();
2178
2179 // check if the splitter starts too early
2180 if (endTime < this->m_values[index_current_log].time()) {
2181 continue; // skip to the next splitter
2182 }
2183
2184 // cache values to reduce number of method calls
2185 const auto beginTime = splitter.start();
2186
2187 // find the first log that should be added
2188 if (this->m_values.back().time() < beginTime) {
2189 // skip directly to the end if the filter starts after the last log
2190 index_current_log = this->m_values.size() - 1;
2191 } else {
2192 // search for the right starting point
2193 while ((this->m_values[index_current_log].time() <= beginTime)) {
2194 if (index_current_log + 1 > this->m_values.size())
2195 break;
2196 index_current_log++;
2197 }
2198 // need to back up by one
2199 if (index_current_log > 0)
2200 index_current_log--;
2201 // go backwards more while times are equal to the one being started at
2202 while (index_current_log > 0 &&
2203 this->m_values[index_current_log].time() == this->m_values[index_current_log - 1].time()) {
2204 index_current_log--;
2205 }
2206 }
2207
2208 // add everything up to the end time
2209 for (; index_current_log < this->m_values.size(); ++index_current_log) {
2210 if (this->m_values[index_current_log].time() >= endTime)
2211 break;
2212
2213 // the current value goes into the filter
2214 filteredValues.emplace_back(this->m_values[index_current_log].value());
2215 }
2216 // go back one so the next splitter can add a value
2217 if (index_current_log > 0)
2218 index_current_log--;
2219 }
2220 }
2221 return filteredValues;
2222 } else {
2223 return this->valuesAsVector();
2224 }
2225}
2226
2227template <typename TYPE> std::vector<TYPE> TimeSeriesProperty<TYPE>::filteredValuesAsVector() const {
2228 return this->valuesAsVector();
2229}
2230
2239template <typename TYPE> std::vector<TimeInterval> TimeSeriesProperty<TYPE>::getTimeIntervals() const {
2240 std::vector<TimeInterval> intervals;
2241 auto lastInterval = this->nthInterval(this->size() - 1);
2242 intervals.emplace_back(firstTime(), lastInterval.stop());
2243 return intervals;
2244}
2245
2247// -------------------------- Macro to instantiation concrete types
2248// --------------------------------
2249#define INSTANTIATE(TYPE) template class TimeSeriesProperty<TYPE>;
2250
2251// -------------------------- Concrete instantiation
2252// -----------------------------------------------
2253INSTANTIATE(int32_t)
2254INSTANTIATE(int64_t)
2255INSTANTIATE(uint32_t)
2256INSTANTIATE(uint64_t)
2257INSTANTIATE(float)
2258INSTANTIATE(double)
2259INSTANTIATE(std::string)
2260INSTANTIATE(bool)
2261
2262
2263
2264} // namespace Kernel
2265} // namespace Mantid
std::string name
Definition Run.cpp:60
const std::vector< double > & rhs
size_t m_size
Maximum size of the store.
double value
The value of the point.
Definition FitMW.cpp:51
double error
std::map< DeltaEMode::Type, std::string > index
#define INSTANTIATE(TYPE)
double right
boost::python::list getTimeIntervals(const TimeROI &self)
Definition TimeROI.cpp:19
#define UNUSED_ARG(x)
Function arguments are sometimes unused in certain implmentations but are required for documentation ...
Definition System.h:44
Marks code as not implemented yet.
Definition Exception.h:138
void debug(const std::string &msg)
Logs at debug level.
Definition Logger.cpp:145
void notice(const std::string &msg)
Logs at notice level.
Definition Logger.cpp:126
void warning(const std::string &msg)
Logs at warning level.
Definition Logger.cpp:117
void information(const std::string &msg)
Logs at information level.
Definition Logger.cpp:136
Base class for properties.
Definition Property.h:94
const std::string & name() const
Get the property's name.
Definition Property.cpp:63
Represents a time interval.
Definition DateAndTime.h:25
const Types::Core::DateAndTime & start() const
Beginning of the interval.
Definition DateAndTime.h:34
const Types::Core::DateAndTime & stop() const
End of the interval.
Definition DateAndTime.h:36
TimeROI : Object that holds information about when the time measurement was active.
Definition TimeROI.h:18
double durationInSeconds() const
Duration of the whole TimeROI.
Definition TimeROI.cpp:624
void addMask(const std::string &startTime, const std::string &stopTime)
Definition TimeROI.cpp:151
const std::vector< Kernel::TimeInterval > toTimeIntervals() const
This method is to lend itself to helping with transition.
Definition TimeROI.cpp:557
const std::vector< Types::Core::DateAndTime > & getAllTimes() const
Definition TimeROI.h:49
Types::Core::DateAndTime firstTime() const
Definition TimeROI.cpp:328
void addROI(const std::string &startTime, const std::string &stopTime)
Definition TimeROI.cpp:76
bool useNone() const
TimeROI selects no time to be used as all is invalid.
Definition TimeROI.cpp:695
bool useAll() const
TimeROI selects all time to be used.
Definition TimeROI.cpp:693
void update_intersection(const TimeROI &other)
Updates the TimeROI values with the intersection with another TimeROI.
Definition TimeROI.cpp:503
Types::Core::DateAndTime lastTime() const
Definition TimeROI.cpp:335
static const TimeROI USE_NONE
Constant for TimeROI where no time is used.
Definition TimeROI.h:22
A specialised Property class for holding a series of time-value pairs.
void replaceValues(const std::vector< Types::Core::DateAndTime > &times, const std::vector< TYPE > &values)
Replaces the time series with new values time series values.
TimeSeriesProperty< TYPE > * clone() const override
"Virtual" copy constructor
std::string getDefault() const override
Returns the default value.
std::map< Types::Core::DateAndTime, TYPE > valueAsMap() const
Return the time series as a C++ map<DateAndTime, TYPE>
size_t getMemorySize() const override
Return the memory used by the property, in bytes.
Json::Value valueAsJson() const override
std::string toString() const
Stringize the property.
int size() const override
Returns the number of values at UNIQUE time intervals in the time series.
TYPE minValue() const
Returns the minimum value found in the series.
TYPE maxValue() const
Returns the maximum value found in the series.
std::vector< TYPE > valuesAsVector() const
Return the time series's values (unfiltered) as a vector<TYPE>
Types::Core::DateAndTime firstTime() const
Returns the first time regardless of filter.
double durationInSeconds(const Kernel::TimeROI *roi=nullptr) const
Returns the duration of the time series, possibly restricted by a TimeROI object.
std::pair< double, double > timeAverageValueAndStdDev(const Kernel::TimeROI *timeRoi=nullptr) const override
Returns the calculated time weighted mean and standard deviation values.
void eliminateDuplicates()
Detects whether there are duplicated entries (of time) in property & eliminates them.
std::vector< double > timesAsVectorSeconds() const
Return the series as list of times, where the time is the number of seconds since the start.
std::multimap< Types::Core::DateAndTime, TYPE > valueAsMultiMap() const
Return the time series as a correct C++ multimap<DateAndTime, TYPE>.
TimeSeriesProperty & operator+=(Property const *right) override
Add the value of another property.
TimeSeriesProperty< TYPE > & merge(Property *rhs) override
Merge the given property with this one.
std::string value() const override
Get the time series property as a string of 'time value'.
void createFilteredData(const TimeROI &timeROI, std::vector< TimeValueUnit< TYPE > > &filteredData) const
Fill in the supplied vector of time series data according to the input TimeROI.
Property * cloneInTimeROI(const TimeROI &timeROI) const override
Create a partial copy according to TimeROI.
void makeFilterByValue(std::vector< SplittingInterval > &split, double min, double max, double TimeTolerance=0.0, bool centre=false) const override
Fill a SplittingIntervalVec that will filter the events by matching.
bool isDefault() const override
Returns if the value is at the default.
void sortIfNecessary() const
Sort the property into increasing times, if not already sorted.
double mean() const
Returns the mean value found in the series.
TimeSeriesPropertyStatistics getStatistics(const Kernel::TimeROI *roi=nullptr) const override
Return a TimeSeriesPropertyStatistics object.
int findIndex(Types::Core::DateAndTime t) const
Find the index of the entry of time t in the mP vector (sorted)
TYPE firstValue() const
Returns the first value regardless of filter.
void saveTimeVector(Nexus::File *file)
Saves the time vector has time + start attribute.
virtual bool operator!=(const TimeSeriesProperty< TYPE > &right) const
Deep comparison (not equal).
void histogramData(const Types::Core::DateAndTime &tMin, const Types::Core::DateAndTime &tMax, std::vector< double > &counts) const
generate constant time-step histogram from the property values
void expandFilterToRange(std::vector< SplittingInterval > &split, double min, double max, const TimeInterval &range) const override
Make sure an existing filter covers the full time range given.
~TimeSeriesProperty() override
Virtual destructor.
double timeAverageValue(const TimeROI *timeRoi=nullptr) const override
Returns the calculated time weighted average value.
std::string setValue(const std::string &) override
Set a property from a string.
virtual std::vector< Types::Core::DateAndTime > filteredTimesAsVector() const
void addValue(const Types::Core::DateAndTime &time, const TYPE &value)
Add a value to the map using a DateAndTime object.
int m_size
The number of values (or time intervals) in the time series.
TimeSeriesProperty(const std::string &name)
Constructor.
void removeDataOutsideTimeROI(const TimeROI &timeRoi) override
Remove time values outside of TimeROI regions each defined as [roi_begin,roi_end].
void clear() override
Deletes the series of values in the property.
Types::Core::DateAndTime lastTime() const
Returns the last time.
void setName(const std::string &name)
Set name of property.
TYPE lastValue() const
Returns the last value.
double extractStatistic(Math::StatisticType selection, const TimeROI *roi=nullptr) const override
Calculate a particular statistical quantity from the values of the time series.
std::map< Types::Core::DateAndTime, TYPE > valueAsCorrectMap() const
Return the time series as a correct C++ map<DateAndTime, TYPE>.
std::string setValueFromJson(const Json::Value &) override
Set a property from a string.
std::unique_ptr< TimeSeriesProperty< double > > getDerivative() const
Return time series property, containing time derivative of current property.
std::vector< TimeValueUnit< TYPE > > m_values
Holds the time series data.
void create(const Types::Core::DateAndTime &start_time, const std::vector< double > &time_sec, const std::vector< TYPE > &new_values)
Clears and creates a TimeSeriesProperty from these parameters.
Types::Core::DateAndTime getFakeEndTime() const
Returns an end time that will have the same spacing to the right of the last value as the last non-ze...
virtual bool operator==(const TimeSeriesProperty< TYPE > &right) const
Deep comparison.
double averageValueInFilter(const std::vector< TimeInterval > &filter) const
Calculate the time-weighted average of a property in a filtered range.
virtual std::vector< Mantid::Kernel::TimeInterval > getTimeIntervals() const
If filtering by log, get the time intervals for splitting.
void clearOutdated() override
Deletes all but the 'last entry' in the property.
virtual TYPE nthValue(int n) const
Returns n-th value of n-th interval in an incredibly inefficient way.
std::string setDataItem(const std::shared_ptr< DataItem > &) override
Set a property from a DataItem.
std::string isValid() const override
This doesn't check anything -we assume these are always valid.
std::vector< std::string > time_tValue() const
New method to return time series value pairs as std::vector<std::string>
std::vector< Types::Core::DateAndTime > timesAsVector() const override
Return the time series's times as a vector<DateAndTime>
std::pair< double, double > averageAndStdDevInFilter(const std::vector< TimeInterval > &intervals) const
Calculate the time-weighted average and std-deviation of a property in a filtered range.
virtual Types::Core::DateAndTime nthTime(int n) const
Returns n-th time. NOTE: Complexity is order(n)! regardless of filter.
TYPE getSingleValue(const Types::Core::DateAndTime &t) const
Returns the value at a particular time.
int upperBound(Types::Core::DateAndTime t, int istart, int iend) const
Find the upper_bound of time t in container.
virtual TimeInterval nthInterval(int n) const
Returns n-th valid time interval, in a very inefficient way.
void addValues(const std::vector< Types::Core::DateAndTime > &times, const std::vector< TYPE > &values)
Adds vectors of values to the map.
int realSize() const override
Returns the real size of the time series property map:
static bool isTimeString(const std::string &str)
Check if str has the right time format.
Property * cloneWithTimeShift(const double timeShift) const override
"Virtual" copy constructor with a time shift in seconds
std::string setValueFromProperty(const Property &right) override
Set a value from another property.
virtual std::vector< TYPE > filteredValuesAsVector() const
void saveProperty(Nexus::File *file) override
Class to hold unit value (DateAndTime, T)
static unsigned short constexpr CHAR
MatrixWorkspace_sptr MANTID_API_DLL operator+=(const MatrixWorkspace_sptr &lhs, const MatrixWorkspace_sptr &rhs)
Adds two workspaces.
void split(const int A, int &S, int &V)
Split a number into the sign and positive value.
Definition Acomp.cpp:42
StatisticType
Maps a "statistic" to a number.
Definition Statistics.h:18
Statistics getStatistics(const std::vector< TYPE > &data, const unsigned int flags=StatOptions::AllStats)
Return a statistics object for the given data set.
std::vector< SplittingInterval > SplittingIntervalVec
A typedef for splitting events according their pulse time.
Definition LogManager.h:27
MANTID_KERNEL_DLL bool operator==(const Mantid::Kernel::Property &lhs, const Mantid::Kernel::Property &rhs)
Compares this to another property for equality.
Definition Property.cpp:246
std::vector< dimsize_t > DimVector
Helper class which provides the Collimation Length for SANS instruments.
constexpr double EMPTY_DBL() noexcept
Returns what we consider an "empty" double within a property.
Definition EmptyValues.h:42
STL namespace.
Simple struct to store statistics.
Definition Statistics.h:35
double mean
Mean value.
Definition Statistics.h:41
double median
Median value.
Definition Statistics.h:43
double minimum
Minimum value.
Definition Statistics.h:37
double maximum
Maximum value.
Definition Statistics.h:39
double standard_deviation
standard_deviation of the values
Definition Statistics.h:45
Struct holding some useful statistics for a TimeSeriesProperty.
double time_standard_deviation
time weighted standard deviation
double standard_deviation
standard_deviation of the values