BornAgain  1.19.0
Simulate and fit neutron and x-ray scattering at grazing incidence
FormFactorCone.cpp
Go to the documentation of this file.
1 // ************************************************************************************************
2 //
3 // BornAgain: simulate and fit reflection and scattering
4 //
5 //! @file Sample/HardParticle/FormFactorCone.cpp
6 //! @brief Implements class FormFactorCone.
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"
19 #include "Base/Math/IntegratorGK.h"
21 #include <limits>
22 
23 FormFactorCone::FormFactorCone(const std::vector<double> P)
24  : IBornFF({"Cone",
25  "frustum with circular base",
26  {{"Radius", "nm", "radius of base", 0, +INF, 0},
27  {"Height", "nm", "height", 0, +INF, 0},
28  {"Alpha", "rad", "angle between base and side", 0., M_PI_2, 0}}},
29  P)
30  , m_radius(m_P[0])
31  , m_height(m_P[1])
32  , m_alpha(m_P[2])
33 {
34  m_cot_alpha = Math::cot(m_alpha);
35  if (!std::isfinite(m_cot_alpha) || m_cot_alpha < 0)
36  throw std::runtime_error("pyramid angle alpha out of bounds");
37  if (m_cot_alpha * m_height > m_radius) {
38  std::ostringstream ostr;
39  ostr << "FormFactorCone() -> Error in class initialization ";
40  ostr << "with parameters radius:" << m_radius;
41  ostr << " m_height:" << m_height;
42  ostr << " alpha[rad]:" << m_alpha << "\n\n";
43  ostr << "Check for 'height <= radius*tan(alpha)' failed.";
44  throw std::runtime_error(ostr.str());
45  }
46  onChange();
47 }
48 
49 FormFactorCone::FormFactorCone(double radius, double height, double alpha)
50  : FormFactorCone(std::vector<double>{radius, height, alpha})
51 {
52 }
53 
54 //! Integrand for complex form factor.
56 {
57  double Rz = m_radius - Z * m_cot_alpha;
58  complex_t q_p = std::sqrt(m_q.x() * m_q.x() + m_q.y() * m_q.y()); // sqrt(x*x + y*y)
59  return Rz * Rz * Math::Bessel::J1c(q_p * Rz) * exp_I(m_q.z() * Z);
60 }
61 
63 {
64  m_q = q;
65  if (std::abs(m_q.mag()) < std::numeric_limits<double>::epsilon()) {
66  double R = m_radius;
67  double H = m_height;
68  if (m_cot_alpha == 0.0)
69  return M_PI * R * R * H; // cylinder case
70  double R2 = R - H * m_cot_alpha;
71  double apex_height = R / m_cot_alpha;
72  return M_PI / 3. * (R * R * H + (R * R - R2 * R2) * (apex_height - H));
73  } else {
74  complex_t integral =
75  ComplexIntegrator().integrate([&](double Z) { return Integrand(Z); }, 0., m_height);
76  return M_TWOPI * integral;
77  }
78 }
79 
81  kvector_t translation) const
82 {
83  auto effects = computeSlicingEffects(limits, translation, m_height);
84  double dradius = effects.dz_bottom * m_cot_alpha;
85  FormFactorCone slicedff(m_radius - dradius, m_height - effects.dz_bottom - effects.dz_top,
86  m_alpha);
87  return createTransformedFormFactor(slicedff, rot, effects.position);
88 }
89 
91 {
93  double radius2 = m_radius - m_height * m_cot_alpha;
94  m_shape3D = std::make_unique<DoubleEllipse>(m_radius, m_radius, m_height, radius2, radius2);
95 }
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_2
Definition: Constants.h:45
#define M_PI
Definition: Constants.h:44
Defines class DoubleEllipse.
Defines class FormFactorCone.
Defines functions in namespace Math.
const double INF
Definition: INode.h:25
Defines classes RealIntegrator, ComplexIntegrator.
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)
A conical frustum (cone truncated parallel to the base) with circular base.
const double & m_alpha
FormFactorCone(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.
const double & m_radius
IFormFactor * sliceFormFactor(ZLimits limits, const IRotation &rot, kvector_t translation) const final
Actually slices the form factor or throws an exception.
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.
const double & m_height
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 cot(double x)
cotangent function:
Definition: Functions.cpp:48
Definition: filesystem.h:81