BornAgain  1.19.79
Open-source research software to simulate and fit neutron and x-ray reflectometry and grazing-incidence small-angle scattering
HemiEllipsoid.cpp
Go to the documentation of this file.
1 // ************************************************************************************************
2 //
3 // BornAgain: simulate and fit reflection and scattering
4 //
5 //! @file Sample/HardParticle/HemiEllipsoid.cpp
6 //! @brief Implements class HemiEllipsoid.
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 
22 HemiEllipsoid::HemiEllipsoid(const std::vector<double> P)
23  : IFormFactor(P)
24  , m_radius_x(m_P[0])
25  , m_radius_y(m_P[1])
26  , m_height(m_P[2])
27 {
28  checkNodeArgs();
29  m_shape3D =
30  std::make_unique<TruncatedEllipsoidNet>(m_radius_x, m_radius_x, m_height, m_height, 0.0);
31 }
32 
33 HemiEllipsoid::HemiEllipsoid(double radius_x, double radius_y, double height)
34  : HemiEllipsoid(std::vector<double>{radius_x, radius_y, height})
35 {
36 }
37 
39 {
40  return (m_radius_x + m_radius_y) / 2.0;
41 }
42 
44 {
45  double R = m_radius_x;
46  double W = m_radius_y;
47  double H = m_height;
48 
49  if (std::abs(q.mag()) <= std::numeric_limits<double>::epsilon())
50  return M_TWOPI * R * W * H / 3.;
51 
52  return M_TWOPI
54  [=](double Z) {
55  double R = m_radius_x;
56  double W = m_radius_y;
57  double H = m_height;
58 
59  double Rz = R * std::sqrt(1.0 - Z * Z / (H * H));
60  double Wz = W * std::sqrt(1.0 - Z * Z / (H * H));
61 
62  complex_t qxRz = q.x() * Rz;
63  complex_t qyWz = q.y() * Wz;
64 
65  complex_t gamma = std::sqrt(qxRz * qxRz + qyWz * qyWz);
66  complex_t J1_gamma_div_gamma = Math::Bessel::J1c(gamma);
67 
68  return Rz * Wz * J1_gamma_div_gamma * exp_I(q.z() * Z);
69  },
70  0., H);
71 }
Defines Bessel functions in namespace Math.
Defines M_PI and some more mathematical constants.
#define M_TWOPI
Definition: Constants.h:54
Defines class HemiEllipsoid.
Defines classes RealIntegrator, ComplexIntegrator.
Defines class TruncatedEllipsoidNet.
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...
Definition: HemiEllipsoid.h:24
const double & m_height
Definition: HemiEllipsoid.h:55
const double & m_radius_x
Definition: HemiEllipsoid.h:53
complex_t formfactor_at_bottom(C3 q) const override
double height() const
Definition: HemiEllipsoid.h:44
const double & m_radius_y
Definition: HemiEllipsoid.h:54
double radialExtension() const override
Returns the (approximate in some cases) radial size of the particle of this form factor's shape....
HemiEllipsoid(double radius_x, double radius_y, double height)
Abstract base class for Born form factors.
Definition: IFormFactor.h:36
std::unique_ptr< IShape3D > m_shape3D
IShape3D object, used to retrieve vertices (which may be approximate in the case of round shapes)....
Definition: IFormFactor.h:74
void checkNodeArgs() const
Raises exception if a parameter value is invalid.
Definition: INode.cpp:27
double J1c(double x)
Bessel function J1(x)/x.
Definition: Bessel.cpp:172
double gamma(double x)
constexpr Double_t H()
Definition: TMath.h:140
constexpr Double_t R()
Definition: TMath.h:213