BornAgain  1.18.0
Simulate and fit neutron and x-ray scattering at grazing incidence
FormFactorBox.cpp
Go to the documentation of this file.
1 // ************************************************************************** //
2 //
3 // BornAgain: simulate and fit scattering at grazing incidence
4 //
5 //! @file Sample/HardParticle/FormFactorBox.cpp
6 //! @brief Implements class FormFactorBox.
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 
18 FormFactorBox::FormFactorBox(const std::vector<double> P)
19  : IFormFactorPrism({"Box",
20  "rectangular cuboid",
21  {{"Length", "nm", "side length in x direction", 0, +INF, 0},
22  {"Width", "nm", "side length in y direction", 0, +INF, 0},
23  {"Height", "nm", "side length in z direction", 0, +INF, 0}}},
24  P),
25  m_length(m_P[0]), m_width(m_P[1]), m_height(m_P[2])
26 {
27  onChange();
28 }
29 
31  : FormFactorBox(std::vector<double>{length, width, height})
32 {
33 }
34 
36 {
37  complex_t qzHdiv2 = m_height / 2 * q.z();
38  return m_length * m_width * m_height * MathFunctions::sinc(m_length / 2 * q.x())
39  * MathFunctions::sinc(m_width / 2 * q.y()) * MathFunctions::sinc(qzHdiv2)
40  * exp_I(qzHdiv2);
41 }
42 
44  kvector_t translation) const
45 {
46  auto effects = computeSlicingEffects(limits, translation, m_height);
47  FormFactorBox slicedff(m_length, m_width, m_height - effects.dz_bottom - effects.dz_top);
48  return createTransformedFormFactor(slicedff, rot, effects.position);
49 }
50 
52 {
53  double a = m_length / 2;
54  double b = m_width / 2;
55  std::vector<kvector_t> V{{a, b, 0.}, {-a, b, 0.}, {-a, -b, 0.}, {a, -b, 0}};
56  setPrism(true, V);
57 }
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 FormFactorBox.
IFormFactor * createTransformedFormFactor(const IFormFactor &formfactor, const IRotation &rot, kvector_t translation)
Definition: IFormFactor.cpp:77
const double INF
Definition: INode.h:24
Defines namespace MathFunctions.
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 rectangular prism (parallelepiped).
Definition: FormFactorBox.h:24
const double & m_width
Definition: FormFactorBox.h:52
void onChange() override final
Action to be taken in inherited class when a parameter has changed.
complex_t evaluate_for_q(cvector_t q) const override final
Returns scattering amplitude for complex scattering wavevector q=k_i-k_f.
IFormFactor * sliceFormFactor(ZLimits limits, const IRotation &rot, kvector_t translation) const override final
Actually slices the form factor or throws an exception.
const double & m_length
Definition: FormFactorBox.h:51
FormFactorBox(const std::vector< double > P)
const double & m_height
Definition: FormFactorBox.h:53
double height() const final
Definition: FormFactorBox.h:48
SlicingEffects computeSlicingEffects(ZLimits limits, const kvector_t &position, double height) const
Helper method for slicing.
A prism with a polygonal base, for form factor computation.
void setPrism(bool symmetry_Ci, const std::vector< kvector_t > &vertices)
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 sinc(double x)
sinc function: