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 const std::size_t num_values = values.size();
35 assert(num_values > 1);
36 const auto &first_value = values.front().value();
37 for (std::size_t i = 1; i < num_values; ++i) {
38 if (first_value != values[i].value())
39 return false;
40 }
41 return true;
42}
43} // namespace
44
49template <typename TYPE>
51 : Property(name, typeid(std::vector<TimeValueUnit<TYPE>>)), m_values(), m_size(), m_propSortedFlag() {}
52
59template <typename TYPE>
61 const std::vector<Types::Core::DateAndTime> &times,
62 const std::vector<TYPE> &values)
64 addValues(times, values);
65}
66
68template <typename TYPE> TimeSeriesProperty<TYPE>::~TimeSeriesProperty() = default;
69
74 return new TimeSeriesProperty<TYPE>(*this);
75}
76
81template <typename TYPE>
83 : Property(*p), m_values(), m_size(), m_propSortedFlag() {}
84
92template <typename TYPE> Property *TimeSeriesProperty<TYPE>::cloneInTimeROI(const TimeROI &timeROI) const {
93 auto filteredTS = new TimeSeriesProperty<TYPE>(this);
94
95 createFilteredData(timeROI, filteredTS->m_values);
96
97 filteredTS->m_size = static_cast<int>(filteredTS->m_values.size());
98
99 return filteredTS;
100}
101
106template <typename TYPE> Property *TimeSeriesProperty<TYPE>::cloneWithTimeShift(const double timeShift) const {
107 auto timeSeriesProperty = this->clone();
108 auto values = timeSeriesProperty->valuesAsVector();
109 auto times = timeSeriesProperty->timesAsVector();
110 // Shift the time
111 for (auto it = times.begin(); it != times.end(); ++it) {
112 // There is a known issue which can cause cloneWithTimeShift to be called
113 // with a large (~9e+9 s) shift. Actual shifting is capped to be ~4.6e+19
114 // seconds in DateAndTime::operator+=
115 (*it) += timeShift;
116 }
117 timeSeriesProperty->clear();
118 timeSeriesProperty->addValues(times, values);
119 return timeSeriesProperty;
120}
121
129template <typename TYPE> std::unique_ptr<TimeSeriesProperty<double>> TimeSeriesProperty<TYPE>::getDerivative() const {
130
131 if (this->m_values.size() < 2) {
132 throw std::runtime_error("Derivative is not defined for a time-series "
133 "property with less then two values");
134 }
135
136 this->sortIfNecessary();
137 auto it = this->m_values.begin();
138 int64_t t0 = it->time().totalNanoseconds();
139 TYPE v0 = it->value();
141 it++;
142 auto timeSeriesDeriv = std::make_unique<TimeSeriesProperty<double>>(this->name() + "_derivative");
143 timeSeriesDeriv->reserve(this->m_values.size() - 1);
144 for (; it != m_values.end(); it++) {
145 TYPE v1 = it->value();
146 int64_t t1 = it->time().totalNanoseconds();
147 if (t1 != t0) {
148 double deriv = 1.e+9 * (double(v1 - v0) / double(t1 - t0));
149 auto tm = static_cast<int64_t>((t1 + t0) / 2);
150 timeSeriesDeriv->addValue(Types::Core::DateAndTime(tm), deriv);
151 }
152 t0 = t1;
153 v0 = v1;
154 }
155 return timeSeriesDeriv;
158template <> std::unique_ptr<TimeSeriesProperty<double>> TimeSeriesProperty<std::string>::getDerivative() const {
159 throw std::runtime_error("Time series property derivative is not defined for strings");
161
165template <typename TYPE> size_t TimeSeriesProperty<TYPE>::getMemorySize() const {
166 // Rough estimate
167 return m_values.size() * (sizeof(TYPE) + sizeof(DateAndTime));
169
177}
185 auto const *rhs = dynamic_cast<TimeSeriesProperty<TYPE> const *>(right);
186
187 if (rhs) {
188 if (this->operator!=(*rhs)) {
189 m_values.insert(m_values.end(), rhs->m_values.begin(), rhs->m_values.end());
190 m_propSortedFlag = TimeSeriesSortStatus::TSUNKNOWN;
191 } else {
192 // Do nothing if appending yourself to yourself. The net result would be
193 // the same anyway
195 }
196
197 // Count the REAL size.
198 m_size = static_cast<int>(m_values.size());
199
200 } else
201 g_log.warning() << "TimeSeriesProperty " << this->name()
202 << " could not be added to another property of the same "
203 "name but incompatible type.\n";
204
205 return *this;
206}
213template <typename TYPE> bool TimeSeriesProperty<TYPE>::operator==(const TimeSeriesProperty<TYPE> &right) const {
214 sortIfNecessary();
215
216 if (this->name() != right.name()) // should this be done?
218 return false;
220
221 if (this->m_size != right.m_size) {
222 return false;
223 }
225 if (this->realSize() != right.realSize()) {
226 return false;
227 } else {
228 const std::vector<DateAndTime> lhsTimes = this->timesAsVector();
229 const std::vector<DateAndTime> rhsTimes = right.timesAsVector();
230 if (!std::equal(lhsTimes.begin(), lhsTimes.end(), rhsTimes.begin())) {
231 return false;
233
234 const std::vector<TYPE> lhsValues = this->valuesAsVector();
235 const std::vector<TYPE> rhsValues = right.valuesAsVector();
236 if (!std::equal(lhsValues.begin(), lhsValues.end(), rhsValues.begin())) {
237 return false;
239 }
241 return true;
242}
249template <typename TYPE> bool TimeSeriesProperty<TYPE>::operator==(const Property &right) const {
250 auto rhs_tsp = dynamic_cast<const TimeSeriesProperty<TYPE> *>(&right);
251 if (!rhs_tsp)
252 return false;
253 return this->operator==(*rhs_tsp);
255
261template <typename TYPE> bool TimeSeriesProperty<TYPE>::operator!=(const TimeSeriesProperty<TYPE> &right) const {
262 return !(*this == right);
264
270template <typename TYPE> bool TimeSeriesProperty<TYPE>::operator!=(const Property &right) const {
271 return !(*this == right);
273
277template <typename TYPE> void TimeSeriesProperty<TYPE>::setName(const std::string &name) { m_name = name; }
278
287template <typename TYPE>
289 std::vector<TimeValueUnit<TYPE>> &filteredData) const {
290 filteredData.clear();
291
292 // Expediently treat a few special cases
293
294 // Nothing to copy
295 if (m_values.empty()) {
296 return;
297 }
298
299 // Copy the only value
300 if (m_values.size() == 1) {
301 filteredData.push_back(m_values.front());
302 return;
303 }
305 // Copy the first value only, if all values are the same
306 // Exclude "proton_charge" logs from consideration, because in a real measurement those values can't be the same,
307 // Removing some of them just because they are equal will cause wrong total proton charge results.
308 if (allValuesAreSame(m_values) && this->name() != "proton_charge") {
309 filteredData.push_back(m_values.front());
310 return;
312
313 // Copy everything
314 if (timeROI.useAll()) {
315 std::copy(m_values.cbegin(), m_values.cend(), std::back_inserter(filteredData));
316 return;
318
319 // Copy the first value only
320 if (timeROI.useNone()) {
321 filteredData.push_back(m_values.front());
322 return;
324
325 // Now treat the general case
326
327 // Get all ROI time boundaries. Every other value is start/stop of an ROI "use" region.
328 const std::vector<Types::Core::DateAndTime> &roiTimes = timeROI.getAllTimes();
329 auto itROI = roiTimes.cbegin();
330 const auto itROIEnd = roiTimes.cend();
331
332 auto itValue = m_values.cbegin();
333 const auto itValueEnd = m_values.cend();
334 auto itLastValueUsed = itValue; // last value used up to the moment
335
336 while (itROI != itROIEnd && itValue != itValueEnd) {
337 // Try fast-forwarding the current ROI "use" region towards the current time value. Note, the current value might
338 // be in an ROI "ignore" region together with one or more following values.
339 while (std::distance(itROI, itROIEnd) > 2 && *(std::next(itROI, 2)) <= itValue->time())
340 std::advance(itROI, 2);
341 // Try finding the first value equal or past the beginning of the current ROI "use" region.
342 itValue = std::lower_bound(itValue, itValueEnd, *itROI,
343 [](const auto &value, const auto &roi_time) { return value.time() < roi_time; });
344 // Calculate a [begin,end) range for the values to use
345 auto itBeginUseValue = itValue;
346 auto itEndUseValue = itValue;
347 // If there are no values past the current ROI "use" region, get the previous value
348 if (itValue == itValueEnd) {
349 itBeginUseValue =
350 std::prev(itValue); // std::prev is safe here, because "m_values is empty" case has already been treated above
351 itEndUseValue = itValueEnd;
353 // If the value is inside the current ROI "use" region, look for other values in the same ROI "use" region
354 else if (itValue->time() <= *(std::next(itROI))) {
355 // First, try including a value immediately preceding the first value in the ROI "use" region.
356 itBeginUseValue = itValue == m_values.begin() ? itValue : std::prev(itValue);
357 // Now try finding the first value past the end of the current ROI "use" region.
358 while (itValue != itValueEnd && itValue->time() <= *(std::next(itROI)))
359 itValue++;
360 // Include the current value, therefore, advance itEndUseValue, because std::copy works as [begin,end).
361 itEndUseValue = itValue == itValueEnd ? itValue : std::next(itValue);
362 }
363 // If we are at the last ROI "use" region or the value is not past the beginning of the next ROI "use" region, keep
364 // it for the current ROI "use" region.
365 else if (std::distance(itROI, itROIEnd) == 2 ||
366 (std::distance(itROI, itROIEnd) > 2 && itValue->time() < *(std::next(itROI, 2)))) {
367 // Try including the value immediately preceding the current value
368 itBeginUseValue = itValue == m_values.begin() ? itValue : std::prev(itValue);
369 itEndUseValue = std::next(itValue);
371 // Do not use a value already copied for the previous ROI
372 if (!filteredData.empty()) {
373 itBeginUseValue = std::max(itBeginUseValue, std::next(itLastValueUsed));
374 }
376 // Copy all [begin,end) values and mark the last value copied
377 if (itBeginUseValue < itEndUseValue) {
378 std::copy(itBeginUseValue, itEndUseValue, std::back_inserter(filteredData));
379 itLastValueUsed = std::prev(itEndUseValue);
380 }
381
382 // Move to the next ROI "use" region
383 std::advance(itROI, 2);
384 }
385}
386
392template <typename TYPE> void TimeSeriesProperty<TYPE>::removeDataOutsideTimeROI(const TimeROI &timeROI) {
393 std::vector<TimeValueUnit<TYPE>> mp_copy;
394 createFilteredData(timeROI, mp_copy);
395
396 m_values.clear();
397 m_values = mp_copy;
398 mp_copy.clear();
399
400 m_size = static_cast<int>(m_values.size());
401}
402
403// The makeFilterByValue & expandFilterToRange methods generate a bunch of
404// warnings when the template type is the wider integer types
405// (when it's being assigned back to a double such as in a call to minValue or
406// firstValue)
407// However, in reality these methods are only used for TYPE=int or double (they
408// are only called from FilterByLogValue) so suppress the warnings
409#ifdef _WIN32
410#pragma warning(push)
411#pragma warning(disable : 4244)
412#pragma warning(disable : 4804) // This one comes about for TYPE=bool - again
413 // the method is never called for this type
414#endif
415#if defined(__GNUC__) && !(defined(__INTEL_COMPILER))
416#pragma GCC diagnostic ignored "-Wconversion"
417#endif
418
432template <typename TYPE>
433void TimeSeriesProperty<TYPE>::makeFilterByValue(std::vector<SplittingInterval> &split, double min, double max,
434 double TimeTolerance, bool centre) const {
435 const bool emptyMin = (min == EMPTY_DBL());
436 const bool emptyMax = (max == EMPTY_DBL());
437
438 if (!emptyMin && !emptyMax && max < min) {
439 std::stringstream ss;
440 ss << "TimeSeriesProperty::makeFilterByValue: 'max' argument must be "
441 "greater than 'min' "
442 << "(got min=" << min << " max=" << max << ")";
443 throw std::invalid_argument(ss.str());
444 }
445
446 // If min or max were unset ("empty") in the algorithm, set to the min or max
447 // value of the log
448 if (emptyMin)
449 min = static_cast<double>(minValue());
450 if (emptyMax)
451 max = static_cast<double>(maxValue());
452
453 // Make sure the splitter starts out empty
454 split.clear();
455
456 // Do nothing if the log is empty.
457 if (m_values.empty())
458 return;
459
460 // 1. Sort
461 sortIfNecessary();
462
463 // 2. Do the rest
464 bool lastGood(false);
465 time_duration tol = DateAndTime::durationFromSeconds(TimeTolerance);
466 int numgood = 0;
467 DateAndTime t;
468 DateAndTime start, stop;
469
470 for (size_t i = 0; i < m_values.size(); ++i) {
471 const DateAndTime lastGoodTime = t;
472 // The new entry
473 t = m_values[i].time();
474 TYPE val = m_values[i].value();
475
476 // A good value?
477 const bool isGood = ((val >= min) && (val <= max));
478 if (isGood)
479 numgood++;
480
481 if (isGood != lastGood) {
482 // We switched from bad to good or good to bad
483
484 if (isGood) {
485 // Start of a good section. Subtract tolerance from the time if
486 // boundaries are centred.
487 start = centre ? t - tol : t;
488 } else {
489 // End of the good section. Add tolerance to the LAST GOOD time if
490 // boundaries are centred.
491 // Otherwise, use the first 'bad' time.
492 stop = centre ? lastGoodTime + tol : t;
493 split.emplace_back(start, stop, 0);
494 // Reset the number of good ones, for next time
495 numgood = 0;
496 }
497 lastGood = isGood;
498 }
499 }
500
501 if (numgood > 0) {
502 // The log ended on "good" so we need to close it using the last time we
503 // found
504 stop = t + tol;
505 split.emplace_back(start, stop, 0);
506 }
507}
508
512template <>
513void TimeSeriesProperty<std::string>::makeFilterByValue(std::vector<SplittingInterval> & /*split*/, double /*min*/,
514 double /*max*/, double /*TimeTolerance*/,
515 bool /*centre*/) const {
516 throw Exception::NotImplementedError("TimeSeriesProperty::makeFilterByValue "
517 "is not implemented for string "
518 "properties");
519}
520
537template <typename TYPE>
538TimeROI TimeSeriesProperty<TYPE>::makeFilterByValue(double min, double max, bool expand,
539 const TimeInterval &expandRange, double TimeTolerance, bool centre,
540 const TimeROI *existingROI) const {
541 const bool emptyMin = (min == EMPTY_DBL());
542 const bool emptyMax = (max == EMPTY_DBL());
543
544 if (!emptyMin && !emptyMax && max < min) {
545 std::stringstream ss;
546 ss << "TimeSeriesProperty::makeFilterByValue: 'max' argument must be "
547 "greater than 'min' "
548 << "(got min=" << min << " max=" << max << ")";
549 throw std::invalid_argument(ss.str());
550 }
551
552 // If min or max were unset ("empty") in the algorithm, set to the min or max
553 // value of the log
554 if (emptyMin)
555 min = static_cast<double>(minValue());
556 if (emptyMax)
557 max = static_cast<double>(maxValue());
558
559 TimeROI newROI;
560
561 // Do nothing if the log is empty.
562 if (m_values.empty())
563 return newROI;
564
565 // 1. Sort
566 sortIfNecessary();
567
568 // 2. Do the rest
569 const time_duration tol = DateAndTime::durationFromSeconds(TimeTolerance);
570 DateAndTime stop_t;
571 DateAndTime start, stop;
572
573 bool isGood = false;
574 for (size_t i = 0; i < m_values.size(); ++i) {
575 TYPE val = m_values[i].value();
576
577 if ((val >= min) && (val <= max)) {
578 if (isGood) {
579 stop_t = m_values[i].time();
580 } else {
581 isGood = true;
582 stop_t = m_values[i].time();
583 start = centre ? m_values[i].time() - tol : m_values[i].time();
584 }
585 } else if (isGood) {
586 stop = centre ? stop_t + tol : m_values[i].time();
587 if (start < stop)
588 newROI.addROI(start, stop);
589 isGood = false;
590 }
591 }
592 if (isGood) {
593 stop = centre ? stop_t + tol : stop_t;
594 if (start < stop)
595 newROI.addROI(start, stop);
596 }
597
598 if (expand) {
599 if (expandRange.start() < firstTime()) {
600 auto val = static_cast<double>(firstValue());
601 if ((val >= min) && (val <= max)) {
602 newROI.addROI(expandRange.start(), firstTime());
603 }
604 }
605 if (expandRange.stop() > lastTime()) {
606 auto val = static_cast<double>(lastValue());
607 if ((val >= min) && (val <= max)) {
608 newROI.addROI(lastTime(), expandRange.stop());
609 }
610 }
611 }
612
613 // If the TimeROI is empty there are no values inside the filter
614 // so we should return USE_NONE
615 if (newROI.useAll()) {
616 return TimeROI::USE_NONE;
617 }
618
619 if (existingROI != nullptr && !existingROI->useAll()) {
620 newROI.update_intersection(*existingROI);
621 }
622 return newROI;
623}
624
628template <>
629TimeROI TimeSeriesProperty<std::string>::makeFilterByValue(double /*min*/, double /*max*/, bool /*expand*/,
630 const TimeInterval & /*expandRange*/,
631 double /*TimeTolerance*/, bool /*centre*/,
632 const TimeROI * /*existingROI*/) const {
633 throw Exception::NotImplementedError("TimeSeriesProperty::makeFilterByValue "
634 "is not implemented for string "
635 "properties");
636}
637
648template <typename TYPE>
649void TimeSeriesProperty<TYPE>::expandFilterToRange(std::vector<SplittingInterval> &split, double min, double max,
650 const TimeInterval &range) const {
651 const bool emptyMin = (min == EMPTY_DBL());
652 const bool emptyMax = (max == EMPTY_DBL());
653
654 if (!emptyMin && !emptyMax && max < min) {
655 std::stringstream ss;
656 ss << "TimeSeriesProperty::expandFilterToRange: 'max' argument must be "
657 "greater than 'min' "
658 << "(got min=" << min << " max=" << max << ")";
659 throw std::invalid_argument(ss.str());
660 }
661
662 // If min or max were unset ("empty") in the algorithm, set to the min or max
663 // value of the log
664 if (emptyMin)
665 min = static_cast<double>(minValue());
666 if (emptyMax)
667 max = static_cast<double>(maxValue());
668
669 if (range.start() < firstTime()) {
670 // Assume everything before the 1st value is constant
671 double val = static_cast<double>(firstValue());
672 if ((val >= min) && (val <= max)) {
673 SplittingIntervalVec extraFilter;
674 extraFilter.emplace_back(range.start(), firstTime(), 0);
675 // Include everything from the start of the run to the first time measured
676 // (which may be a null time interval; this'll be ignored)
677 split = split | extraFilter;
678 }
679 }
680
681 if (lastTime() < range.stop()) {
682 // Assume everything after the LAST value is constant
683 double val = static_cast<double>(lastValue());
684 if ((val >= min) && (val <= max)) {
685 SplittingIntervalVec extraFilter;
686 extraFilter.emplace_back(lastTime(), range.stop(), 0);
687 // Include everything from the start of the run to the first time measured
688 // (which may be a null time interval; this'll be ignored)
689 split = split | extraFilter;
690 }
691 }
692}
693
697template <>
698void TimeSeriesProperty<std::string>::expandFilterToRange(std::vector<SplittingInterval> & /*split*/, double /*min*/,
699 double /*max*/, const TimeInterval & /*range*/) const {
700 throw Exception::NotImplementedError("TimeSeriesProperty::makeFilterByValue "
701 "is not implemented for string "
702 "properties");
703}
704
709template <typename TYPE> double TimeSeriesProperty<TYPE>::timeAverageValue(const TimeROI *timeRoi) const {
710 double retVal = 0.0;
711 try {
712 if ((timeRoi == nullptr) || (timeRoi->useAll())) {
713 const auto &intervals = getTimeIntervals();
714 retVal = this->averageValueInFilter(intervals);
715 } else if (timeRoi->useNone()) {
716 // if TimeROI bans everything, use the simple mean
717 const auto stats =
719 return stats.mean;
720 } else {
721 const auto &filter = timeRoi->toTimeIntervals();
722 retVal = this->averageValueInFilter(filter);
723 g_log.warning("Calls to TimeSeriesProperty::timeAverageValue should be replaced with "
724 "Run::getTimeAveragedValue");
725 }
726 } catch (std::exception &) {
727 // just return nan
728 retVal = std::numeric_limits<double>::quiet_NaN();
729 }
730 return retVal;
731}
732
733template <> double TimeSeriesProperty<std::string>::timeAverageValue(const TimeROI * /*timeRoi*/) const {
734 throw Exception::NotImplementedError("TimeSeriesProperty::timeAverageValue is not implemented for string properties");
735}
736
743template <typename TYPE>
744double TimeSeriesProperty<TYPE>::averageValueInFilter(const std::vector<TimeInterval> &filter) const {
745 // TODO: Consider logs that aren't giving starting values.
746
747 // If there's just a single value in the log, return that.
748 if (size() == 1) {
749 return static_cast<double>(this->firstValue());
750 }
751
752 // First of all, if the log or the filter is empty, return NaN
753 if (realSize() == 0 || filter.empty()) {
754 return std::numeric_limits<double>::quiet_NaN();
755 }
756
757 sortIfNecessary();
758
759 double numerator(0.0), totalTime(0.0);
760 // Loop through the filter ranges
761 for (const auto &time : filter) {
762 // Calculate the total time duration (in seconds) within by the filter
763 totalTime += time.duration();
764
765 // Get the log value and index at the start time of the filter
766 int index;
767 double currentValue = static_cast<double>(getSingleValue(time.start(), index));
768 DateAndTime startTime = time.start();
769
770 while (index < realSize() - 1 && m_values[index + 1].time() < time.stop()) {
771 ++index;
772 numerator += DateAndTime::secondsFromDuration(m_values[index].time() - startTime) * currentValue;
773 startTime = m_values[index].time();
774 currentValue = static_cast<double>(m_values[index].value());
775 }
776
777 // Now close off with the end of the current filter range
778 numerator += DateAndTime::secondsFromDuration(time.stop() - startTime) * currentValue;
779 }
780
781 if (totalTime > 0) {
782 // 'Normalise' by the total time
783 return numerator / totalTime;
784 } else {
785 // give simple mean
786 const auto stats = Mantid::Kernel::getStatistics(this->valuesAsVector(), Mantid::Kernel::Math::StatisticType::Mean);
787 return stats.mean;
788 }
789}
790
794template <>
795double TimeSeriesProperty<std::string>::averageValueInFilter(const std::vector<TimeInterval> & /*filter*/) const {
796 throw Exception::NotImplementedError("TimeSeriesProperty::"
797 "averageValueInFilter is not "
798 "implemented for string properties");
799}
800
801template <typename TYPE>
802std::pair<double, double>
803TimeSeriesProperty<TYPE>::averageAndStdDevInFilter(const std::vector<TimeInterval> &intervals) const {
804 double mean_prev, mean_current(0.0), s(0.0), variance, duration, weighted_sum(0.0);
805
806 // First of all, if the log or the intervals are empty or is a single value,
807 // return NaN for the uncertainty
808 if (realSize() <= 1 || intervals.empty()) {
809 return std::pair<double, double>{this->averageValueInFilter(intervals), std::numeric_limits<double>::quiet_NaN()};
810 }
811 auto real_size = realSize();
812 for (const auto &time : intervals) {
813 int index;
814 auto currentValue = static_cast<double>(getSingleValue(time.start(), index));
815 DateAndTime startTime = time.start();
816 while (index < realSize() - 1 && m_values[index + 1].time() < time.stop()) {
817 index++;
818 if (index == real_size) {
819 duration = DateAndTime::secondsFromDuration(time.stop() - startTime);
820 } else {
821 duration = DateAndTime::secondsFromDuration(m_values[index].time() - startTime);
822 startTime = m_values[index].time();
823 }
824 mean_prev = mean_current;
825 if (duration > 0.) {
826 weighted_sum += duration;
827
828 mean_current = mean_prev + (duration / weighted_sum) * (currentValue - mean_prev);
829 s += duration * (currentValue - mean_prev) * (currentValue - mean_current);
830 }
831 currentValue = static_cast<double>(m_values[index].value());
832 }
833
834 // Now close off with the end of the current filter range
835 duration = DateAndTime::secondsFromDuration(time.stop() - startTime);
836 if (duration > 0.) {
837 weighted_sum += duration;
838 mean_prev = mean_current;
839
840 mean_current = mean_prev + (duration / weighted_sum) * (currentValue - mean_prev);
841 s += duration * (currentValue - mean_prev) * (currentValue - mean_current);
842 }
843 }
844 variance = s / weighted_sum;
845 // Normalise by the total time
846 return std::pair<double, double>{mean_current, std::sqrt(variance)};
847}
848
852template <>
853std::pair<double, double>
854TimeSeriesProperty<std::string>::averageAndStdDevInFilter(const std::vector<TimeInterval> & /*filter*/) const {
855 throw Exception::NotImplementedError("TimeSeriesProperty::"
856 "averageAndStdDevInFilter is not "
857 "implemented for string properties");
858}
859
860template <typename TYPE>
861std::pair<double, double> TimeSeriesProperty<TYPE>::timeAverageValueAndStdDev(const Kernel::TimeROI *timeRoi) const {
862 // time series with less than two entries are conner cases
863 if (this->realSize() == 0)
864 return std::pair<double, double>{std::numeric_limits<double>::quiet_NaN(),
865 std::numeric_limits<double>::quiet_NaN()};
866 else if (this->realSize() == 1)
867 return std::pair<double, double>(static_cast<double>(this->firstValue()), 0.0);
868
869 // Derive splitting intervals from either the roi or from the first/last entries in the time series
870 std::vector<TimeInterval> intervals;
871 if (timeRoi && !timeRoi->useAll()) {
872 intervals = timeRoi->toTimeIntervals(this->firstTime());
873 } else {
874 intervals = this->getTimeIntervals();
875 }
876
877 return this->averageAndStdDevInFilter(intervals);
878}
879
883template <>
884std::pair<double, double>
887 "TimeSeriesProperty::timeAverageValueAndStdDev is not implemented for string properties");
888}
889
890// Re-enable the warnings disabled before makeFilterByValue
891#ifdef _WIN32
892#pragma warning(pop)
893#endif
894#if defined(__GNUC__) && !(defined(__INTEL_COMPILER))
895#pragma GCC diagnostic warning "-Wconversion"
896#endif
897
904template <typename TYPE> std::map<DateAndTime, TYPE> TimeSeriesProperty<TYPE>::valueAsCorrectMap() const {
905 // 1. Sort if necessary
906 sortIfNecessary();
907
908 // 2. Data Strcture
909 std::map<DateAndTime, TYPE> asMap;
910
911 if (!m_values.empty()) {
912 for (size_t i = 0; i < m_values.size(); i++)
913 asMap[m_values[i].time()] = m_values[i].value();
914 }
915
916 return asMap;
917}
918
923template <typename TYPE> std::vector<TYPE> TimeSeriesProperty<TYPE>::valuesAsVector() const {
924 sortIfNecessary();
925
926 std::vector<TYPE> out;
927 out.reserve(m_values.size());
928
929 for (size_t i = 0; i < m_values.size(); i++)
930 out.emplace_back(m_values[i].value());
931
932 return out;
933}
934
941template <typename TYPE> std::multimap<DateAndTime, TYPE> TimeSeriesProperty<TYPE>::valueAsMultiMap() const {
942 std::multimap<DateAndTime, TYPE> asMultiMap;
943
944 if (!m_values.empty()) {
945 for (size_t i = 0; i < m_values.size(); i++)
946 asMultiMap.insert(std::make_pair(m_values[i].time(), m_values[i].value()));
947 }
948
949 return asMultiMap;
950}
951
956template <typename TYPE> std::vector<DateAndTime> TimeSeriesProperty<TYPE>::timesAsVector() const {
957 sortIfNecessary();
958
959 std::vector<DateAndTime> out;
960 out.reserve(m_values.size());
961
962 for (size_t i = 0; i < m_values.size(); i++) {
963 out.emplace_back(m_values[i].time());
964 }
965
966 return out;
967}
968
978template <typename TYPE>
979std::vector<DateAndTime> TimeSeriesProperty<TYPE>::filteredTimesAsVector(const Kernel::TimeROI *roi) const {
980 if (roi && !roi->useAll()) {
981 this->sortIfNecessary();
982 std::vector<DateAndTime> filteredTimes;
983 if (roi->firstTime() > this->m_values.back().time()) {
984 // Since the ROI starts after everything, just return the last time in the log
985 filteredTimes.emplace_back(roi->firstTime());
986 } else { // only use the times in the filter - this is very similar to FilteredTimeSeriesProperty::applyFilter
987 // the index into the m_values array of the time, or -1 (before) or m_values.size() (after)
988 std::size_t index_current_log{0};
989
990 for (const auto &splitter : roi->toTimeIntervals()) {
991 const auto endTime = splitter.stop();
992
993 // check if the splitter starts too early
994 if (endTime < this->m_values[index_current_log].time()) {
995 continue; // skip to the next splitter
996 }
997
998 // cache values to reduce number of method calls
999 const auto beginTime = splitter.start();
1000
1001 // find the first log that should be added
1002 if (this->m_values.back().time() < beginTime) {
1003 // skip directly to the end if the filter starts after the last log
1004 index_current_log = this->m_values.size() - 1;
1005 } else {
1006 // search for the right starting point
1007 while ((this->m_values[index_current_log].time() <= beginTime)) {
1008 if (index_current_log + 1 > this->m_values.size())
1009 break;
1010 index_current_log++;
1011 }
1012 // need to back up by one
1013 if (index_current_log > 0)
1014 index_current_log--;
1015 // go backwards more while times are equal to the one being started at
1016 while (index_current_log > 0 &&
1017 this->m_values[index_current_log].time() == this->m_values[index_current_log - 1].time()) {
1018 index_current_log--;
1019 }
1020 }
1021
1022 // add everything up to the end time
1023 for (; index_current_log < this->m_values.size(); ++index_current_log) {
1024 if (this->m_values[index_current_log].time() >= endTime)
1025 break;
1026
1027 // start time is when this value was created or when the filter started
1028 filteredTimes.emplace_back(std::max(beginTime, this->m_values[index_current_log].time()));
1029 }
1030 // go back one so the next splitter can add a value
1031 if (index_current_log > 0)
1032 index_current_log--;
1033 }
1034 }
1035
1036 return filteredTimes;
1037 } else {
1038 return this->timesAsVector();
1039 }
1040}
1041
1042template <typename TYPE> std::vector<DateAndTime> TimeSeriesProperty<TYPE>::filteredTimesAsVector() const {
1043 return this->timesAsVector();
1044}
1045
1050template <typename TYPE> std::vector<double> TimeSeriesProperty<TYPE>::timesAsVectorSeconds() const {
1051 // 1. Sort if necessary
1052 sortIfNecessary();
1053
1054 // 2. Output data structure
1055 std::vector<double> out;
1056 if (!m_values.empty()) {
1057 out.reserve(m_values.size());
1058
1059 Types::Core::DateAndTime start = m_values[0].time();
1060 for (size_t i = 0; i < m_values.size(); i++) {
1061 out.emplace_back(DateAndTime::secondsFromDuration(m_values[i].time() - start));
1062 }
1063 }
1064
1065 return out;
1066}
1067
1073template <typename TYPE>
1074std::vector<double> TimeSeriesProperty<TYPE>::timesAsVectorSeconds(Types::Core::DateAndTime start) const {
1075 // 1. Sort if necessary
1076 sortIfNecessary();
1077
1078 // 2. Output data structure
1079 std::vector<double> out;
1080 if (!m_values.empty()) {
1081 out.reserve(m_values.size());
1082 for (size_t i = 0; i < m_values.size(); i++) {
1083 out.emplace_back(DateAndTime::secondsFromDuration(m_values[i].time() - start));
1084 }
1085 }
1086
1087 return out;
1088}
1089
1095template <typename TYPE>
1096void TimeSeriesProperty<TYPE>::addValue(const Types::Core::DateAndTime &time, const TYPE &value) {
1097 TimeValueUnit<TYPE> newvalue(time, value);
1098 // Add the value to the back of the vector
1099 m_values.emplace_back(newvalue);
1100 // Increment the separate record of the property's size
1101 m_size++;
1102
1103 // Toggle the sorted flag if necessary
1104 // (i.e. if the flag says we're sorted and the added time is before the prior
1105 // last time)
1106 if (m_size == 1) {
1107 // First item, must be sorted.
1108 m_propSortedFlag = TimeSeriesSortStatus::TSSORTED;
1109 } else if (m_propSortedFlag == TimeSeriesSortStatus::TSUNKNOWN && m_values.back() < *(m_values.rbegin() + 1)) {
1110 // Previously unknown and still unknown
1111 m_propSortedFlag = TimeSeriesSortStatus::TSUNSORTED;
1112 } else if (m_propSortedFlag == TimeSeriesSortStatus::TSSORTED && m_values.back() < *(m_values.rbegin() + 1)) {
1113 // Previously sorted but last added is not in order
1114 m_propSortedFlag = TimeSeriesSortStatus::TSUNSORTED;
1115 }
1116
1117 // m_filterApplied = false;
1118}
1119
1125template <typename TYPE> void TimeSeriesProperty<TYPE>::addValue(const std::string &time, const TYPE &value) {
1126 return addValue(Types::Core::DateAndTime(time), value);
1127}
1128
1134template <typename TYPE> void TimeSeriesProperty<TYPE>::addValue(const std::time_t &time, const TYPE &value) {
1135 Types::Core::DateAndTime dt;
1136 dt.set_from_time_t(time);
1137 return addValue(dt, value);
1138}
1139
1145template <typename TYPE>
1146void TimeSeriesProperty<TYPE>::addValues(const std::vector<Types::Core::DateAndTime> &times,
1147 const std::vector<TYPE> &values) {
1148 size_t length = std::min(times.size(), values.size());
1149 m_size += static_cast<int>(length);
1150 for (size_t i = 0; i < length; ++i) {
1151 m_values.emplace_back(times[i], values[i]);
1152 }
1153
1154 if (!values.empty())
1155 m_propSortedFlag = TimeSeriesSortStatus::TSUNKNOWN;
1156}
1157
1163template <typename TYPE>
1164void TimeSeriesProperty<TYPE>::replaceValues(const std::vector<Types::Core::DateAndTime> &times,
1165 const std::vector<TYPE> &values) {
1166 clear();
1167 addValues(times, values);
1168}
1169
1174template <typename TYPE> DateAndTime TimeSeriesProperty<TYPE>::lastTime() const {
1175 if (m_values.empty()) {
1176 const std::string error("lastTime(): TimeSeriesProperty '" + name() + "' is empty");
1177 g_log.debug(error);
1178 throw std::runtime_error(error);
1179 }
1180
1181 sortIfNecessary();
1182
1183 return m_values.rbegin()->time();
1184}
1185
1189template <typename TYPE> TYPE TimeSeriesProperty<TYPE>::firstValue() const {
1190 if (m_values.empty()) {
1191 const std::string error("firstValue(): TimeSeriesProperty '" + name() + "' is empty");
1192 g_log.debug(error);
1193 throw std::runtime_error(error);
1194 }
1195
1196 sortIfNecessary();
1197
1198 return m_values[0].value();
1199}
1200
1201template <typename TYPE> TYPE TimeSeriesProperty<TYPE>::firstValue(const Kernel::TimeROI &roi) const {
1202 const auto startTime = roi.firstTime();
1203 if (startTime <= this->firstTime()) {
1204 return this->firstValue();
1205 } else if (startTime >= this->lastTime()) {
1206 return this->lastValue();
1207 } else {
1208 const auto times = this->timesAsVector();
1209 auto iter = std::lower_bound(times.cbegin(), times.cend(), startTime);
1210 if (*iter > startTime)
1211 iter--;
1212 const auto index = std::size_t(std::distance(times.cbegin(), iter));
1213
1214 const auto values = this->valuesAsVector();
1215 const TYPE ret = values[index];
1216 return ret;
1217 }
1218}
1219
1223template <typename TYPE> DateAndTime TimeSeriesProperty<TYPE>::firstTime() const {
1224 if (m_values.empty()) {
1225 const std::string error("firstTime(): TimeSeriesProperty '" + name() + "' is empty");
1226 g_log.debug(error);
1227 throw std::runtime_error(error);
1228 }
1229
1230 sortIfNecessary();
1231
1232 return m_values[0].time();
1233}
1234
1239template <typename TYPE> TYPE TimeSeriesProperty<TYPE>::lastValue() const {
1240 if (m_values.empty()) {
1241 const std::string error("lastValue(): TimeSeriesProperty '" + name() + "' is empty");
1242 g_log.debug(error);
1243 throw std::runtime_error(error);
1244 }
1245
1246 sortIfNecessary();
1247
1248 return m_values.rbegin()->value();
1249}
1250
1251template <typename TYPE> TYPE TimeSeriesProperty<TYPE>::lastValue(const Kernel::TimeROI &roi) const {
1252 const auto stopTime = roi.lastTime();
1253 const auto times = this->timesAsVector();
1254 if (stopTime <= times.front()) {
1255 return this->firstValue();
1256 } else if (stopTime >= times.back()) {
1257 return this->lastValue();
1258 } else {
1259 auto iter = std::lower_bound(times.cbegin(), times.cend(), stopTime);
1260 if ((iter != times.cbegin()) && (*iter > stopTime))
1261 --iter;
1262 const auto index = std::size_t(std::distance(times.cbegin(), iter));
1263
1264 const auto values = this->valuesAsVector();
1265 const TYPE ret = values[index];
1266 return ret;
1267 }
1268}
1269
1279template <typename TYPE> double TimeSeriesProperty<TYPE>::durationInSeconds(const Kernel::TimeROI *roi) const {
1280 if (this->size() == 0)
1281 return std::numeric_limits<double>::quiet_NaN();
1282 if (roi && !roi->useAll()) {
1283 Kernel::TimeROI seriesSpan(*roi);
1284 const auto thisFirstTime = this->firstTime();
1285 // remove everything before the start time
1286 if (thisFirstTime > DateAndTime::GPS_EPOCH) {
1287 seriesSpan.addMask(DateAndTime::GPS_EPOCH, thisFirstTime);
1288 }
1289 return seriesSpan.durationInSeconds();
1290 } else {
1291 const auto &intervals = this->getTimeIntervals();
1292 const double duration_sec =
1293 std::accumulate(intervals.cbegin(), intervals.cend(), 0.,
1294 [](double sum, const auto &interval) { return sum + interval.duration(); });
1295 return duration_sec;
1296 }
1297}
1298
1299template <typename TYPE> TYPE TimeSeriesProperty<TYPE>::minValue() const {
1300 return std::min_element(m_values.begin(), m_values.end(), TimeValueUnit<TYPE>::valueCmp)->value();
1301}
1302
1303template <typename TYPE> TYPE TimeSeriesProperty<TYPE>::maxValue() const {
1304 return std::max_element(m_values.begin(), m_values.end(), TimeValueUnit<TYPE>::valueCmp)->value();
1305}
1306
1307template <typename TYPE> double TimeSeriesProperty<TYPE>::mean() const {
1308 Mantid::Kernel::Statistics raw_stats =
1309 Mantid::Kernel::getStatistics(this->filteredValuesAsVector(), StatOptions::Mean);
1310 return raw_stats.mean;
1311}
1312
1315template <typename TYPE> int TimeSeriesProperty<TYPE>::size() const { return m_size; }
1316
1321template <typename TYPE> int TimeSeriesProperty<TYPE>::realSize() const { return static_cast<int>(m_values.size()); }
1322
1323/*
1324 * Get the time series property as a string of 'time value'
1325 * @return time series property as a string
1326 */
1327template <typename TYPE> std::string TimeSeriesProperty<TYPE>::value() const {
1328 sortIfNecessary();
1329
1330 std::stringstream ins;
1331 for (size_t i = 0; i < m_values.size(); i++) {
1332 try {
1333 ins << m_values[i].time().toSimpleString();
1334 ins << " " << m_values[i].value() << "\n";
1335 } catch (...) {
1336 // Some kind of error; for example, invalid year, can occur when
1337 // converting boost time.
1338 ins << "Error Error"
1339 << "\n";
1340 }
1341 }
1342
1343 return ins.str();
1344}
1345
1350template <typename TYPE> std::vector<std::string> TimeSeriesProperty<TYPE>::time_tValue() const {
1351 sortIfNecessary();
1352
1353 std::vector<std::string> values;
1354 values.reserve(m_values.size());
1355
1356 for (size_t i = 0; i < m_values.size(); i++) {
1357 std::stringstream line;
1358 line << m_values[i].time().toSimpleString() << " " << m_values[i].value();
1359 values.emplace_back(line.str());
1360 }
1361
1362 return values;
1363}
1364
1372template <typename TYPE> std::map<DateAndTime, TYPE> TimeSeriesProperty<TYPE>::valueAsMap() const {
1373 // 1. Sort if necessary
1374 sortIfNecessary();
1375
1376 // 2. Build map
1377
1378 std::map<DateAndTime, TYPE> asMap;
1379 if (m_values.empty())
1380 return asMap;
1381
1382 TYPE d = m_values[0].value();
1383 asMap[m_values[0].time()] = d;
1384
1385 for (size_t i = 1; i < m_values.size(); i++) {
1386 if (m_values[i].value() != d) {
1387 // Only put entry with different value from last entry to map
1388 asMap[m_values[i].time()] = m_values[i].value();
1389 d = m_values[i].value();
1390 }
1391 }
1392 return asMap;
1393}
1394
1400template <typename TYPE> std::string TimeSeriesProperty<TYPE>::setValue(const std::string & /*unused*/) {
1401 throw Exception::NotImplementedError("TimeSeriesProperty<TYPE>::setValue - "
1402 "Cannot extract TimeSeries from a "
1403 "std::string");
1404}
1405
1411template <typename TYPE> std::string TimeSeriesProperty<TYPE>::setValueFromJson(const Json::Value & /*unused*/) {
1412 throw Exception::NotImplementedError("TimeSeriesProperty<TYPE>::setValue - "
1413 "Cannot extract TimeSeries from a "
1414 "Json::Value");
1415}
1416
1421template <typename TYPE>
1422std::string TimeSeriesProperty<TYPE>::setDataItem(const std::shared_ptr<DataItem> & /*unused*/) {
1423 throw Exception::NotImplementedError("TimeSeriesProperty<TYPE>::setValue - "
1424 "Cannot extract TimeSeries from "
1425 "DataItem");
1426}
1427
1430template <typename TYPE> void TimeSeriesProperty<TYPE>::clear() {
1431 m_size = 0;
1432 m_values.clear();
1433
1434 m_propSortedFlag = TimeSeriesSortStatus::TSSORTED;
1435 // m_filterApplied = false;
1436}
1437
1445template <typename TYPE> void TimeSeriesProperty<TYPE>::clearOutdated() {
1446 if (realSize() > 1) {
1447 auto lastValueInVec = m_values.back();
1448 clear();
1449 m_values.emplace_back(lastValueInVec);
1450 m_size = 1;
1451 }
1452}
1453
1454//--------------------------------------------------------------------------------------------
1464template <typename TYPE>
1465void TimeSeriesProperty<TYPE>::create(const Types::Core::DateAndTime &start_time, const std::vector<double> &time_sec,
1466 const std::vector<TYPE> &new_values) {
1467 if (time_sec.size() != new_values.size()) {
1468 std::stringstream msg;
1469 msg << "TimeSeriesProperty \"" << name() << "\" create: mismatched size "
1470 << "for the time and values vectors.";
1471 throw std::invalid_argument(msg.str());
1472 }
1473
1474 clear();
1475 const std::size_t num = new_values.size();
1476 m_values.reserve(num);
1477
1478 // set the sorted flag
1479 if (std::is_sorted(time_sec.cbegin(), time_sec.cend()))
1480 m_propSortedFlag = TimeSeriesSortStatus::TSSORTED;
1481 else
1482 m_propSortedFlag = TimeSeriesSortStatus::TSUNSORTED;
1483 // set the values
1484 constexpr double SEC_TO_NANO{1000000000.0};
1485 const int64_t start_time_ns = start_time.totalNanoseconds();
1486 for (std::size_t i = 0; i < num; i++) {
1487 m_values.emplace_back(start_time_ns + static_cast<int64_t>(time_sec[i] * SEC_TO_NANO), new_values[i]);
1488 }
1489
1490 // reset the size
1491 m_size = static_cast<int>(m_values.size());
1492}
1493
1494//--------------------------------------------------------------------------------------------
1503template <typename TYPE>
1504void TimeSeriesProperty<TYPE>::create(const std::vector<Types::Core::DateAndTime> &new_times,
1505 const std::vector<TYPE> &new_values) {
1506 if (new_times.size() != new_values.size())
1507 throw std::invalid_argument("TimeSeriesProperty::create: mismatched size "
1508 "for the time and values vectors.");
1509
1510 clear();
1511 // nothing to do without values
1512 if (new_times.empty()) {
1513 m_propSortedFlag = TimeSeriesSortStatus::TSSORTED;
1514 return;
1515 }
1516
1517 const std::size_t num = new_values.size();
1518 m_values.reserve(num);
1519
1520 // set the sorted flag
1521 if (std::is_sorted(new_times.cbegin(), new_times.cend()))
1522 m_propSortedFlag = TimeSeriesSortStatus::TSSORTED;
1523 else
1524 m_propSortedFlag = TimeSeriesSortStatus::TSUNSORTED;
1525 // add the values
1526 for (std::size_t i = 0; i < num; i++) {
1527 m_values.emplace_back(new_times[i], new_values[i]);
1528 }
1529
1530 // reset the size
1531 m_size = static_cast<int>(m_values.size());
1532}
1533
1538template <typename TYPE> TYPE TimeSeriesProperty<TYPE>::getSingleValue(const Types::Core::DateAndTime &t) const {
1539 if (m_values.empty()) {
1540 const std::string error("getSingleValue(): TimeSeriesProperty '" + name() + "' is empty");
1541 g_log.debug(error);
1542 throw std::runtime_error(error);
1543 }
1544
1545 // 1. Get sorted
1546 sortIfNecessary();
1547
1548 // 2.
1549 TYPE valueAtTime;
1550 if (t < m_values[0].time()) {
1551 // 1. Out side of lower bound
1552 valueAtTime = m_values[0].value();
1553 } else if (t >= m_values.back().time()) {
1554 // 2. Out side of upper bound
1555 valueAtTime = m_values.back().value();
1556 } else {
1557 // 3. Within boundary
1558 int index = this->findIndex(t);
1559
1560 if (index < 0) {
1561 // If query time "t" is earlier than the begin time of the series
1562 index = 0;
1563 } else if (index == int(m_values.size())) {
1564 // If query time "t" is later than the end time of the series
1565 index = static_cast<int>(m_values.size()) - 1;
1566 } else if (index > int(m_values.size())) {
1567 std::stringstream errss;
1568 errss << "TimeSeriesProperty.findIndex() returns index (" << index << " ) > maximum defined value "
1569 << m_values.size();
1570 throw std::logic_error(errss.str());
1571 }
1572
1573 valueAtTime = m_values[static_cast<size_t>(index)].value();
1574 }
1575
1576 return valueAtTime;
1577} // END-DEF getSinglevalue()
1578
1584template <typename TYPE>
1585TYPE TimeSeriesProperty<TYPE>::getSingleValue(const Types::Core::DateAndTime &t, int &index) const {
1586 if (m_values.empty()) {
1587 const std::string error("getSingleValue(): TimeSeriesProperty '" + name() + "' is empty");
1588 g_log.debug(error);
1589 throw std::runtime_error(error);
1590 }
1591
1592 // 1. Get sorted
1593 sortIfNecessary();
1594
1595 // 2.
1596 TYPE valueAtTime;
1597 if (t < m_values[0].time()) {
1598 // 1. Out side of lower bound
1599 valueAtTime = m_values[0].value();
1600 index = 0;
1601 } else if (t >= m_values.back().time()) {
1602 // 2. Out side of upper bound
1603 valueAtTime = m_values.back().value();
1604 index = int(m_values.size()) - 1;
1605 } else {
1606 // 3. Within boundary
1607 index = this->findIndex(t);
1608
1609 if (index < 0) {
1610 // If query time "t" is earlier than the begin time of the series
1611 index = 0;
1612 } else if (index == int(m_values.size())) {
1613 // If query time "t" is later than the end time of the series
1614 index = static_cast<int>(m_values.size()) - 1;
1615 } else if (index > int(m_values.size())) {
1616 std::stringstream errss;
1617 errss << "TimeSeriesProperty.findIndex() returns index (" << index << " ) > maximum defined value "
1618 << m_values.size();
1619 throw std::logic_error(errss.str());
1620 }
1621
1622 valueAtTime = m_values[static_cast<size_t>(index)].value();
1623 }
1624
1625 return valueAtTime;
1626} // END-DEF getSinglevalue()
1627
1638template <typename TYPE> TimeInterval TimeSeriesProperty<TYPE>::nthInterval(int n) const {
1639 // 0. Throw exception
1640 if (m_values.empty()) {
1641 const std::string error("nthInterval(): TimeSeriesProperty '" + name() + "' is empty");
1642 g_log.debug(error);
1643 throw std::runtime_error(error);
1644 }
1645
1646 // 1. Sort
1647 sortIfNecessary();
1648
1649 // 2. Calculate time interval
1650
1651 Kernel::TimeInterval deltaT;
1652
1653 // No filter
1654 if (n >= static_cast<int>(m_values.size()) || (n == static_cast<int>(m_values.size()) - 1 && m_values.size() == 1)) {
1655 // Out of bound
1656 ;
1657 } else if (n == static_cast<int>(m_values.size()) - 1) {
1658 // Last one by making up an end time.
1659 DateAndTime endTime = getFakeEndTime();
1660
1661 deltaT = Kernel::TimeInterval(m_values.rbegin()->time(), endTime);
1662 } else {
1663 // Regular
1664 DateAndTime startT = m_values[static_cast<std::size_t>(n)].time();
1665 DateAndTime endT = m_values[static_cast<std::size_t>(n) + 1].time();
1666 TimeInterval dt(startT, endT);
1667 deltaT = dt;
1668 }
1669
1670 return deltaT;
1671}
1672
1673template <typename TYPE> Types::Core::DateAndTime TimeSeriesProperty<TYPE>::getFakeEndTime() const {
1674 sortIfNecessary();
1675
1676 // the last time is the last thing known
1677 const auto ultimate = m_values.rbegin()->time();
1678
1679 // go backwards from the time before it that is different
1680 int counter = 0;
1681 while (DateAndTime::secondsFromDuration(ultimate - (m_values.rbegin() + counter)->time()) == 0.) {
1682 counter += 1;
1683 }
1684
1685 // get the last time that is different
1686 time_duration lastDuration = m_values.rbegin()->time() - (m_values.rbegin() + counter)->time();
1687
1688 // the last duration is equal to the previous, non-zero, duration
1689 return m_values.rbegin()->time() + lastDuration;
1690}
1691
1692//-----------------------------------------------------------------------------------------------
1698template <typename TYPE> TYPE TimeSeriesProperty<TYPE>::nthValue(int n) const {
1699
1700 // 1. Throw error if property is empty
1701 if (m_values.empty()) {
1702 const std::string error("nthValue(): TimeSeriesProperty '" + name() + "' is empty");
1703 g_log.debug(error);
1704 throw std::runtime_error(error);
1705 }
1706
1707 // 2. Sort and apply filter
1708 sortIfNecessary();
1709
1710 TYPE nthValue;
1711
1712 // 3. Situation 1: No filter
1713 if (static_cast<size_t>(n) < m_values.size()) {
1714 const auto entry = m_values[static_cast<std::size_t>(n)];
1715 nthValue = entry.value();
1716 } else {
1717 const auto entry = m_values[static_cast<std::size_t>(m_size) - 1];
1718 nthValue = entry.value();
1719 }
1720
1721 return nthValue;
1722}
1723
1729template <typename TYPE> Types::Core::DateAndTime TimeSeriesProperty<TYPE>::nthTime(int n) const {
1730 sortIfNecessary();
1731
1732 if (m_values.empty()) {
1733 const std::string error("nthTime(): TimeSeriesProperty '" + name() + "' is empty");
1734 g_log.debug(error);
1735 throw std::runtime_error(error);
1736 }
1737
1738 if (n < 0 || n >= static_cast<int>(m_values.size()))
1739 n = static_cast<int>(m_values.size()) - 1;
1740
1741 return m_values[static_cast<size_t>(n)].time();
1742}
1743
1747template <typename TYPE> void TimeSeriesProperty<TYPE>::countSize() const {
1748 // 1. Not filter
1749 m_size = int(m_values.size());
1750}
1751
1756template <typename TYPE> bool TimeSeriesProperty<TYPE>::isTimeString(const std::string &str) {
1757 static const boost::regex re("^[0-9]{4}.[0-9]{2}.[0-9]{2}.[0-9]{2}.[0-9]{2}.[0-9]{2}");
1758 return boost::regex_search(str.begin(), str.end(), re);
1759}
1760
1765template <typename TYPE> std::string TimeSeriesProperty<TYPE>::isValid() const { return ""; }
1766
1767/*
1768 * A TimeSeriesProperty never has a default, so return empty string
1769 * @returns Empty string as no defaults can be provided
1770 */
1771template <typename TYPE> std::string TimeSeriesProperty<TYPE>::getDefault() const {
1772 return ""; // No defaults can be provided=empty string
1773}
1774
1778template <typename TYPE> bool TimeSeriesProperty<TYPE>::isDefault() const { return false; }
1779
1787template <typename TYPE>
1789 // Start with statistics that are not time-weighted
1790 TimeSeriesPropertyStatistics out(Mantid::Kernel::getStatistics(this->filteredValuesAsVector(roi)));
1791 out.duration = this->durationInSeconds(roi);
1792 if (out.standard_deviation == 0.) {
1793 // if the simple std-dev is zero, just copy the simple mean
1794 out.time_mean = out.mean;
1795 out.time_standard_deviation = 0.;
1796 } else {
1797 // follow with time-weighted statistics
1798 auto avAndDev = this->timeAverageValueAndStdDev(roi);
1799 out.time_mean = avAndDev.first;
1800 out.time_standard_deviation = avAndDev.second;
1801 }
1802 return out;
1803}
1804
1805template <>
1807 // statistics of a string property doesn't make sense
1809 out.setAllToNan();
1810
1811 return out;
1812}
1813
1819template <typename TYPE>
1821 using namespace Kernel::Math;
1822 double singleValue = 0;
1823 switch (selection) {
1824 case FirstValue:
1825 if (roi && !roi->useAll())
1826 singleValue = double(this->firstValue(*roi));
1827 else
1828 singleValue = double(this->nthValue(0));
1829 break;
1830 case LastValue:
1831 if (roi && !roi->useAll())
1832 singleValue = double(this->lastValue(*roi));
1833 else
1834 singleValue = double(this->nthValue(this->size() - 1));
1835 break;
1836 case Minimum:
1837 singleValue = static_cast<double>(this->getStatistics(roi).minimum);
1838 break;
1839 case Maximum:
1840 singleValue = static_cast<double>(this->getStatistics(roi).maximum);
1841 break;
1842 case Mean:
1843 singleValue = this->getStatistics(roi).mean;
1844 break;
1845 case Median:
1846 singleValue = this->getStatistics(roi).median;
1847 break;
1848 case TimeAveragedMean:
1849 singleValue = this->getStatistics(roi).time_mean;
1850 break;
1851 case StdDev:
1852 singleValue = this->getStatistics(roi).standard_deviation;
1853 break;
1854 case TimeAverageStdDev:
1855 singleValue = this->getStatistics(roi).time_standard_deviation;
1856 break;
1857 default:
1858 throw std::invalid_argument("extractStatistic - Unknown statistic type: " + boost::lexical_cast<std::string>(this));
1859 };
1860 return singleValue;
1861}
1862
1866template <>
1868 UNUSED_ARG(selection);
1869 UNUSED_ARG(roi);
1870 throw Exception::NotImplementedError("TimeSeriesProperty::"
1871 "extractStatistic is not "
1872 "implemented for string properties");
1873}
1874
1875/*
1876 * Detects whether there are duplicated entries (of time) in property
1877 * If there is any, keep one of them
1878 */
1879template <typename TYPE> void TimeSeriesProperty<TYPE>::eliminateDuplicates() {
1880 // ensure that the values are sorted
1881 sortIfNecessary();
1882
1883 // cache the original size so the number removed can be reported
1884 const auto origSize{m_size};
1885
1886 // remove the first n-repeats
1887 // taken from
1888 // https://stackoverflow.com/questions/21060636/using-stdunique-and-vector-erase-to-remove-all-but-last-occurrence-of-duplicat
1889 auto it = std::unique(m_values.rbegin(), m_values.rend(),
1890 [](const auto &a, const auto &b) { return a.time() == b.time(); });
1891 m_values.erase(m_values.begin(), it.base());
1892
1893 // update m_size
1894 countSize();
1895
1896 // log how many values were removed
1897 const auto numremoved = origSize - m_size;
1898 if (numremoved > 0)
1899 g_log.notice() << "Log \"" << this->name() << "\" has " << numremoved << " entries removed due to duplicated time. "
1900 << "\n";
1901}
1902
1903/*
1904 * Print the content to string
1905 */
1906template <typename TYPE> std::string TimeSeriesProperty<TYPE>::toString() const {
1907 std::stringstream ss;
1908 for (size_t i = 0; i < m_values.size(); ++i)
1909 ss << m_values[i].time() << "\t\t" << m_values[i].value() << "\n";
1910
1911 return ss.str();
1912}
1913
1914//-------------------------------------------------------------------------
1915// Private methods
1916//-------------------------------------------------------------------------
1917
1918//----------------------------------------------------------------------------------
1919/*
1920 * Sort vector mP and set the flag. Only sorts if the values are not already
1921 * sorted.
1922 */
1923template <typename TYPE> void TimeSeriesProperty<TYPE>::sortIfNecessary() const {
1924 if (m_propSortedFlag == TimeSeriesSortStatus::TSUNKNOWN) {
1925 bool sorted = is_sorted(m_values.begin(), m_values.end());
1926 if (sorted)
1927 m_propSortedFlag = TimeSeriesSortStatus::TSSORTED;
1928 else
1929 m_propSortedFlag = TimeSeriesSortStatus::TSUNSORTED;
1930 }
1931
1932 if (m_propSortedFlag == TimeSeriesSortStatus::TSUNSORTED) {
1933 g_log.information() << "TimeSeriesProperty \"" << this->name()
1934 << "\" is not sorted. Sorting is operated on it. \n";
1935 std::stable_sort(m_values.begin(), m_values.end());
1936 m_propSortedFlag = TimeSeriesSortStatus::TSSORTED;
1937 }
1938}
1939
1946template <typename TYPE> int TimeSeriesProperty<TYPE>::findIndex(Types::Core::DateAndTime t) const {
1947 // 0. Return with an empty container
1948 if (m_values.empty())
1949 return 0;
1950
1951 // 1. Sort
1952 sortIfNecessary();
1953
1954 // 2. Extreme value
1955 if (t <= m_values[0].time()) {
1956 return -1;
1957 } else if (t >= m_values.back().time()) {
1958 return (int(m_values.size()));
1959 }
1960
1961 // 3. Find by lower_bound()
1962 typename std::vector<TimeValueUnit<TYPE>>::const_iterator fid;
1963 TimeValueUnit<TYPE> temp(t, m_values[0].value());
1964 fid = std::lower_bound(m_values.begin(), m_values.end(), temp);
1965
1966 int newindex = int(fid - m_values.begin());
1967 if (fid->time() > t)
1968 newindex--;
1969
1970 return newindex;
1971}
1972
1979template <typename TYPE>
1980int TimeSeriesProperty<TYPE>::upperBound(Types::Core::DateAndTime t, int istart, int iend) const {
1981 // 0. Check validity
1982 if (istart < 0) {
1983 throw std::invalid_argument("Start Index cannot be less than 0");
1984 }
1985 if (iend >= static_cast<int>(m_values.size())) {
1986 throw std::invalid_argument("End Index cannot exceed the boundary");
1987 }
1988 if (istart > iend) {
1989 throw std::invalid_argument("Start index cannot be greater than end index");
1990 }
1991
1992 // 1. Return instantly if it is out of boundary
1993 if (t < (m_values.begin() + istart)->time()) {
1994 return -1;
1995 }
1996 if (t > (m_values.begin() + iend)->time()) {
1997 return static_cast<int>(m_values.size());
1998 }
1999
2000 // 2. Sort
2001 sortIfNecessary();
2002
2003 // 3. Construct the pair for comparison and do lower_bound()
2004 const TimeValueUnit<TYPE> temppair(t, m_values[0].value());
2005 const auto first = m_values.cbegin() + istart;
2006 const auto last = m_values.cbegin() + iend + 1;
2007 const auto iter = std::lower_bound(first, last, temppair);
2008
2009 // 4. Calculate return value
2010 if (iter == last)
2011 throw std::runtime_error("Cannot find data");
2012 return static_cast<int>(std::distance(m_values.cbegin(), iter));
2013}
2014
2022template <typename TYPE> std::string TimeSeriesProperty<TYPE>::setValueFromProperty(const Property &right) {
2023 auto prop = dynamic_cast<const TimeSeriesProperty<TYPE> *>(&right);
2024 if (!prop) {
2025 return "Could not set value: properties have different type.";
2026 }
2027 m_values = prop->m_values;
2028 m_size = prop->m_size;
2029 m_propSortedFlag = prop->m_propSortedFlag;
2030 // m_filter = prop->m_filter;
2031 // m_filterQuickRef = prop->m_filterQuickRef;
2032 // m_filterApplied = prop->m_filterApplied;
2033 return "";
2034}
2035
2036//----------------------------------------------------------------------------------------------
2038template <typename TYPE> void TimeSeriesProperty<TYPE>::saveTimeVector(Nexus::File *file) {
2039 std::vector<DateAndTime> times = this->timesAsVector();
2040 const DateAndTime &start = times.front();
2041 std::vector<double> timeSec(times.size());
2042 for (size_t i = 0; i < times.size(); i++)
2043 timeSec[i] = static_cast<double>(times[i].totalNanoseconds() - start.totalNanoseconds()) * 1e-9;
2044 file->writeData("time", timeSec);
2045 file->openData("time");
2046 file->putAttr("start", start.toISO8601String());
2047 file->closeData();
2048}
2049
2050//----------------------------------------------------------------------------------------------
2052template <> void TimeSeriesProperty<std::string>::saveProperty(Nexus::File *file) {
2053 std::vector<std::string> values = this->valuesAsVector();
2054 if (values.empty())
2055 return;
2056 file->makeGroup(this->name(), "NXlog", true);
2057
2058 // Find the max length of any string
2059 auto max_it = std::max_element(values.begin(), values.end(),
2060 [](const std::string &a, const std::string &b) { return a.size() < b.size(); });
2061 // Increment by 1 to have the 0 terminator
2062 size_t maxlen = max_it->size() + 1;
2063 // Copy into one array
2064 std::vector<char> strs(values.size() * maxlen);
2065 size_t index = 0;
2066 for (const auto &prop : values) {
2067 std::copy(prop.begin(), prop.end(), &strs[index]);
2068 index += maxlen;
2069 }
2070
2071 const Mantid::Nexus::DimVector dims{values.size(), maxlen};
2072 file->makeData("value", NXnumtype::CHAR, dims, true);
2073 file->putData(strs.data());
2074 file->closeData();
2075 saveTimeVector(file);
2076 file->closeGroup();
2077}
2078
2085template <> void TimeSeriesProperty<bool>::saveProperty(Nexus::File *file) {
2086 std::vector<bool> value = this->valuesAsVector();
2087 if (value.empty())
2088 return;
2089 std::vector<uint8_t> asUint(value.begin(), value.end());
2090 file->makeGroup(this->name(), "NXlog", true);
2091 file->writeData("value", asUint);
2092 file->putAttr("boolean", "1");
2093 saveTimeVector(file);
2094 file->closeGroup();
2095}
2096
2097template <typename TYPE> void TimeSeriesProperty<TYPE>::saveProperty(Nexus::File *file) {
2098 auto values = this->valuesAsVector();
2099 if (values.empty())
2100 return;
2101 file->makeGroup(this->name(), "NXlog", 1);
2102 file->writeData("value", values);
2103 file->openData("value");
2104 file->putAttr("units", this->units());
2105 file->closeData();
2106 saveTimeVector(file);
2107 file->closeGroup();
2108}
2109
2114template <typename TYPE> Json::Value TimeSeriesProperty<TYPE>::valueAsJson() const { return Json::Value(value()); }
2115
2126template <typename TYPE>
2127void TimeSeriesProperty<TYPE>::histogramData(const Types::Core::DateAndTime &tMin, const Types::Core::DateAndTime &tMax,
2128 std::vector<double> &counts) const {
2129
2130 size_t nPoints = counts.size();
2131 if (nPoints == 0)
2132 return; // nothing to do
2133
2134 auto t0 = static_cast<double>(tMin.totalNanoseconds());
2135 auto t1 = static_cast<double>(tMax.totalNanoseconds());
2136 if (t0 > t1)
2137 throw std::invalid_argument("invalid arguments for histogramData; tMax<tMin");
2138
2139 double dt = (t1 - t0) / static_cast<double>(nPoints);
2140
2141 for (auto const &ev : m_values) {
2142 auto time = static_cast<double>(ev.time().totalNanoseconds());
2143 if (time < t0 || time >= t1)
2144 continue;
2145 auto ind = static_cast<size_t>((time - t0) / dt);
2146 counts[ind] += static_cast<double>(ev.value());
2147 }
2148}
2149
2150template <>
2151void TimeSeriesProperty<std::string>::histogramData(const Types::Core::DateAndTime &tMin,
2152 const Types::Core::DateAndTime &tMax,
2153 std::vector<double> &counts) const {
2154 UNUSED_ARG(tMin);
2155 UNUSED_ARG(tMax);
2156 UNUSED_ARG(counts);
2157 throw std::runtime_error("histogramData is not implememnted for time series "
2158 "properties containing strings");
2159}
2160
2170template <typename TYPE> std::vector<TYPE> TimeSeriesProperty<TYPE>::filteredValuesAsVector(const TimeROI *roi) const {
2171 if (roi && !roi->useAll()) {
2172 this->sortIfNecessary();
2173 std::vector<TYPE> filteredValues;
2174 if (roi->firstTime() > this->m_values.back().time()) {
2175 // Since the ROI starts after everything, just return the last value in the log
2176 filteredValues.emplace_back(this->m_values.back().value());
2177 } else { // only use the values in the filter - this is very similar to FilteredTimeSeriesProperty::applyFilter
2178 // the index into the m_values array of the time, or -1 (before) or m_values.size() (after)
2179 std::size_t index_current_log{0};
2180 for (const auto &splitter : roi->toTimeIntervals()) {
2181 const auto endTime = splitter.stop();
2182
2183 // check if the splitter starts too early
2184 if (endTime < this->m_values[index_current_log].time()) {
2185 continue; // skip to the next splitter
2186 }
2187
2188 // cache values to reduce number of method calls
2189 const auto beginTime = splitter.start();
2190
2191 // find the first log that should be added
2192 if (this->m_values.back().time() < beginTime) {
2193 // skip directly to the end if the filter starts after the last log
2194 index_current_log = this->m_values.size() - 1;
2195 } else {
2196 // search for the right starting point
2197 while ((this->m_values[index_current_log].time() <= beginTime)) {
2198 if (index_current_log + 1 > this->m_values.size())
2199 break;
2200 index_current_log++;
2201 }
2202 // need to back up by one
2203 if (index_current_log > 0)
2204 index_current_log--;
2205 // go backwards more while times are equal to the one being started at
2206 while (index_current_log > 0 &&
2207 this->m_values[index_current_log].time() == this->m_values[index_current_log - 1].time()) {
2208 index_current_log--;
2209 }
2210 }
2211
2212 // add everything up to the end time
2213 for (; index_current_log < this->m_values.size(); ++index_current_log) {
2214 if (this->m_values[index_current_log].time() >= endTime)
2215 break;
2216
2217 // the current value goes into the filter
2218 filteredValues.emplace_back(this->m_values[index_current_log].value());
2219 }
2220 // go back one so the next splitter can add a value
2221 if (index_current_log > 0)
2222 index_current_log--;
2223 }
2224 }
2225 return filteredValues;
2226 } else {
2227 return this->valuesAsVector();
2228 }
2229}
2230
2231template <typename TYPE> std::vector<TYPE> TimeSeriesProperty<TYPE>::filteredValuesAsVector() const {
2232 return this->valuesAsVector();
2233}
2234
2243template <typename TYPE> std::vector<TimeInterval> TimeSeriesProperty<TYPE>::getTimeIntervals() const {
2244 std::vector<TimeInterval> intervals;
2245 auto lastInterval = this->nthInterval(this->size() - 1);
2246 intervals.emplace_back(firstTime(), lastInterval.stop());
2247 return intervals;
2248}
2249
2251// -------------------------- Macro to instantiation concrete types
2252// --------------------------------
2253#define INSTANTIATE(TYPE) template class TimeSeriesProperty<TYPE>;
2254
2255// -------------------------- Concrete instantiation
2256// -----------------------------------------------
2257INSTANTIATE(int32_t)
2258INSTANTIATE(int64_t)
2259INSTANTIATE(uint32_t)
2260INSTANTIATE(uint64_t)
2261INSTANTIATE(float)
2262INSTANTIATE(double)
2263INSTANTIATE(std::string)
2264INSTANTIATE(bool)
2265
2266
2267
2268} // namespace Kernel
2269} // 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:48
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:61
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:28
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:244
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