BornAgain  1.19.0
Simulate and fit neutron and x-ray scattering at grazing incidence
FormFactorEllipsoidalCylinder.cpp
Go to the documentation of this file.
1 // ************************************************************************************************
2 //
3 // BornAgain: simulate and fit reflection and scattering
4 //
5 //! @file Sample/HardParticle/FormFactorEllipsoidalCylinder.cpp
6 //! @brief Implements class FormFactorEllipsoidalCylinder.
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/Functions.h"
20 
22  : IBornFF({"EllipsoidalCylinder",
23  "elliptical cylinder",
24  {{"RadiusX", "nm", "radius in x direction", 0, +INF, 0},
25  {"RadiusY", "nm", "radius in y direction", 0, +INF, 0},
26  {"Height", "nm", "height", 0, +INF, 0}}},
27  P)
28  , m_radius_x(m_P[0])
29  , m_radius_y(m_P[1])
30  , m_height(m_P[2])
31 {
32  onChange();
33 }
34 
36  double height)
37  : FormFactorEllipsoidalCylinder(std::vector<double>{radius_x, radius_y, height})
38 {
39 }
40 
42 {
43  return (m_radius_x + m_radius_y) / 2.0;
44 }
45 
47 {
48  complex_t qxRa = q.x() * m_radius_x;
49  complex_t qyRb = q.y() * m_radius_y;
50  complex_t qzHdiv2 = m_height / 2 * q.z();
51 
52  complex_t Fz = exp_I(qzHdiv2) * Math::sinc(qzHdiv2);
53  complex_t gamma = std::sqrt((qxRa) * (qxRa) + (qyRb) * (qyRb));
54  complex_t J1_gamma_div_gamma = Math::Bessel::J1c(gamma);
55 
56  return M_TWOPI * m_radius_x * m_radius_y * m_height * Fz * J1_gamma_div_gamma;
57 }
58 
60  kvector_t translation) const
61 {
62  auto effects = computeSlicingEffects(limits, translation, m_height);
64  m_height - effects.dz_bottom - effects.dz_top);
65  return createTransformedFormFactor(slicedff, rot, effects.position);
66 }
67 
69 {
70  m_shape3D =
71  std::make_unique<DoubleEllipse>(m_radius_x, m_radius_y, m_height, m_radius_x, m_radius_y);
72 }
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 DoubleEllipse.
Defines class FormFactorEllipsoidalCylinder.
Defines functions in namespace Math.
const double INF
Definition: INode.h:25
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
A cylinder with elliptical base.
double radialExtension() const final
Returns the (approximate in some cases) radial size of the particle of this form factor's shape.
IFormFactor * sliceFormFactor(ZLimits limits, const IRotation &rot, kvector_t translation) const final
Actually slices the form factor or throws an exception.
complex_t evaluate_for_q(cvector_t q) const final
Returns scattering amplitude for complex scattering wavevector q=k_i-k_f.
FormFactorEllipsoidalCylinder(const std::vector< double > P)
void onChange() final
Action to be taken in inherited class when a parameter has changed.
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
static SlicingEffects computeSlicingEffects(ZLimits limits, const kvector_t &position, double height)
Helper method for slicing.
Definition: IBornFF.cpp:66
Abstract base class for all form factors.
Definition: IFormFactor.h:36
static IFormFactor * createTransformedFormFactor(const IFormFactor &formfactor, const IRotation &rot, kvector_t translation)
Definition: IFormFactor.cpp:99
Abstract base class for rotations.
Definition: Rotations.h:28
Class that contains upper and lower limits of the z-coordinate for the slicing of form factors.
Definition: ZLimits.h:45
double J1c(double x)
Bessel function J1(x)/x.
Definition: Bessel.cpp:172
double sinc(double x)
sinc function:
Definition: Functions.cpp:53
Definition: filesystem.h:81