BornAgain  1.19.0
Simulate and fit neutron and x-ray scattering at grazing incidence
FormFactorTruncatedSphere.cpp
Go to the documentation of this file.
1 // ************************************************************************************************
2 //
3 // BornAgain: simulate and fit reflection and scattering
4 //
5 //! @file Sample/HardParticle/FormFactorTruncatedSphere.cpp
6 //! @brief Implements class FormFactorTruncatedSphere.
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/IntegratorGK.h"
20 #include <limits>
21 
23  : IBornFF({"TruncatedSphere",
24  "class_tooltip",
25  {{"Radius", "nm", "radius", 0, +INF, 0},
26  {"Height", "nm", "height before removal of cap", 0, +INF, 0},
27  {"DeltaHeight", "nm", "height of removed cap", 0, +INF, 0}}},
28  P)
29  , m_radius(m_P[0])
30  , m_height(m_P[1])
31  , m_dh(m_P[2])
32 {
33  check_initialization();
34  onChange();
35 }
36 
37 FormFactorTruncatedSphere::FormFactorTruncatedSphere(double radius, double height, double dh)
38  : FormFactorTruncatedSphere(std::vector<double>{radius, height, dh})
39 {
40 }
41 
43 {
44  bool result(true);
45  if (m_height > 2. * m_radius || m_dh > m_height) {
46  std::ostringstream ostr;
47  ostr << "::FormFactorTruncatedSphere() -> Error in class initialization ";
48  ostr << "with parameters 'radius':" << m_radius << " 'height':" << m_height
49  << " 'delta_height':" << m_dh << "\n\n";
50  ostr << "Check for height <= 2.*radius AND delta_height < height failed.";
51  throw std::runtime_error(ostr.str());
52  }
53  return result;
54 }
55 
56 //! Integrand for complex form factor.
58 {
59  double Rz = std::sqrt(m_radius * m_radius - Z * Z);
60  complex_t qx = m_q.x();
61  complex_t qy = m_q.y();
62  complex_t q_p = std::sqrt(qx * qx + qy * qy); // NOT the modulus!
63  return Rz * Rz * Math::Bessel::J1c(q_p * Rz) * exp_I(m_q.z() * Z);
64 }
65 
66 //! Complex form factor.
68 {
69  m_q = q;
70  if (std::abs(q.mag()) < std::numeric_limits<double>::epsilon()) {
71  return M_PI / 3.
72  * (m_height * m_height * (3. * m_radius - m_height)
73  - m_dh * m_dh * (3. * m_radius - m_dh));
74  }
75  // else
76  complex_t integral = ComplexIntegrator().integrate([&](double Z) { return Integrand(Z); },
78  return M_TWOPI * integral * exp_I(q.z() * (m_height - m_radius));
79 }
80 
82  kvector_t translation) const
83 {
84  double height = m_height - m_dh;
85  auto effects = computeSlicingEffects(limits, translation, height);
86  FormFactorTruncatedSphere slicedff(m_radius, m_height - effects.dz_bottom,
87  effects.dz_top + m_dh);
88  return createTransformedFormFactor(slicedff, rot, effects.position);
89 }
90 
92 {
93  m_shape3D = std::make_unique<TruncatedEllipsoid>(m_radius, m_radius, m_radius, m_height, m_dh);
94 }
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
#define M_PI
Definition: Constants.h:44
Defines class FormFactorTruncatedSphere.
const double INF
Definition: INode.h:25
Defines classes RealIntegrator, ComplexIntegrator.
Defines class TruncatedEllipsoid.
double mag() const
Returns magnitude of the vector.
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
To integrate a complex function of a real variable.
Definition: IntegratorGK.h:44
complex_t integrate(const std::function< complex_t(double)> &f, double lmin, double lmax)
complex_t evaluate_for_q(cvector_t q) const final
Complex form factor.
complex_t Integrand(double Z) const
Integrand for complex form factor.
void onChange() final
Action to be taken in inherited class when a parameter has changed.
FormFactorTruncatedSphere(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.
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
Definition: filesystem.h:81