BornAgain  1.19.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 reflection and scattering
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 
16 #include "Base/Math/Constants.h"
17 #include "Base/Math/Functions.h"
20 #include <limits>
21 
23  : IBornFF(
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])
30  , m_height(m_P[1])
31 {
32  onChange();
33 }
34 
36  : FormFactorFullSpheroid(std::vector<double>{radius, height})
37 {
38 }
39 
41 {
42  double h = m_height / 2;
43  double R = m_radius;
44 
45  // complex length of q (not a sesquilinear dot product!),
46  // xy components multiplied with R, z component multiplied with h
47  complex_t qR = sqrt(R * R * (q.x() * q.x() + q.y() * q.y()) + h * h * q.z() * q.z());
48 
49  complex_t zFactor = exp_I(h * q.z());
50 
51  if (std::abs(qR) < 1e-4)
52  // expand sin(qR)-qR*cos(qR) up to qR^5
53  return 4 * M_PI / 3 * R * R * h * (1. - 0.1 * pow(qR, 2)) * zFactor;
54 
55  return 4 * M_PI / pow(qR, 3) * R * R * h * (sin(qR) - qR * cos(qR)) * zFactor;
56 }
57 
59  kvector_t translation) const
60 {
61  double flattening = m_height / (2.0 * m_radius);
62  auto effects = computeSlicingEffects(limits, translation, m_height);
63  FormFactorTruncatedSpheroid slicedff(m_radius, m_height - effects.dz_bottom, flattening,
64  effects.dz_top);
65  return createTransformedFormFactor(slicedff, rot, effects.position);
66 }
67 
69 {
70  m_shape3D =
71  std::make_unique<TruncatedEllipsoid>(m_radius, m_radius, m_height / 2.0, m_height, 0.0);
72 }
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_PI
Definition: Constants.h:44
Defines class FormFactorFullSpheroid.
Defines class FormFactorTruncatedSpheroid.
Defines functions in namespace Math.
const double INF
Definition: INode.h:25
Defines class TruncatedEllipsoid.
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 full spheroid (an ellipsoid with two equal axes, hence with circular cross section)
complex_t evaluate_for_q(cvector_t q) const 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 final
Actually slices the form factor or throws an exception.
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
Definition: filesystem.h:81