BornAgain  1.18.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 scattering at grazing incidence
4 //
5 //! @file Sample/Scattering/IFormFactor.h
6 //! @brief Defines and implements pure virtual 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 BORNAGAIN_CORE_SCATTERING_IFORMFACTOR_H
16 #define BORNAGAIN_CORE_SCATTERING_IFORMFACTOR_H
17 
18 #include "Base/Types/Complex.h"
21 #include <Eigen/Core>
22 
24 class IRotation;
25 class Material;
26 class WavevectorInfo;
27 
28 //! Pure virtual 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 IBornFormFactor.
33 //!
34 //! Other children besides IBornFormFactor are IFormFactorDecorator, FormFactorWeighted,
35 //! FormFactorDWBA, FormFactorDWBAPol and FormFactorCrystal.
36 
37 //! @ingroup formfactors_internal
38 
39 class IFormFactor : public ISample
40 {
41 public:
42  IFormFactor() = default;
43  IFormFactor(const NodeMeta& meta, const std::vector<double>& PValues);
44 
45  ~IFormFactor() = default;
46  IFormFactor* clone() const override = 0;
47 
48  //! Creates a (possibly sliced) form factor with the given rotation and translation
50  kvector_t translation) const;
51 
52  //! Passes the material in which this particle is embedded.
53  virtual void setAmbientMaterial(const Material&) = 0;
54 
55  //! Returns scattering amplitude for complex wavevectors ki, kf.
56  virtual complex_t evaluate(const WavevectorInfo& wavevectors) const = 0;
57 
58 #ifndef SWIG
59  //! Returns scattering amplitude for matrix interactions
60  virtual Eigen::Matrix2cd evaluatePol(const WavevectorInfo& wavevectors) const;
61 #endif
62 
63  //! Returns the total volume of the particle of this form factor's shape
64  virtual double volume() const;
65 
66  //! Returns the (approximate in some cases) radial size of the particle of this
67  //! form factor's shape. This is used for SSCA calculations
68  virtual double radialExtension() const = 0;
69 
70  //! Returns the z-coordinate of the lowest point in this shape after a given rotation
71  virtual double bottomZ(const IRotation& rotation) const = 0;
72 
73  //! Returns the z-coordinate of the lowest point in this shape after a given rotation
74  virtual double topZ(const IRotation& rotation) const = 0;
75 
76 #ifndef SWIG
77  //! Sets reflection/transmission info
78  virtual void setSpecularInfo(std::unique_ptr<const ILayerRTCoefficients>,
79  std::unique_ptr<const ILayerRTCoefficients>);
80 #endif
81 
82 protected:
83  //! Checks if slicing has a fast analytical solution
84  virtual bool canSliceAnalytically(const IRotation& rot) const;
85 
86  //! Actually slices the form factor or throws an exception
87  virtual IFormFactor* sliceFormFactor(ZLimits limits, const IRotation& rot,
88  kvector_t translation) const;
89 };
90 
91 IFormFactor* createTransformedFormFactor(const IFormFactor& formfactor, const IRotation& rot,
92  kvector_t translation);
93 
94 #endif // BORNAGAIN_CORE_SCATTERING_IFORMFACTOR_H
Defines complex_t, and a few elementary functions.
std::complex< double > complex_t
Definition: Complex.h:20
IFormFactor * createTransformedFormFactor(const IFormFactor &formfactor, const IRotation &rot, kvector_t translation)
Definition: IFormFactor.cpp:77
Defines interface class ISample.
Defines class ZLimits.
Pure virtual base class for all form factors.
Definition: IFormFactor.h:40
~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:72
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:67
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:35
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 ISample object.
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:49
virtual void setSpecularInfo(std::unique_ptr< const ILayerRTCoefficients >, std::unique_ptr< const ILayerRTCoefficients >)
Sets reflection/transmission info.
Definition: IFormFactor.cpp:62
virtual double volume() const
Returns the total volume of the particle of this form factor's shape.
Definition: IFormFactor.cpp:56
Interface to access reflection/transmission coefficients.
Pure virtual interface for rotations.
Definition: Rotations.h:27
Pure virtual base class for sample components and properties related to scattering.
Definition: ISample.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:41
Metadata of one model node.
Definition: INode.h:37