BornAgain  1.18.0
Simulate and fit neutron and x-ray scattering at grazing incidence
FormFactorFullSpheroid.cpp
Go to the documentation of this file.
1 // ************************************************************************** //
2 //
3 // BornAgain: simulate and fit scattering at grazing incidence
4 //
5 //! @file Sample/HardParticle/FormFactorFullSpheroid.cpp
6 //! @brief Implements class FormFactorFullSpheroid.
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 
20 #include <limits>
21 
24  {"FullSpheroid",
25  "ellipsoid of revolution",
26  {{"Radius", "nm", "revolution radius", 0, +INF, 0},
27  {"Height", "nm", "height = twice the radius in non-revolution direction", 0, +INF, 0}}},
28  P),
29  m_radius(m_P[0]), m_height(m_P[1])
30 {
31  onChange();
32 }
33 
35  : FormFactorFullSpheroid(std::vector<double>{radius, height})
36 {
37 }
38 
40 {
41  double h = m_height / 2;
42  double R = m_radius;
43 
44  // complex length of q (not a sesquilinear dot product!),
45  // xy components multiplied with R, z component multiplied with h
46  complex_t qR = sqrt(R * R * (q.x() * q.x() + q.y() * q.y()) + h * h * q.z() * q.z());
47 
48  complex_t zFactor = exp_I(h * q.z());
49 
50  if (std::abs(qR) < 1e-4)
51  // expand sin(qR)-qR*cos(qR) up to qR^5
52  return 4 * M_PI / 3 * R * R * h * (1. - 0.1 * pow(qR, 2)) * zFactor;
53 
54  return 4 * M_PI / pow(qR, 3) * R * R * h * (sin(qR) - qR * cos(qR)) * zFactor;
55 }
56 
58  kvector_t translation) const
59 {
60  double flattening = m_height / (2.0 * m_radius);
61  auto effects = computeSlicingEffects(limits, translation, m_height);
62  FormFactorTruncatedSpheroid slicedff(m_radius, m_height - effects.dz_bottom, flattening,
63  effects.dz_top);
64  return createTransformedFormFactor(slicedff, rot, effects.position);
65 }
66 
68 {
69  mP_shape =
70  std::make_unique<TruncatedEllipsoid>(m_radius, m_radius, m_height / 2.0, m_height, 0.0);
71 }
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 FormFactorFullSpheroid.
Defines class FormFactorTruncatedSpheroid.
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_PI
Definition: MathConstants.h:39
Defines namespace MathFunctions.
Defines class TruncatedEllipsoid.
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 full spheroid (an ellipsoid with two equal axes, hence with circular cross section)
complex_t evaluate_for_q(cvector_t q) const override final
Returns scattering amplitude for complex scattering wavevector q=k_i-k_f.
FormFactorFullSpheroid(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.
void onChange() override final
Action to be taken in inherited class when a parameter has changed.
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
const double radius(5 *Units::nanometer)