BornAgain  1.18.0
Simulate and fit neutron and x-ray scattering at grazing incidence
FormFactorTruncatedSpheroid.cpp
Go to the documentation of this file.
1 // ************************************************************************** //
2 //
3 // BornAgain: simulate and fit scattering at grazing incidence
4 //
5 //! @file Sample/HardParticle/FormFactorTruncatedSpheroid.cpp
6 //! @brief Implements class FormFactorTruncatedSpheroid.
7 //!
8 //! @homepage http://www.bornagainproject.org
9 //! @license GNU General Public License v3 or higher (see COPYING)
10 //! @copyright Forschungszentrum Jülich GmbH 2018
11 //! @authors Scientific Computing Group at MLZ (see CITATION, AUTHORS)
12 //
13 // ************************************************************************** //
14 
17 #include "Base/Types/Exceptions.h"
18 #include "Base/Utils/Integrator.h"
21 #include <limits>
22 
23 FormFactorTruncatedSpheroid::FormFactorTruncatedSpheroid(const std::vector<double> P)
25  {"TruncatedSpheroid",
26  "class_tooltip",
27  {{"Radius", "nm", "horizontal radius", 0, +INF, 0},
28  {"Height", "nm", "height before removal of cap", 0, +INF, 0},
29  {"HeightFlattening", "", "ratio of vertical to horizontal radius", 0, +INF, 0},
30  {"DeltaHeight", "nm", "height of removed cap", 0, +INF, 0}}},
31  P),
32  m_radius(m_P[0]), m_height(m_P[1]), m_height_flattening(m_P[2]), m_dh(m_P[3])
33 {
34  check_initialization();
35  onChange();
36 }
37 
38 FormFactorTruncatedSpheroid::FormFactorTruncatedSpheroid(double radius, double height,
39  double height_flattening, double dh)
40  : FormFactorTruncatedSpheroid(std::vector<double>{radius, height, height_flattening, dh})
41 {
42 }
43 
44 bool FormFactorTruncatedSpheroid::check_initialization() const
45 {
46  bool result(true);
47  if (m_height > 2. * m_radius * m_height_flattening || m_dh > m_height) {
48  std::ostringstream ostr;
49  ostr << "::FormFactorTruncatedSpheroid() -> Error in class initialization with parameters ";
50  ostr << " radius:" << m_radius;
51  ostr << " height:" << m_height;
52  ostr << " height_flattening:" << m_height_flattening << "\n\n";
53  ostr << "Check for 'height <= 2.*radius*height_flattening' failed.";
55  }
56  return result;
57 }
58 
59 //! Integrand for complex form factor.
60 complex_t FormFactorTruncatedSpheroid::Integrand(double Z) const
61 {
62  double R = m_radius;
63  double fp = m_height_flattening;
64 
65  double Rz = std::sqrt(R * R - Z * Z / (fp * fp));
66  complex_t qrRz = std::sqrt(m_q.x() * m_q.x() + m_q.y() * m_q.y()) * Rz;
67  complex_t J1_qrRz_div_qrRz = MathFunctions::Bessel_J1c(qrRz);
68 
69  return Rz * Rz * J1_qrRz_div_qrRz * exp_I(m_q.z() * Z);
70 }
71 
73 {
74  double H = m_height;
75  double R = m_radius;
76  double fp = m_height_flattening;
77  m_q = q;
78 
79  if (std::abs(m_q.mag()) <= std::numeric_limits<double>::epsilon())
80  return M_PI / 3. / fp * (H * H * (3. * R - H / fp) - m_dh * m_dh * (3. * R - m_dh / fp));
81  complex_t z_part = exp_I(m_q.z() * (H - fp * R));
82  return M_TWOPI * z_part
83  * ComplexIntegrator().integrate([&](double Z) { return Integrand(Z); }, fp * R - H,
84  fp * R - m_dh);
85 }
86 
88  kvector_t translation) const
89 {
90  double height = m_height - m_dh;
91  auto effects = computeSlicingEffects(limits, translation, height);
92  FormFactorTruncatedSpheroid slicedff(m_radius, height - effects.dz_bottom, m_height_flattening,
93  effects.dz_top + m_dh);
94  return createTransformedFormFactor(slicedff, rot, effects.position);
95 }
96 
98 {
99  mP_shape.reset(
100  new TruncatedEllipsoid(m_radius, m_radius, m_height_flattening * m_radius, m_height, m_dh));
101 }
complex_t exp_I(complex_t z)
Returns exp(I*z), where I is the imaginary unit.
Definition: Complex.h:30
Defines many exception classes in namespace Exceptionss.
Defines class FormFactorTruncatedSpheroid.
Defines classes RealIntegrator, ComplexIntegrator.
Defines M_PI and some more mathematical constants.
Defines namespace MathFunctions.
Defines class TruncatedEllipsoid.
double mag() const
Returns magnitude of the vector.
T z() const
Returns z-component in cartesian coordinate system.
Definition: BasicVector3D.h:68
T y() const
Returns y-component in cartesian coordinate system.
Definition: BasicVector3D.h:66
T x() const
Returns x-component in cartesian coordinate system.
Definition: BasicVector3D.h:64
To integrate a complex function of a real variable.
Definition: Integrator.h:41
IFormFactor * sliceFormFactor(ZLimits limits, const IRotation &rot, kvector_t translation) const override final
Actually slices the form factor or throws an exception.
void onChange() override final
Action to be taken in inherited class when a parameter has changed.
complex_t evaluate_for_q(cvector_t q) const override final
Returns scattering amplitude for complex scattering wavevector q=k_i-k_f.
Pure virtual base class for Born form factors.
SlicingEffects computeSlicingEffects(ZLimits limits, const kvector_t &position, double height) const
Helper method for slicing.
std::unique_ptr< IShape > mP_shape
IShape object, used to retrieve vertices (which may be approximate in the case of round shapes).
Pure virtual base class for all form factors.
Definition: IFormFactor.h:40
Pure virtual interface for rotations.
Definition: Rotations.h:27
Class that contains upper and lower limits of the z-coordinate for the slicing of form factors.
Definition: ZLimits.h:41
double Bessel_J1c(double x)
Bessel function Bessel_J1(x)/x.