BornAgain  1.19.0
Simulate and fit neutron and x-ray scattering at grazing incidence
IFormFactor.h
Go to the documentation of this file.
1 // ************************************************************************************************
2 //
3 // BornAgain: simulate and fit reflection and scattering
4 //
5 //! @file Sample/Scattering/IFormFactor.h
6 //! @brief Defines and implements interface IFormFactor.
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_IFORMFACTOR_H
17 #define BORNAGAIN_SAMPLE_SCATTERING_IFORMFACTOR_H
18 
21 #include <Eigen/Core>
22 
24 class IRotation;
25 class Material;
26 class WavevectorInfo;
27 
28 //! Abstract base class for all form factors.
29 //!
30 //! The actual form factor is returned by the complex valued function IFormFactor::evaluate,
31 //! which depends on the incoming and outgoing wave vectors ki and kf.
32 //! If it only depends on the scattering vector q=ki-kf, then it is a IBornFF.
33 
34 //! @ingroup formfactors_internal
35 
36 class IFormFactor : public ISampleNode {
37 public:
38  IFormFactor() = default;
39  IFormFactor(const NodeMeta& meta, const std::vector<double>& PValues);
40 
41  ~IFormFactor() = default;
42  IFormFactor* clone() const override = 0;
43 
44  //! Creates a (possibly sliced) form factor with the given rotation and translation
46  kvector_t translation) const;
47 
48  //! Passes the material in which this particle is embedded.
49  virtual void setAmbientMaterial(const Material&) = 0;
50 
51  //! Returns scattering amplitude for complex wavevectors ki, kf.
52  virtual complex_t evaluate(const WavevectorInfo& wavevectors) const = 0;
53 
54 #ifndef SWIG
55  //! Returns scattering amplitude for matrix interactions
56  virtual Eigen::Matrix2cd evaluatePol(const WavevectorInfo& wavevectors) const;
57 #endif
58 
59  //! Returns the total volume of the particle of this form factor's shape
60  virtual double volume() const;
61 
62  //! Returns the (approximate in some cases) radial size of the particle of this
63  //! form factor's shape. This is used for SSCA calculations
64  virtual double radialExtension() const = 0;
65 
66  //! Returns the z-coordinate of the lowest point in this shape after a given rotation
67  virtual double bottomZ(const IRotation& rotation) const = 0;
68 
69  //! Returns the z-coordinate of the lowest point in this shape after a given rotation
70  virtual double topZ(const IRotation& rotation) const = 0;
71 
72 #ifndef SWIG
73  //! Sets reflection/transmission info
74  virtual void setSpecularInfo(std::unique_ptr<const ILayerRTCoefficients>,
75  std::unique_ptr<const ILayerRTCoefficients>);
76 #endif
77 
78 protected:
79  //! Checks if slicing has a fast analytical solution
80  virtual bool canSliceAnalytically(const IRotation& rot) const;
81 
82  //! Actually slices the form factor or throws an exception
83  virtual IFormFactor* sliceFormFactor(ZLimits limits, const IRotation& rot,
84  kvector_t translation) const;
85 
86  static IFormFactor* createTransformedFormFactor(const IFormFactor& formfactor,
87  const IRotation& rot, kvector_t translation);
88 };
89 
90 #endif // BORNAGAIN_SAMPLE_SCATTERING_IFORMFACTOR_H
91 #endif // USER_API
std::complex< double > complex_t
Definition: Complex.h:20
Defines interface class ISampleNode.
Defines class ZLimits.
Abstract base class for all form factors.
Definition: IFormFactor.h:36
~IFormFactor()=default
virtual double radialExtension() const =0
Returns the (approximate in some cases) radial size of the particle of this form factor's shape.
virtual IFormFactor * sliceFormFactor(ZLimits limits, const IRotation &rot, kvector_t translation) const
Actually slices the form factor or throws an exception.
Definition: IFormFactor.cpp:94
virtual complex_t evaluate(const WavevectorInfo &wavevectors) const =0
Returns scattering amplitude for complex wavevectors ki, kf.
virtual bool canSliceAnalytically(const IRotation &rot) const
Checks if slicing has a fast analytical solution.
Definition: IFormFactor.cpp:89
IFormFactor()=default
IFormFactor * createSlicedFormFactor(ZLimits limits, const IRotation &rot, kvector_t translation) const
Creates a (possibly sliced) form factor with the given rotation and translation.
Definition: IFormFactor.cpp:58
virtual double topZ(const IRotation &rotation) const =0
Returns the z-coordinate of the lowest point in this shape after a given rotation.
IFormFactor * clone() const override=0
Returns a clone of this ISampleNode object.
static IFormFactor * createTransformedFormFactor(const IFormFactor &formfactor, const IRotation &rot, kvector_t translation)
Definition: IFormFactor.cpp:99
virtual double bottomZ(const IRotation &rotation) const =0
Returns the z-coordinate of the lowest point in this shape after a given rotation.
virtual void setAmbientMaterial(const Material &)=0
Passes the material in which this particle is embedded.
virtual Eigen::Matrix2cd evaluatePol(const WavevectorInfo &wavevectors) const
Returns scattering amplitude for matrix interactions.
Definition: IFormFactor.cpp:72
virtual void setSpecularInfo(std::unique_ptr< const ILayerRTCoefficients >, std::unique_ptr< const ILayerRTCoefficients >)
Sets reflection/transmission info.
Definition: IFormFactor.cpp:84
virtual double volume() const
Returns the total volume of the particle of this form factor's shape.
Definition: IFormFactor.cpp:78
Interface to access reflection/transmission coefficients.
Abstract base class for rotations.
Definition: Rotations.h:28
Abstract base class for sample components and properties related to scattering.
Definition: ISampleNode.h:28
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