BornAgain  1.19.0
Simulate and fit neutron and x-ray scattering at grazing incidence
FormFactorDecoratorPositionFactor.cpp
Go to the documentation of this file.
1 // ************************************************************************************************
2 //
3 // BornAgain: simulate and fit reflection and scattering
4 //
5 //! @file Sample/Scattering/FormFactorDecoratorPositionFactor.cpp
6 //! @brief Implements class FormFactorDecoratorPositionFactor.
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 
18 
20  const kvector_t& position)
21  : IFormFactorDecorator(ff), m_position(position)
22 {
23  setName("FormFactorDecoratorPositionFactor");
24 }
25 
27 {
28  kvector_t rotated_translation = rotation.transformed(m_position);
29  return m_ff->bottomZ(rotation) + rotated_translation.z();
30 }
31 
33 {
34  kvector_t rotated_translation = rotation.transformed(m_position);
35  return m_ff->topZ(rotation) + rotated_translation.z();
36 }
37 
39 {
40  return getPositionFactor(wavevectors) * m_ff->evaluate(wavevectors);
41 }
42 
43 Eigen::Matrix2cd
45 {
46  return getPositionFactor(wavevectors) * m_ff->evaluatePol(wavevectors);
47 }
48 
51 {
52  cvector_t q = wavevectors.getQ();
53  return exp_I(m_position.dot(q));
54 }
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 FormFactorDecoratorPositionFactor.
Defines IRotation classes.
Defines WavevectorInfo.
auto dot(const BasicVector3D< U > &v) const
Returns dot product of vectors (antilinear in the first [=self] argument).
T z() const
Returns z-component in cartesian coordinate system.
Definition: BasicVector3D.h:67
double bottomZ(const IRotation &rotation) const final
Returns the z-coordinate of the lowest point in this shape after a given rotation.
complex_t evaluate(const WavevectorInfo &wavevectors) const final
Returns scattering amplitude for complex wavevectors ki, kf.
FormFactorDecoratorPositionFactor(const IFormFactor &ff, const kvector_t &position)
complex_t getPositionFactor(const WavevectorInfo &wavevectors) const
Eigen::Matrix2cd evaluatePol(const WavevectorInfo &wavevectors) const final
Returns scattering amplitude for matrix interactions.
double topZ(const IRotation &rotation) const final
Returns the z-coordinate of the lowest point in this shape after a given rotation.
Encapsulates another formfactor and adds extra functionality (a scalar factor, a position-dependent p...
Abstract base class for all form factors.
Definition: IFormFactor.h:36
virtual complex_t evaluate(const WavevectorInfo &wavevectors) const =0
Returns scattering amplitude for complex wavevectors ki, kf.
virtual double topZ(const IRotation &rotation) const =0
Returns the z-coordinate of the lowest point in this shape after a given rotation.
virtual double bottomZ(const IRotation &rotation) const =0
Returns the z-coordinate of the lowest point in this shape after a given rotation.
virtual Eigen::Matrix2cd evaluatePol(const WavevectorInfo &wavevectors) const
Returns scattering amplitude for matrix interactions.
Definition: IFormFactor.cpp:72
void setName(const std::string &name)
Abstract base class for rotations.
Definition: Rotations.h:28
kvector_t transformed(const kvector_t &v) const
Definition: Rotations.cpp:53
Holds all wavevector information relevant for calculating form factors.
cvector_t getQ() const