Mantid
Loading...
Searching...
No Matches
MDPlane.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 +
9#include "MantidKernel/VMD.h"
10#include <gsl/gsl_linalg.h>
11#include <stdexcept>
12
13using namespace Mantid::Kernel;
14
15namespace Mantid::Geometry {
16
17//----------------------------------------------------------------------------------------------
24MDPlane::MDPlane(const std::vector<coord_t> &normal, const std::vector<coord_t> &point) {
25 m_nd = normal.size();
26 if ((m_nd < 1) || (m_nd > 100))
27 throw std::invalid_argument("MDPlane::ctor(): Invalid number of dimensions in the normal vector !");
28 if (point.size() != normal.size())
29 throw std::invalid_argument("MDPlane::ctor(): Inconsistent number of "
30 "dimensions in the normal/point vectors!");
31 construct(normal, point);
32}
33
34//----------------------------------------------------------------------------------------------
42 m_nd = normal.getNumDims();
43 if ((m_nd < 1) || (m_nd > 100))
44 throw std::invalid_argument("MDPlane::ctor(): Invalid number of dimensions in the normal vector !");
45 if (point.getNumDims() != normal.getNumDims())
46 throw std::invalid_argument("MDPlane::ctor(): Inconsistent number of "
47 "dimensions in the normal/point vectors!");
48 construct(normal, point);
49}
50
51//----------------------------------------------------------------------------------------------
59MDPlane::MDPlane(const size_t nd, const float *normal, const float *point) : m_nd(nd) {
60 if ((m_nd < 1) || (m_nd > 100))
61 throw std::invalid_argument("MDPlane::ctor(): Invalid number of dimensions in the workspace!");
62 construct(normal, point);
63}
64
65//----------------------------------------------------------------------------------------------
73MDPlane::MDPlane(const size_t nd, const double *normal, const double *point) : m_nd(nd) {
74 if ((m_nd < 1) || (m_nd > 100))
75 throw std::invalid_argument("MDPlane::ctor(): Invalid number of dimensions in the workspace!");
76 construct(normal, point);
77}
78
79//----------------------------------------------------------------------------------------------
88MDPlane::MDPlane(const std::vector<Mantid::Kernel::VMD> &vectors, const Mantid::Kernel::VMD &origin,
89 const Mantid::Kernel::VMD &insidePoint) {
90 // Get the normal vector by the determinant method
91 VMD normal = VMD::getNormalVector(vectors);
92
93 // The dimensionality of the plane
94 m_nd = normal.getNumDims();
95
96 // Whew. We have a normal, and a point on the plane. We can construct
97 construct(normal, origin);
98
99 // Did we get the wrong sign of the normal?
100 if (!this->isPointBounded(insidePoint)) {
101 // Flip the normal over
102 delete[] this->m_normal;
103 for (size_t d = 0; d < normal.getNumDims(); d++)
104 normal[d] = -1.0f * normal[d];
105 // And re-construct
106 construct(normal, origin);
107 }
108}
109
110// //----------------------------------------------------------------------------------------------
111// MDPlane::MDPlane(const std::vector<Mantid::Kernel::VMD> & points, const VMD
112// & insidePoint)
113// {
114// if (points.size() <= 0)
115// throw std::invalid_argument("MDPlane::ctor(): Must give at least 1
116// point");
117// VMD origin = points[0];
118// m_nd = origin.getNumDims();
119// if (insidePoint.getNumDims() != m_nd)
120// throw std::invalid_argument("MDPlane::ctor(): The insidePoint parameter
121// must match the dimensions of the other points!");
122// if (m_nd < 1)
123// throw std::invalid_argument("MDPlane::ctor(): Must have at least 1
124// dimension!");
125// if (points.size() != m_nd)
126// throw std::invalid_argument("MDPlane::ctor(): Must have as many points
127// as there are dimensions!");
128// for (size_t i=0; i < points.size(); i++)
129// if (points[i].getNumDims() != m_nd)
130// throw std::invalid_argument("MDPlane::ctor(): Inconsistent number of
131// dimensions in the points given!");
132//
133// // Special case for 1D (no need to use the linear equation system)
134// if (m_nd == 1)
135// {
136// VMD origin = points[0];
137// VMD normal(1);
138// normal[0] = (insidePoint[0] > origin[0]) ? +1.0 : -1.0;
139// construct(normal, origin);
140// return;
141// }
142//
143// // Make the A and B matrices as described above
144// double * a_data = new double[(m_nd-1) * (m_nd-1)];;
145// double * b_data = new double[(m_nd-1)];
146// for (size_t row=0; row<(m_nd-1); row++)
147// {
148// // Offset the point
149// VMD point = points[row+1] - origin;
150// // Put the 1th+ coords
151// for (size_t col=0; col<(m_nd-1); col++)
152// a_data[row*(m_nd-1) + col] = point[col+1];
153// // The row in B = - the 0th coordinate of the offset point.
154// b_data[row] = -point[0];
155// }
156// // To avoid aborting
157// gsl_set_error_handler_off();
158//
159// int ret = 0;
160//
161// gsl_matrix_view m = gsl_matrix_view_array (a_data, m_nd-1, m_nd-1);
162// gsl_vector_view b = gsl_vector_view_array (b_data, m_nd-1);
163// gsl_vector *x = gsl_vector_alloc (m_nd-1);
164//
165// // Use GSL to solve the linear algebra problem
171//
172// gsl_matrix * V = gsl_matrix_alloc(m_nd-1, m_nd-1);
173// gsl_vector * S = gsl_vector_alloc(m_nd-1);
174// gsl_vector * work = gsl_vector_alloc(m_nd-1);
175// ret = gsl_linalg_SV_decomp(&m.matrix, V, S, work);
176// if (ret == 0)
177// ret = gsl_linalg_SV_solve(&m.matrix, V, S, &b.vector, x);
178//
179// if (ret != 0)
180// {
181// std::string pointsStr;
182// for (size_t i=0; i < points.size(); i++)
183// pointsStr += points[i].toString(",") + ((i != points.size()-1) ? "; "
184// : "");
185// // Something failed
186// throw std::runtime_error("MDPlane::ctor(): unable to solve the system of
187// equations.\n"
188// "The points given likely do not form a plane (some may be
189// collinear), meaning the plane cannot be constructed.\n"
190// "Points were: " + pointsStr);
191// }
192//
193// // The normal vector
194// VMD normal(m_nd);
195//
196// // Recall that we fixed a1 to 1.0
197// normal[0] = 1.0;
198// for (size_t d=0; d<(m_nd-1); d++)
199// normal[d+1] = x->data[d];
200//
201// gsl_vector_free (x);
202// gsl_vector_free (S);
203// gsl_matrix_free (V);
204// delete [] a_data;
205// delete [] b_data;
206//
207// // Whew. We have a normal, and a point on the plane. We can construct
208// construct(normal, origin);
209//
210// // Did we get the wrong normal?
211// if (!this->isPointBounded(insidePoint))
212// {
213// // Flip the normal over
214// delete [] this->m_normal;
215// for (size_t d=0; d<normal.getNumDims(); d++)
216// normal[d] = -1.0 * normal[d];
217// // And re-construct
218// construct(normal, origin);
219// }
220// }
221
222//----------------------------------------------------------------------------------------------
227MDPlane::MDPlane(const MDPlane &other) : m_nd(other.m_nd), m_inequality(other.m_inequality) {
228 m_normal = new coord_t[m_nd];
229 for (size_t d = 0; d < m_nd; d++)
230 m_normal[d] = other.m_normal[d];
231}
232
233//----------------------------------------------------------------------------------------------
239 if (this != &other) // protect against invalid self-assignment
240 {
241 m_nd = other.m_nd;
242 m_inequality = other.m_inequality;
243 // Replace the array
244 delete[] m_normal;
245 m_normal = new coord_t[m_nd];
246 for (size_t d = 0; d < m_nd; d++)
247 m_normal[d] = other.m_normal[d];
248 }
249 return *this;
250}
251
252//----------------------------------------------------------------------------------------------
256
257} // namespace Mantid::Geometry
A generalized description of a N-dimensional hyperplane.
Definition MDPlane.h:41
MDPlane(const Mantid::Kernel::VMD &normal, const Mantid::Kernel::VMD &point)
Constructor with normal and point.
Definition MDPlane.cpp:41
coord_t m_inequality
Right-hand side of the linear equation. aka b in : a1*x1 + a2*x2 + ... < b.
Definition MDPlane.h:189
void construct(IterableType normal, IterableType point)
Templated method for constructing the MDPlane.
Definition MDPlane.h:170
size_t m_nd
Number of dimensions in the MDEventWorkspace.
Definition MDPlane.h:181
bool isPointBounded(const coord_t *coords) const
Is a point in MDimensions bounded by this hyperplane, that is, is (a1*x1 + a2*x2 + ....
Definition MDPlane.h:71
coord_t * m_normal
Coefficients to multiply each coordinate, sized m_nd.
Definition MDPlane.h:186
MDPlane & operator=(const MDPlane &other)
Assignment operator.
Definition MDPlane.cpp:238
size_t getNumDims() const
Definition VMD.cpp:236
static VMDBase getNormalVector(const std::vector< VMDBase > &vectors)
Given N-1 vectors defining a N-1 dimensional hyperplane in N dimensions, returns a vector that is nor...
Definition VMD.cpp:494
float coord_t
Typedef for the data type to use for coordinate axes in MD objects such as MDBox, MDEventWorkspace,...
Definition MDTypes.h:27