BornAgain  1.18.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 scattering at grazing incidence
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 
17 #include "Base/Types/Exceptions.h"
18 #include "Base/Utils/Integrator.h"
21 #include <limits>
22 
23 FormFactorCone::FormFactorCone(const std::vector<double> P)
24  : IFormFactorBorn({"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]), m_height(m_P[1]), m_alpha(m_P[2])
31 {
32  m_cot_alpha = MathFunctions::cot(m_alpha);
33  if (!std::isfinite(m_cot_alpha) || m_cot_alpha < 0)
34  throw Exceptions::OutOfBoundsException("pyramid angle alpha out of bounds");
35  if (m_cot_alpha * m_height > m_radius) {
36  std::ostringstream ostr;
37  ostr << "FormFactorCone() -> Error in class initialization ";
38  ostr << "with parameters radius:" << m_radius;
39  ostr << " m_height:" << m_height;
40  ostr << " alpha[rad]:" << m_alpha << "\n\n";
41  ostr << "Check for 'height <= radius*tan(alpha)' failed.";
43  }
44  onChange();
45 }
46 
47 FormFactorCone::FormFactorCone(double radius, double height, double alpha)
48  : FormFactorCone(std::vector<double>{radius, height, alpha})
49 {
50 }
51 
52 //! Integrand for complex form factor.
54 {
55  double Rz = m_radius - Z * m_cot_alpha;
56  complex_t q_p = std::sqrt(m_q.x() * m_q.x() + m_q.y() * m_q.y()); // sqrt(x*x + y*y)
57  return Rz * Rz * MathFunctions::Bessel_J1c(q_p * Rz) * exp_I(m_q.z() * Z);
58 }
59 
61 {
62  m_q = q;
63  if (std::abs(m_q.mag()) < std::numeric_limits<double>::epsilon()) {
64  double R = m_radius;
65  double H = m_height;
66  if (m_cot_alpha == 0.0)
67  return M_PI * R * R * H; // cylinder case
68  double R2 = R - H * m_cot_alpha;
69  double apex_height = R / m_cot_alpha;
70  return M_PI / 3. * (R * R * H + (R * R - R2 * R2) * (apex_height - H));
71  } else {
72  complex_t integral =
73  ComplexIntegrator().integrate([&](double Z) { return Integrand(Z); }, 0., m_height);
74  return M_TWOPI * integral;
75  }
76 }
77 
79  kvector_t translation) const
80 {
81  auto effects = computeSlicingEffects(limits, translation, m_height);
82  double dradius = effects.dz_bottom * m_cot_alpha;
83  FormFactorCone slicedff(m_radius - dradius, m_height - effects.dz_bottom - effects.dz_top,
84  m_alpha);
85  return createTransformedFormFactor(slicedff, rot, effects.position);
86 }
87 
89 {
91  double radius2 = m_radius - m_height * m_cot_alpha;
92  mP_shape = std::make_unique<DoubleEllipse>(m_radius, m_radius, m_height, radius2, radius2);
93 }
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 many exception classes in namespace Exceptionss.
Defines class FormFactorCone.
IFormFactor * createTransformedFormFactor(const IFormFactor &formfactor, const IRotation &rot, kvector_t translation)
Definition: IFormFactor.cpp:77
const double INF
Definition: INode.h:24
Defines classes RealIntegrator, ComplexIntegrator.
Defines M_PI and some more mathematical constants.
#define M_TWOPI
Definition: MathConstants.h:49
#define M_PI_2
Definition: MathConstants.h:40
#define M_PI
Definition: MathConstants.h:39
Defines namespace MathFunctions.
double mag() const
Returns magnitude of the vector.
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
To integrate a complex function of a real variable.
Definition: Integrator.h:41
complex_t integrate(const std::function< complex_t(double)> &f, double lmin, double lmax)
Definition: Integrator.cpp:36
A conical frustum (cone truncated parallel to the base) with circular base.
const double & m_alpha
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.
FormFactorCone(const std::vector< double > P)
const double & m_radius
complex_t Integrand(double Z) const
Integrand for complex form factor.
complex_t evaluate_for_q(cvector_t q) const override final
Returns scattering amplitude for complex scattering wavevector q=k_i-k_f.
const double & m_height
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 cot(double x)
cotangent function:
double Bessel_J1c(double x)
Bessel function Bessel_J1(x)/x.
const double radius(5 *Units::nanometer)