Mantid
Loading...
Searching...
No Matches
EigenVector.h
Go to the documentation of this file.
1// Mantid Repository : https://github.com/mantidproject/mantid
2//
3// Copyright © 2010 ISIS Rutherford Appleton Laboratory UKRI,
4// NScD Oak Ridge National Laboratory, European Spallation Source,
5// Institut Laue - Langevin & CSNS, Institute of High Energy Physics, CAS
6// SPDX - License - Identifier: GPL - 3.0 +
7#pragma once
8
9#include "MantidCurveFitting/DllConfig.h"
10
11#include "Eigen/Core"
12#include "Eigen/Dense"
14
15#include <ostream>
16#include <vector>
17
18namespace Mantid {
19namespace CurveFitting {
27class MANTID_CURVEFITTING_DLL EigenVector {
28public:
32 explicit EigenVector(const size_t n);
34 explicit EigenVector(std::vector<double> v);
36 explicit EigenVector(const Eigen::VectorXd *v);
38 EigenVector(std::initializer_list<double> ilist);
40 EigenVector(const EigenVector &v);
42 EigenVector(EigenVector &&v) noexcept;
44 EigenVector &operator=(const EigenVector &v);
46 EigenVector &operator=(const std::vector<double> &v);
48 EigenVector &operator=(const Eigen::VectorXd v);
49
51 inline vec_map_type &mutator() { return m_view.vector_mutator(); }
53 const vec_map_type inspector() const { return m_view.vector_inspector(); }
55 vec_map_type copy_view() const { return m_view.vector_copy(); }
56
58 void resize(const size_t n);
60 size_t size() const;
61
63 void set(const size_t i, const double value);
65 double get(const size_t i) const;
67 const double &operator[](const size_t i) const { return m_data[i]; }
69 double &operator[](const size_t i) { return m_data[i]; }
70 // Set all elements to zero
71 void zero();
73 void normalize();
75 double norm() const;
77 double norm2() const;
79 double dot(const EigenVector &v) const;
81 size_t indexOfMinElement() const;
83 size_t indexOfMaxElement() const;
85 std::pair<size_t, size_t> indicesOfMinMaxElements() const;
87 std::vector<size_t> sortIndices(bool ascending = true) const;
89 void sort(const std::vector<size_t> &indices);
91 std::vector<double> toStdVector() const;
93 std::vector<double> &StdVectorRef();
94
102 EigenVector &operator*=(const double d);
104 EigenVector &operator+=(const double d);
105
106private:
108 std::vector<double> m_data;
111};
112
114MANTID_CURVEFITTING_DLL std::ostream &operator<<(std::ostream &ostr, const EigenVector &v);
115
116} // namespace CurveFitting
117} // namespace Mantid
double value
The value of the point.
Definition: FitMW.cpp:51
const std::vector< Type > & m_data
Definition: TableColumn.h:417
A wrapper around Eigen::Vector.
Definition: EigenVector.h:27
std::vector< double > m_data
Default element storage.
Definition: EigenVector.h:108
vec_map_type & mutator()
Get the map of the eigen vector.
Definition: EigenVector.h:51
double & operator[](const size_t i)
Get a reference to an element.
Definition: EigenVector.h:69
EigenVector_View m_view
The map to the Eigen vector.
Definition: EigenVector.h:110
const double & operator[](const size_t i) const
Get a const reference to an element.
Definition: EigenVector.h:67
vec_map_type copy_view() const
Get a copy of the eigen vector.
Definition: EigenVector.h:55
const vec_map_type inspector() const
Get the const map of the eigen vector.
Definition: EigenVector.h:53
MatrixWorkspace_sptr MANTID_API_DLL operator+=(const MatrixWorkspace_sptr &lhs, const MatrixWorkspace_sptr &rhs)
Adds two workspaces.
MatrixWorkspace_sptr MANTID_API_DLL operator-=(const MatrixWorkspace_sptr &lhs, const MatrixWorkspace_sptr &rhs)
Subtracts two workspaces.
MANTID_API_DLL std::ostream & operator<<(std::ostream &, const AlgorithmHistory &)
Prints a text representation.
MatrixWorkspace_sptr MANTID_API_DLL operator*=(const MatrixWorkspace_sptr &lhs, const MatrixWorkspace_sptr &rhs)
Multiply two workspaces.
Eigen::Map< Eigen::VectorXd, 0, dynamic_stride > vec_map_type
MANTID_KERNEL_DLL V3D normalize(V3D v)
Normalizes a V3D.
Definition: V3D.h:341
Helper class which provides the Collimation Length for SANS instruments.