BornAgain  1.19.0
Simulate and fit neutron and x-ray scattering at grazing incidence
IBornFF.h
Go to the documentation of this file.
1 // ************************************************************************************************
2 //
3 // BornAgain: simulate and fit reflection and scattering
4 //
5 //! @file Sample/Scattering/IBornFF.h
6 //! @brief Defines interface IBornFF.
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 #ifndef USER_API
16 #ifndef BORNAGAIN_SAMPLE_SCATTERING_IBORNFF_H
17 #define BORNAGAIN_SAMPLE_SCATTERING_IBORNFF_H
18 
20 
21 class IShape3D;
22 
23 //! Nested structure that holds slicing effects on position and removed parts.
24 
25 //! @ingroup formfactors_internal
26 
29  double dz_bottom;
30  double dz_top;
31 };
32 
33 //! Abstract base class for Born form factors.
34 //!
35 //! In contrast to the generic IFormFactor, a Born form factor does not depend
36 //! on the incoming and outgoing wave vectors ki and kf, except through their
37 //! difference, the scattering vector q=ki-kf.
38 
39 //! @ingroup formfactors_internal
40 
41 class IBornFF : public IFormFactor {
42 public:
44  IBornFF(const NodeMeta& meta, const std::vector<double>& PValues);
46 
47  IBornFF* clone() const override = 0;
48 
49  void setAmbientMaterial(const Material&) override {}
50 
51  complex_t evaluate(const WavevectorInfo& wavevectors) const override;
52 
53 #ifndef SWIG
54  Eigen::Matrix2cd evaluatePol(const WavevectorInfo& wavevectors) const override;
55 #endif
56 
57  virtual double bottomZ(const IRotation& rotation) const override;
58  virtual double topZ(const IRotation& rotation) const override;
59 
60  //! Returns scattering amplitude for complex scattering wavevector q=k_i-k_f.
61  //! This method is public only for convenience of plotting form factors in Python.
62  virtual complex_t evaluate_for_q(cvector_t q) const = 0;
63 
64 protected:
65  //! Default implementation only allows rotations along z-axis
66  bool canSliceAnalytically(const IRotation& rot) const override;
67 
68 #ifndef SWIG
69  //! Returns scattering amplitude for complex scattering wavevector q=k_i-k_f in case
70  //! of matrix interactions. Default implementation calls evaluate_for_q(q) and
71  //! multiplies with the unit matrix.
72  virtual Eigen::Matrix2cd evaluate_for_q_pol(cvector_t q) const;
73 #endif
74 
75  //! IShape3D object, used to retrieve vertices (which may be approximate in the case
76  //! of round shapes). For soft particles, this will be a hard mean shape.
77  std::unique_ptr<IShape3D> m_shape3D;
78 
79  //! Helper method for slicing
80  static SlicingEffects computeSlicingEffects(ZLimits limits, const kvector_t& position,
81  double height);
82 
83  //! Calculates the z-coordinate of the lowest vertex after rotation
84  static double BottomZ(const std::vector<kvector_t>& vertices, const IRotation& rotation);
85 
86  //! Calculates the z-coordinate of the highest vertex after rotation
87  static double TopZ(const std::vector<kvector_t>& vertices, const IRotation& rotation);
88 };
89 
90 #endif // BORNAGAIN_SAMPLE_SCATTERING_IBORNFF_H
91 #endif // USER_API
std::complex< double > complex_t
Definition: Complex.h:20
Defines and implements interface IFormFactor.
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
void setAmbientMaterial(const Material &) override
Passes the material in which this particle is embedded.
Definition: IBornFF.h:49
virtual Eigen::Matrix2cd evaluate_for_q_pol(cvector_t q) const
Returns scattering amplitude for complex scattering wavevector q=k_i-k_f in case of matrix interactio...
Definition: IBornFF.cpp:61
Eigen::Matrix2cd evaluatePol(const WavevectorInfo &wavevectors) const override
Returns scattering amplitude for matrix interactions.
Definition: IBornFF.cpp:35
static double TopZ(const std::vector< kvector_t > &vertices, const IRotation &rotation)
Calculates the z-coordinate of the highest vertex after rotation.
Definition: IBornFF.cpp:98
bool canSliceAnalytically(const IRotation &rot) const override
Default implementation only allows rotations along z-axis.
Definition: IBornFF.cpp:54
virtual double bottomZ(const IRotation &rotation) const override
Returns the z-coordinate of the lowest point in this shape after a given rotation.
Definition: IBornFF.cpp:40
IBornFF * clone() const override=0
Returns a clone of this ISampleNode object.
virtual double topZ(const IRotation &rotation) const override
Returns the z-coordinate of the lowest point in this shape after a given rotation.
Definition: IBornFF.cpp:47
virtual complex_t evaluate_for_q(cvector_t q) const =0
Returns scattering amplitude for complex scattering wavevector q=k_i-k_f.
complex_t evaluate(const WavevectorInfo &wavevectors) const override
Returns scattering amplitude for complex wavevectors ki, kf.
Definition: IBornFF.cpp:30
static double BottomZ(const std::vector< kvector_t > &vertices, const IRotation &rotation)
Calculates the z-coordinate of the lowest vertex after rotation.
Definition: IBornFF.cpp:90
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
Abstract base class for rotations.
Definition: Rotations.h:28
Abstract base class for different shapes.
Definition: IShape3D.h:33
A wrapper for underlying material implementation.
Definition: Material.h:29
Holds all wavevector information relevant for calculating form factors.
Class that contains upper and lower limits of the z-coordinate for the slicing of form factors.
Definition: ZLimits.h:45
Metadata of one model node.
Definition: INode.h:38
Nested structure that holds slicing effects on position and removed parts.
Definition: IBornFF.h:27
double dz_top
Definition: IBornFF.h:30
kvector_t position
Definition: IBornFF.h:28
double dz_bottom
Definition: IBornFF.h:29