BornAgain  1.18.0
Simulate and fit neutron and x-ray scattering at grazing incidence
FormFactorHemiEllipsoid.cpp
Go to the documentation of this file.
1 // ************************************************************************** //
2 //
3 // BornAgain: simulate and fit scattering at grazing incidence
4 //
5 //! @file Sample/HardParticle/FormFactorHemiEllipsoid.cpp
6 //! @brief Implements class FormFactorHemiEllipsoid.
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/Utils/Integrator.h"
20 #include <limits>
21 
22 FormFactorHemiEllipsoid::FormFactorHemiEllipsoid(const std::vector<double> P)
23  : IFormFactorBorn({"HemiEllipsoid",
24  "actually a spheroid, truncated at central xy plane",
25  {{"RadiusX", "nm", "radius in x direction", 0, +INF, 0},
26  {"RadiusY", "nm", "radius in y direction", 0, +INF, 0},
27  {"Height", "nm", "height = radius in z direction", 0, +INF, 0}}},
28  P),
29  m_radius_x(m_P[0]), m_radius_y(m_P[1]), m_height(m_P[2])
30 {
31  onChange();
32 }
33 
34 FormFactorHemiEllipsoid::FormFactorHemiEllipsoid(double radius_x, double radius_y, double height)
35  : FormFactorHemiEllipsoid(std::vector<double>{radius_x, radius_y, height})
36 {
37 }
38 
40 {
41  return (m_radius_x + m_radius_y) / 2.0;
42 }
43 
44 //! Integrand for complex form factor.
45 complex_t FormFactorHemiEllipsoid::Integrand(double Z) const
46 {
47  double R = m_radius_x;
48  double W = m_radius_y;
49  double H = m_height;
50 
51  double Rz = R * std::sqrt(1.0 - Z * Z / (H * H));
52  double Wz = W * std::sqrt(1.0 - Z * Z / (H * H));
53 
54  complex_t qxRz = m_q.x() * Rz;
55  complex_t qyWz = m_q.y() * Wz;
56 
57  complex_t gamma = std::sqrt(qxRz * qxRz + qyWz * qyWz);
58  complex_t J1_gamma_div_gamma = MathFunctions::Bessel_J1c(gamma);
59 
60  return Rz * Wz * J1_gamma_div_gamma * exp_I(m_q.z() * Z);
61 }
62 
64 {
65  m_q = q;
66  double R = m_radius_x;
67  double W = m_radius_y;
68  double H = m_height;
69 
70  if (std::abs(m_q.mag()) <= std::numeric_limits<double>::epsilon())
71  return M_TWOPI * R * W * H / 3.;
72  return M_TWOPI * ComplexIntegrator().integrate([&](double Z) { return Integrand(Z); }, 0., H);
73 }
74 
76 {
77  mP_shape =
78  std::make_unique<TruncatedEllipsoid>(m_radius_x, m_radius_x, m_height, m_height, 0.0);
79 }
complex_t exp_I(complex_t z)
Returns exp(I*z), where I is the imaginary unit.
Definition: Complex.h:30
Defines class FormFactorHemiEllipsoid.
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
An hemi ellipsoid, obtained by truncating a full ellipsoid in the middle plane spanned by two princip...
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.
double radialExtension() const override final
Returns the (approximate in some cases) radial size of the particle of this form factor's shape.
Pure virtual base class for Born form factors.
std::unique_ptr< IShape > mP_shape
IShape object, used to retrieve vertices (which may be approximate in the case of round shapes).
double Bessel_J1c(double x)
Bessel function Bessel_J1(x)/x.