BornAgain  1.19.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 reflection and scattering
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 
16 #include "Base/Math/Bessel.h"
17 #include "Base/Math/Constants.h"
18 #include "Base/Math/IntegratorGK.h"
20 #include <limits>
21 
23  : IBornFF({"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])
30  , m_radius_y(m_P[1])
31  , m_height(m_P[2])
32 {
33  onChange();
34 }
35 
36 FormFactorHemiEllipsoid::FormFactorHemiEllipsoid(double radius_x, double radius_y, double height)
37  : FormFactorHemiEllipsoid(std::vector<double>{radius_x, radius_y, height})
38 {
39 }
40 
42 {
43  return (m_radius_x + m_radius_y) / 2.0;
44 }
45 
46 //! Integrand for complex form factor.
48 {
49  double R = m_radius_x;
50  double W = m_radius_y;
51  double H = m_height;
52 
53  double Rz = R * std::sqrt(1.0 - Z * Z / (H * H));
54  double Wz = W * std::sqrt(1.0 - Z * Z / (H * H));
55 
56  complex_t qxRz = m_q.x() * Rz;
57  complex_t qyWz = m_q.y() * Wz;
58 
59  complex_t gamma = std::sqrt(qxRz * qxRz + qyWz * qyWz);
60  complex_t J1_gamma_div_gamma = Math::Bessel::J1c(gamma);
61 
62  return Rz * Wz * J1_gamma_div_gamma * exp_I(m_q.z() * Z);
63 }
64 
66 {
67  m_q = q;
68  double R = m_radius_x;
69  double W = m_radius_y;
70  double H = m_height;
71 
72  if (std::abs(m_q.mag()) <= std::numeric_limits<double>::epsilon())
73  return M_TWOPI * R * W * H / 3.;
74  return M_TWOPI * ComplexIntegrator().integrate([&](double Z) { return Integrand(Z); }, 0., H);
75 }
76 
78 {
79  m_shape3D =
80  std::make_unique<TruncatedEllipsoid>(m_radius_x, m_radius_x, m_height, m_height, 0.0);
81 }
Defines Bessel functions in namespace Math.
std::complex< double > complex_t
Definition: Complex.h:20
complex_t exp_I(complex_t z)
Returns exp(I*z), where I is the imaginary unit.
Definition: Complex.h:30
Defines M_PI and some more mathematical constants.
#define M_TWOPI
Definition: Constants.h:54
Defines class FormFactorHemiEllipsoid.
const double INF
Definition: INode.h:25
Defines classes RealIntegrator, ComplexIntegrator.
Defines class TruncatedEllipsoid.
double mag() const
Returns magnitude of the vector.
T z() const
Returns z-component in cartesian coordinate system.
Definition: BasicVector3D.h:67
T y() const
Returns y-component in cartesian coordinate system.
Definition: BasicVector3D.h:65
T x() const
Returns x-component in cartesian coordinate system.
Definition: BasicVector3D.h:63
To integrate a complex function of a real variable.
Definition: IntegratorGK.h:44
complex_t integrate(const std::function< complex_t(double)> &f, double lmin, double lmax)
An hemi ellipsoid, obtained by truncating a full ellipsoid in the middle plane spanned by two princip...
double radialExtension() const final
Returns the (approximate in some cases) radial size of the particle of this form factor's shape.
void onChange() final
Action to be taken in inherited class when a parameter has changed.
complex_t Integrand(double Z) const
Integrand for complex form factor.
complex_t evaluate_for_q(cvector_t q) const final
Returns scattering amplitude for complex scattering wavevector q=k_i-k_f.
FormFactorHemiEllipsoid(const std::vector< double > P)
Abstract base class for Born form factors.
Definition: IBornFF.h:41
std::unique_ptr< IShape3D > m_shape3D
IShape3D object, used to retrieve vertices (which may be approximate in the case of round shapes).
Definition: IBornFF.h:77
double J1c(double x)
Bessel function J1(x)/x.
Definition: Bessel.cpp:172
Definition: filesystem.h:81