BornAgain  1.19.0
Simulate and fit neutron and x-ray scattering at grazing incidence
FormFactorDecoratorPositionFactor.h
Go to the documentation of this file.
1 // ************************************************************************************************
2 //
3 // BornAgain: simulate and fit reflection and scattering
4 //
5 //! @file Sample/Scattering/FormFactorDecoratorPositionFactor.h
6 //! @brief Defines 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 
15 #ifdef SWIG
16 #error no need to expose this header to Swig
17 #endif
18 
19 #ifndef USER_API
20 #ifndef BORNAGAIN_SAMPLE_SCATTERING_FORMFACTORDECORATORPOSITIONFACTOR_H
21 #define BORNAGAIN_SAMPLE_SCATTERING_FORMFACTORDECORATORPOSITIONFACTOR_H
22 
24 
25 //! Decorates a form factor with a position dependent phase factor.
26 //! @ingroup formfactors_internal
27 
29 public:
30  FormFactorDecoratorPositionFactor(const IFormFactor& ff, const kvector_t& position);
31 
33  {
35  }
36 
37  void accept(INodeVisitor* visitor) const final { visitor->visit(this); }
38 
39  double bottomZ(const IRotation& rotation) const final;
40 
41  double topZ(const IRotation& rotation) const final;
42 
43  complex_t evaluate(const WavevectorInfo& wavevectors) const final;
44 #ifndef SWIG
45  Eigen::Matrix2cd evaluatePol(const WavevectorInfo& wavevectors) const final;
46 #endif
47 
48 private:
49  complex_t getPositionFactor(const WavevectorInfo& wavevectors) const;
50 
52 };
53 
54 #endif // BORNAGAIN_SAMPLE_SCATTERING_FORMFACTORDECORATORPOSITIONFACTOR_H
55 #endif // USER_API
std::complex< double > complex_t
Definition: Complex.h:20
Defines and implements interface class IFormFactorDecorator.
Decorates a form factor with a position dependent phase factor.
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)
FormFactorDecoratorPositionFactor * clone() const final
Returns a clone of this ISampleNode object.
complex_t getPositionFactor(const WavevectorInfo &wavevectors) const
void accept(INodeVisitor *visitor) const final
Calls the INodeVisitor's visit method.
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
Visitor interface to visit ISampleNode objects.
Definition: INodeVisitor.h:146
Abstract base class for rotations.
Definition: Rotations.h:28
Holds all wavevector information relevant for calculating form factors.