Mantid
Loading...
Searching...
No Matches
VesuvioCalculateMS.h
Go to the documentation of this file.
1// Mantid Repository : https://github.com/mantidproject/mantid
2//
3// Copyright © 2014 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// Includes
10//-----------------------------------------------------------------------------
11#include "MantidAPI/Algorithm.h"
12#include "MantidCurveFitting/DllConfig.h"
15#include "MantidKernel/V3D.h"
16
17namespace Mantid {
18
19namespace API {
20class ISpectrum;
21}
22
23namespace Geometry {
24class IObject;
25}
26
27namespace CurveFitting {
28namespace Functions {
29struct ResolutionParams;
30}
31
32namespace Algorithms {
33struct DetectorParams;
34
39class MANTID_CURVEFITTING_DLL VesuvioCalculateMS final : public API::Algorithm {
40private:
41 // Holds date on the compton scattering properties of an atom
43 ComptonNeutronAtom() : mass(-1.0), sclength(-1.0), profile(-1.0) {}
44 double mass; // in amu
45 double sclength; // 4pi/xsec
46 double profile; // s.d of J(y)
47 };
48 // Holds data about sample as a whole.
50 SampleComptonProperties(const int nprops) : atoms(nprops), density(-1.0), totalxsec(-1.0), mu(-1.0) {}
51
52 std::vector<ComptonNeutronAtom> atoms;
53 double density; // g/cm^3
54 double totalxsec; // total free-scattering cross section
55 double mu; // attenuation factor (1/m)
56 };
57
58public:
61 const std::string name() const override { return "VesuvioCalculateMS"; }
63 int version() const override { return 1; }
65 const std::string category() const override { return "CorrectionFunctions\\SpecialCorrections"; }
67 const std::string summary() const override {
68 return "Calculates the contributions of multiple scattering "
69 "on a flat plate sample for VESUVIO";
70 }
71
72 const std::vector<std::string> seeAlso() const override {
73 return {"MayersSampleCorrection", "MonteCarloAbsorption", "CarpenterSampleCorrection"};
74 }
75
76private:
77 void init() override;
78 void exec() override;
79
80 void cacheInputs();
81
82 void calculateMS(CurveFitting::MSVesuvioHelper::RandomVariateGenerator &rng, const size_t wsIndex,
83 API::ISpectrum &totalsc, API::ISpectrum &multsc) const;
85 const Functions::ResolutionParams &respar, MSVesuvioHelper::Simulation &simulCounts) const;
86 void assignToOutput(const MSVesuvioHelper::SimulationWithErrors &avgCounts, API::ISpectrum &totalsc,
87 API::ISpectrum &multsc) const;
88 double calculateCounts(CurveFitting::MSVesuvioHelper::RandomVariateGenerator &rng, const DetectorParams &detpar,
89 const Functions::ResolutionParams &respar, MSVesuvioHelper::Simulation &simulation) const;
90
91 // single-event helpers
92 Kernel::V3D generateSrcPos(CurveFitting::MSVesuvioHelper::RandomVariateGenerator &rng, const double l1) const;
93 double generateE0(CurveFitting::MSVesuvioHelper::RandomVariateGenerator &rng, const double l1, const double t2,
94 double &weight) const;
95 double generateTOF(CurveFitting::MSVesuvioHelper::RandomVariateGenerator &rng, const double en0, const double dtof,
96 const double dl1) const;
97
98 bool generateScatter(CurveFitting::MSVesuvioHelper::RandomVariateGenerator &rng, const Kernel::V3D &startPos,
99 const Kernel::V3D &direc, double &weight, Kernel::V3D &scatterPt) const;
100 std::pair<double, double> calculateE1Range(const double theta, const double en0) const;
101 double partialDiffXSec(const double en0, const double en1, const double theta) const;
103 const Kernel::V3D &nominalPos, const double energy, const Kernel::V3D &scatterPt,
104 const Kernel::V3D &direcBeforeSc, double &scang, double &distToExit) const;
105 double generateE1(CurveFitting::MSVesuvioHelper::RandomVariateGenerator &rng, const double angle, const double e1nom,
106 const double e1res) const;
107
108 // Member Variables
109 size_t m_acrossIdx, m_upIdx, m_beamIdx; // indices of each direction
110 Kernel::V3D m_beamDir; // Directional vector for beam
111 double m_srcR2; // beam penumbra radius (m)
112 double m_halfSampleHeight, m_halfSampleWidth, m_halfSampleThick; // (m)
113 Geometry::IObject const *m_sampleShape; // sample shape
114 std::unique_ptr<SampleComptonProperties> m_sampleProps; // description of sample properties
115 double m_detHeight, m_detWidth, m_detThick; // (m)
116 double m_tmin, m_tmax, m_delt; // min, max & dt TOF value
117 double m_foilRes; // resolution in energy of foil
118
119 size_t m_nscatters; // highest order of scattering to generate
120 size_t m_nruns; // number of runs per spectrum
121 size_t m_nevents; // number of single events per run
122
123 std::unique_ptr<API::Progress> m_progress;
125};
126
127} // namespace Algorithms
128} // namespace CurveFitting
129} // namespace Mantid
double energy
Definition: GetAllEi.cpp:157
Base class from which all concrete algorithm classes should be derived.
Definition: Algorithm.h:85
A "spectrum" is an object that holds the data for a particular spectrum, in particular:
Definition: ISpectrum.h:39
Calculates the multiple scattering & total scattering contributions for a flat-plate or cylindrical s...
int version() const override
function to return a version of the algorithm, must be overridden in all algorithms
const std::string name() const override
function to return a name of the algorithm, must be overridden in all algorithms
std::unique_ptr< SampleComptonProperties > m_sampleProps
const std::string category() const override
function to return a category of the algorithm.
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 summary() const override
function returns a summary message that will be displayed in the default GUI, and in the help.
IObject : Interface for geometry objects.
Definition: IObject.h:41
Class for 3D vectors.
Definition: V3D.h:34
std::shared_ptr< MatrixWorkspace > MatrixWorkspace_sptr
shared pointer to the matrix workspace base class
Helper class which provides the Collimation Length for SANS instruments.
Simple data structure to store nominal detector values It avoids some functions taking a huge number ...
Simple data structure to store resolution parameter values It avoids some functions taking a huge num...