BornAgain  1.19.0
Simulate and fit neutron and x-ray scattering at grazing incidence
FormFactorCylinder.cpp
Go to the documentation of this file.
1 // ************************************************************************************************
2 //
3 // BornAgain: simulate and fit reflection and scattering
4 //
5 //! @file Sample/HardParticle/FormFactorCylinder.cpp
6 //! @brief Implements class FormFactorCylinder.
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 
21 FormFactorCylinder::FormFactorCylinder(const std::vector<double> P)
22  : IBornFF(
23  {"Cylinder",
24  "circular cylinder",
25  {{"Radius", "nm", "radius of base", 0, +INF, 0}, {"Height", "nm", "height", 0, +INF, 0}}},
26  P)
27  , m_radius(m_P[0])
28  , m_height(m_P[1])
29 {
30  onChange();
31 }
32 
33 FormFactorCylinder::FormFactorCylinder(double radius, double height)
34  : FormFactorCylinder(std::vector<double>{radius, height})
35 {
36 }
37 
39 {
40  double R = m_radius;
41  double H = m_height;
42 
43  complex_t qzH_half = q.z() * H / 2.0;
44  complex_t z_part = H * Math::sinc(qzH_half) * exp_I(qzH_half);
45  complex_t qxy = std::sqrt(q.x() * q.x() + q.y() * q.y());
46  complex_t radial_part = M_TWOPI * R * R * Math::Bessel::J1c(qxy * R);
47  complex_t result = radial_part * z_part;
48 
49  return result;
50 }
51 
53  kvector_t translation) const
54 {
55  auto effects = computeSlicingEffects(limits, translation, m_height);
56  FormFactorCylinder slicedff(m_radius, m_height - effects.dz_bottom - effects.dz_top);
57  return createTransformedFormFactor(slicedff, rot, effects.position);
58 }
59 
61 {
62  m_shape3D = std::make_unique<DoubleEllipse>(m_radius, m_radius, m_height, m_radius, m_radius);
63 }
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 FormFactorCylinder.
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 circular cylinder.
FormFactorCylinder(const std::vector< double > P)
complex_t evaluate_for_q(cvector_t q) const final
Returns scattering amplitude for complex scattering wavevector q=k_i-k_f.
IFormFactor * sliceFormFactor(ZLimits limits, const IRotation &rot, kvector_t translation) const final
Actually slices the form factor or throws an exception.
const double & m_height
const double & m_radius
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