BornAgain  1.18.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 scattering at grazing incidence
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 
19 
21  : IFormFactorBorn({"EllipsoidalCylinder",
22  "elliptical cylinder",
23  {{"RadiusX", "nm", "radius in x direction", 0, +INF, 0},
24  {"RadiusY", "nm", "radius in y direction", 0, +INF, 0},
25  {"Height", "nm", "height", 0, +INF, 0}}},
26  P),
27  m_radius_x(m_P[0]), m_radius_y(m_P[1]), m_height(m_P[2])
28 {
29  onChange();
30 }
31 
33  double height)
34  : FormFactorEllipsoidalCylinder(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  complex_t qxRa = q.x() * m_radius_x;
46  complex_t qyRb = q.y() * m_radius_y;
47  complex_t qzHdiv2 = m_height / 2 * q.z();
48 
49  complex_t Fz = exp_I(qzHdiv2) * MathFunctions::sinc(qzHdiv2);
50  complex_t gamma = std::sqrt((qxRa) * (qxRa) + (qyRb) * (qyRb));
51  complex_t J1_gamma_div_gamma = MathFunctions::Bessel_J1c(gamma);
52 
53  return M_TWOPI * m_radius_x * m_radius_y * m_height * Fz * J1_gamma_div_gamma;
54 }
55 
57  kvector_t translation) const
58 {
59  auto effects = computeSlicingEffects(limits, translation, m_height);
61  m_height - effects.dz_bottom - effects.dz_top);
62  return createTransformedFormFactor(slicedff, rot, effects.position);
63 }
64 
66 {
67  mP_shape =
68  std::make_unique<DoubleEllipse>(m_radius_x, m_radius_y, m_height, m_radius_x, m_radius_y);
69 }
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 class DoubleEllipse.
Defines class FormFactorEllipsoidalCylinder.
IFormFactor * createTransformedFormFactor(const IFormFactor &formfactor, const IRotation &rot, kvector_t translation)
Definition: IFormFactor.cpp:77
const double INF
Definition: INode.h:24
Defines M_PI and some more mathematical constants.
#define M_TWOPI
Definition: MathConstants.h:49
Defines namespace MathFunctions.
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
A cylinder with elliptical base.
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.
double radialExtension() const override final
Returns the (approximate in some cases) radial size of the particle of this form factor's shape.
complex_t evaluate_for_q(cvector_t q) const override final
Returns scattering amplitude for complex scattering wavevector q=k_i-k_f.
FormFactorEllipsoidalCylinder(const std::vector< double > P)
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 sinc(double x)
sinc function:
double Bessel_J1c(double x)
Bessel function Bessel_J1(x)/x.