BornAgain  1.18.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 scattering at grazing incidence
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 
19 
20 FormFactorCylinder::FormFactorCylinder(const std::vector<double> P)
22  {"Cylinder",
23  "circular cylinder",
24  {{"Radius", "nm", "radius of base", 0, +INF, 0}, {"Height", "nm", "height", 0, +INF, 0}}},
25  P),
26  m_radius(m_P[0]), m_height(m_P[1])
27 {
28  onChange();
29 }
30 
32  : FormFactorCylinder(std::vector<double>{radius, height})
33 {
34 }
35 
37 {
38  double R = m_radius;
39  double H = m_height;
40 
41  complex_t qzH_half = q.z() * H / 2.0;
42  complex_t z_part = H * MathFunctions::sinc(qzH_half) * exp_I(qzH_half);
43  complex_t qxy = std::sqrt(q.x() * q.x() + q.y() * q.y());
44  complex_t radial_part = M_TWOPI * R * R * MathFunctions::Bessel_J1c(qxy * R);
45  complex_t result = radial_part * z_part;
46 
47  return result;
48 }
49 
51  kvector_t translation) const
52 {
53  auto effects = computeSlicingEffects(limits, translation, m_height);
54  FormFactorCylinder slicedff(m_radius, m_height - effects.dz_bottom - effects.dz_top);
55  return createTransformedFormFactor(slicedff, rot, effects.position);
56 }
57 
59 {
60  mP_shape = std::make_unique<DoubleEllipse>(m_radius, m_radius, m_height, m_radius, m_radius);
61 }
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 FormFactorCylinder.
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 circular cylinder.
FormFactorCylinder(const std::vector< double > P)
IFormFactor * sliceFormFactor(ZLimits limits, const IRotation &rot, kvector_t translation) const override final
Actually slices the form factor or throws an exception.
const double & m_height
complex_t evaluate_for_q(cvector_t q) const override final
Returns scattering amplitude for complex scattering wavevector q=k_i-k_f.
void onChange() override final
Action to be taken in inherited class when a parameter has changed.
const double & m_radius
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.
const double radius(5 *Units::nanometer)