Mantid
Loading...
Searching...
No Matches
IntegratePeakTimeSlices.h
Go to the documentation of this file.
1// Mantid Repository : https://github.com/mantidproject/mantid
2//
3// Copyright © 2009 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/*
8 * IntegratePeakTimeSlices.h
9 *
10 * Created on: May 5, 2011
11 * Author: ruth
12 */
13#pragma once
14
15#include "MantidAPI/Algorithm.h"
21#include "MantidCrystal/DllConfig.h"
26#include "MantidKernel/V3D.h"
28
29#include <array>
30
31// Forward Declarations
32namespace Mantid {
33namespace HistogramData {
34class HistogramX;
35}
36} // namespace Mantid
37namespace Mantid {
38namespace Crystal {
49public:
51
52 DataModeHandler(const DataModeHandler &handler);
53 DataModeHandler(double baseRCRadius, double lastRCRadius, double lastRow, double lastCol, double CellWidth,
54 double CellHeight, bool CalcVariance, int MinCol, int MaxCol, int MinRow, int MaxRow) {
55 init();
56 this->baseRCRadius = baseRCRadius;
57 this->lastRCRadius = lastRCRadius;
58 this->lastRow = lastRow;
59 this->lastCol = lastCol;
60 this->CellWidth = CellWidth;
61 this->CellHeight = CellHeight;
62 this->CalcVariance = CalcVariance;
63 this->MaxCol = MaxCol;
64 this->MaxRow = MaxRow;
65 this->MinCol = MinCol;
66 this->MinRow = MinRow;
67 }
68
69 void setTime(double time) { this->time = time; }
70
71 // Returns true if "Edge Peak" otherwise returns false
72 bool setStatBase(std::vector<double> const &StatBase);
73
74 bool isEdgePeak(const double *params, int nparams);
75
76 void setHeightHalfWidthInfo(const MantidVec &xvals, const MantidVec &yvals, const MantidVec &counts);
77
78 void setCurrentRadius(double radius) { currentRadius = radius; }
79 void setCurrentCenter(const Kernel::V3D &newCenter) {
80 Kernel::V3D XX(newCenter);
81 currentPosition = XX;
82 }
83
84 double getCurrentRadius() { return currentRadius; }
86 void updateEdgeXsize(double newsize) {
87 if (EdgeX < 0)
88
89 EdgeX = newsize;
90
91 else if (newsize < EdgeX)
92
93 EdgeX = newsize;
94 }
95
96 void updateEdgeYsize(double newsize) {
97 if (EdgeY < 0)
98
99 EdgeY = newsize;
100
101 else if (newsize < EdgeY)
102
103 EdgeY = newsize;
104 }
105
106 void CalcVariancesFromData(double background, double meanx, double meany, double &Varxx, double &Varxy, double &Varyy,
107 const std::vector<double> &StatBase);
108
109 bool IsEnoughData(const double *ParameterValues, Kernel::Logger &);
110
111 double getNewRCRadius();
112
113 double getInitBackground() { return back_calc; }
114
115 double getInitRow() { return row_calc; }
116
117 double getInitCol() { return col_calc; }
118
119 double getInitIntensity() { return Intensity_calc; }
120
121 double getInitVarx() { return Vx_calc; }
122
123 double getInitVary() { return Vy_calc; }
124
125 double getInitVarxy() { return Vxy_calc; }
126
127 std::string CalcConstraints(std::vector<std::pair<double, double>> &Bounds, bool CalcVariances);
128 std::string getTies() { return ""; }
129
130 bool CalcVariances();
131
132 std::vector<double> GetParams(double b);
133
134 double StatBaseVals(int index) {
135 if (index < 0 || index >= (int)StatBase.size())
136 return 0.0;
137 return StatBase[index];
138 }
139
140 double CalcISAWIntensity(const double *params);
141
142 double CalcISAWIntensityVariance(const double *params, const double *errs, double chiSqOvDOF);
143
147 double CalcSampleIntensityMultiplier(const double *params) const;
148
154 std::vector<double> InitValues(double Varx, double Vary, double b);
159
160 double lastRow;
162 double lastCol;
164 double time;
165 double CellWidth;
167
168 double VarxHW, VaryHW;
171 std::vector<double> StatBase;
172
173 double EdgeX, EdgeY;
176 bool case4; // if true result of successful merge of dir =1 chan=0 and chan=1
178
179private:
180 void init() {
181 baseRCRadius = -1;
182 lastRCRadius = -1;
183 lastRow = -1;
184 lastCol = -1;
185 EdgeX = EdgeY = -1;
186 calcNewRCRadius = -1;
187 MaxRow = -1;
188 MaxCol = -1;
189 MinRow = -1;
190 MinCol = -1;
191 time = -1;
192 CalcVariance = true;
193 CellWidth = CellHeight = 0;
194 currentRadius = -1;
196 lastISAWVariance = -1;
199 case4 = false;
200
201 VarxHW = -1;
202 VaryHW = -1;
204 }
205};
206
208public:
211
213 ~IntegratePeakTimeSlices() override;
214
216 const std::string name() const override { return "IntegratePeakTimeSlices"; }
217
219 const std::string summary() const override {
220 return "The algorithm uses CurveFitting::BivariateNormal for fitting a "
221 "time slice";
222 }
223
225 int version() const override { return 1; }
226 const std::vector<std::string> seeAlso() const override { return {"PeakIntegration"}; }
227
229 const std::string category() const override { return "Crystal\\Integration"; }
230
231private:
232 void init() override;
233 void exec() override;
234
236
237 std::string m_AttributeNames[20];
238
239 std::string m_ParameterNames[7];
240
241 std::shared_ptr<DataModeHandler> m_AttributeValues;
242 std::array<double, 7> m_ParameterValues;
243
245
246 int *m_NeighborIDs; // Stores IDs of nearest neighbors
247 double m_R0;
248
252 double m_ROW;
253 double m_COL;
255 double m_cellWidth;
259
260 void SetUpData(API::MatrixWorkspace_sptr &Data, API::MatrixWorkspace_const_sptr const &inpWkSpace,
261 const std::shared_ptr<Geometry::IComponent> &comp, const int chanMin, const int chanMax, double CentX,
262 double CentY, Kernel::V3D &CentNghbr,
263
264 double &neighborRadius, // from CentDet
265 double Radius, std::string &spec_idList);
266
267 bool getNeighborPixIDs(const std::shared_ptr<Geometry::IComponent> &comp, const Kernel::V3D &Center, double &Radius,
268 int *&ArryofID);
269
270 int CalculateTimeChannelSpan(Geometry::IPeak const &peak, const double dQ, const Mantid::HistogramData::HistogramX &X,
271 const int specNum, int &Centerchan);
272
273 double CalculatePositionSpan(DataObjects::Peak const &peak, const double dQ);
274
275 void InitializeColumnNamesInTableWorkspace(DataObjects::TableWorkspace_sptr &TabWS);
276
279 void SetUpData1(API::MatrixWorkspace_sptr &Data, API::MatrixWorkspace_const_sptr const &inpWkSpace, const int chanMin,
280 const int chanMax, double Radius, const Kernel::V3D &CentPos, std::string &spec_idList
281
282 );
283
287 void PreFit(const API::MatrixWorkspace_sptr &Data, double &chisqOverDOF, bool &done, std::vector<std::string> &names,
288 std::vector<double> &params, std::vector<double> &errs, double lastRow, double lastCol,
289 double neighborRadius);
290
291 void Fit(const API::MatrixWorkspace_sptr &Data, double &chisqOverDOF, bool &done, std::vector<std::string> &names,
292 std::vector<double> &params, std::vector<double> &errs, double lastRow, double lastCol,
293 double neighborRadius);
294
295 std::string CalculateFunctionProperty_Fit();
296
297 bool isGoodFit(std::vector<double> const &params, std::vector<double> const &errs,
298 std::vector<std::string> const &names, double chisqOverDOF);
299 // returns last row added
300 int UpdateOutputWS(DataObjects::TableWorkspace_sptr &TabWS, const int dir, const double chan,
301 std::vector<double> const &params, std::vector<double> const &errs,
302 std::vector<std::string> const &names, const double Chisq, const double time,
303 std::string spec_idList);
304
305 void updatePeakInformation(std::vector<double> const &params, std::vector<double> const &errs,
306 std::vector<std::string> const &names, double &TotVariance, double &TotIntensity,
307 double const TotSliceIntensity, double const TotSliceVariance, double const chisqdivDOF,
308 const int ncells);
309
310 void updateStats(const double intensity, const double variance, const double row, const double col,
311 std::vector<double> &StatBase);
312
313 int findNameInVector(std::string const &oneName, std::vector<std::string> const &nameList);
314
315 double CalculateIsawIntegrateError(const double background, const double backError, const double ChiSqOverDOF,
316 const double TotVariance, const int ncells);
317
318 void FindPlane(Kernel::V3D &center, Kernel::V3D &xvec, Kernel::V3D &yvec, double &ROW, double &COL, int &NROWS,
319 int &NCOLS, double &pixWidthx, double &pixHeighty, DataObjects::Peak const &peak) const;
320
321 int findTimeChannel(const Mantid::HistogramData::HistogramX &X, const double time);
322
323 // returns true if Neighborhood list is changed
324 bool updateNeighbors(const std::shared_ptr<Geometry::IComponent> &comp, const Kernel::V3D &CentPos,
325 const Kernel::V3D &oldCenter, double NewRadius, double &neighborRadius);
326};
327} // namespace Crystal
328} // namespace Mantid
Base class from which all concrete algorithm classes should be derived.
Definition Algorithm.h:76
Class for marking algorithms as deprecated.
Integrates each time slice using the BivariateNormal formula, adding the results to the peak object.
void setHeightHalfWidthInfo(const MantidVec &xvals, const MantidVec &yvals, const MantidVec &counts)
For edge peaks, the sample standard deviations do not work.
double CalcISAWIntensity(const double *params)
Calculates the Intensity designed for Edge Peaks.
double CalcISAWIntensityVariance(const double *params, const double *errs, double chiSqOvDOF)
Calculates the Error in the Intensity designed for Edge Peaks.
std::vector< double > InitValues(double Varx, double Vary, double b)
Returns init values with background and variances replaced by arguments.
bool isEdgePeak(const double *params, int nparams)
Determines if a Peak is an edge peak.
bool CalcVariances()
Determines whether the Variances can be calculated.
double getNewRCRadius()
Calculates the new radius for neighborhoods so as to include almost all of a peak.
void CalcVariancesFromData(double background, double meanx, double meany, double &Varxx, double &Varxy, double &Varyy, const std::vector< double > &StatBase)
Utility method to calculate variances from data given background and means.
double CalcSampleIntensityMultiplier(const double *params) const
For Edge Peaks.
std::vector< double > GetParams(double b)
Calculates the initial values of the parameters given background b.
DataModeHandler(double baseRCRadius, double lastRCRadius, double lastRow, double lastCol, double CellWidth, double CellHeight, bool CalcVariance, int MinCol, int MaxCol, int MinRow, int MaxRow)
bool IsEnoughData(const double *ParameterValues, Kernel::Logger &)
Calculates if there is enough data to for there to be a peak.
std::string CalcConstraints(std::vector< std::pair< double, double > > &Bounds, bool CalcVariances)
Calculates the string form of the constraints to be sent to the Fit Algorithm and also saves it in a ...
bool setStatBase(std::vector< double > const &StatBase)
Sets the Accumulated data values into this class, then updates other information like initial values.
void setCurrentCenter(const Kernel::V3D &newCenter)
double m_COL
for Describing the Column(or 0) describing the center of the
Kernel::V3D m_center
for Describing the Plane at the Peak
Kernel::V3D m_yvec
for Describing the Plane at the Peak
Kernel::V3D m_xvec
for Describing the Plane at the Peak
std::shared_ptr< DataModeHandler > m_AttributeValues
const std::vector< std::string > seeAlso() const override
Function to return all of the seeAlso (these are not validated) algorithms related to this algorithm....
const std::string category() const override
Algorithm's category for identification overriding a virtual method.
double m_ROW
for Describing the Row(or 0) describing the center of the Peak
double m_cellHeight
for Describing the Plane at the Peak
int version() const override
Algorithm's version for identification overriding a virtual method.
double m_R0
for Weak Peaks, these can be set using info from close
const std::string summary() const override
Summary of algorithms purpose.
const std::string name() const override
Algorithm's name for identification overriding a virtual method.
A generic fitting algorithm.
Definition Fit.h:78
Structure describing a single-crystal peak.
Definition Peak.h:34
Structure describing a single-crystal peak.
Definition IPeak.h:26
The Logger class is in charge of the publishing messages from the framework through various channels.
Definition Logger.h:51
Class for 3D vectors.
Definition V3D.h:34
std::shared_ptr< const MatrixWorkspace > MatrixWorkspace_const_sptr
shared pointer to the matrix workspace base class (const version)
std::shared_ptr< MatrixWorkspace > MatrixWorkspace_sptr
shared pointer to the matrix workspace base class
std::shared_ptr< TableWorkspace > TableWorkspace_sptr
shared pointer to Mantid::DataObjects::TableWorkspace
Helper class which provides the Collimation Length for SANS instruments.
std::unordered_map< detid_t, size_t > detid2index_map
Map with key = detector ID, value = workspace index.
std::vector< double > MantidVec
typedef for the data storage used in Mantid matrix workspaces
Definition cow_ptr.h:172
Peak indexing algorithm, which works by assigning multiple possible HKL values to each peak and then ...