BornAgain  1.19.79
Open-source research software to simulate and fit neutron and x-ray reflectometry and grazing-incidence small-angle scattering
IFormFactor.h
Go to the documentation of this file.
1 // ************************************************************************************************
2 //
3 // BornAgain: simulate and fit reflection and scattering
4 //
5 //! @file Sample/Particle/IFormFactor.h
6 //! @brief Defines interface IDecoratableBorn.
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_PARTICLE_IFORMFACTOR_H
17 #define BORNAGAIN_SAMPLE_PARTICLE_IFORMFACTOR_H
18 
20 #include <heinz/Complex.h>
21 #include <heinz/Vectors3D.h>
22 
23 #include "Base/Spin/SpinMatrix.h"
24 #include "Base/Spin/Spinor.h"
25 
26 class IRotation;
27 class IShape3D;
28 class WavevectorInfo;
29 
30 //! Abstract base class for Born form factors.
31 //!
32 //! In contrast to the generic IReParticle, a Born form factor does not depend
33 //! on the incoming and outgoing wave vectors ki and kf, except through their
34 //! difference, the scattering vector q=ki-kf.
35 
36 class IFormFactor : public ISampleNode {
37 public:
39  IFormFactor(const std::vector<double>& PValues);
40  ~IFormFactor() override;
41 
42  IFormFactor* clone() const override = 0;
43 
44  virtual complex_t theFF(const WavevectorInfo& wavevectors) const;
45  virtual SpinMatrix thePolFF(const WavevectorInfo& wavevectors) const;
46 
47  std::string shapeName() const;
48 
49  virtual double bottomZ(const IRotation* rotation) const;
50  virtual double topZ(const IRotation* rotation) const;
51 
52  virtual double volume() const;
53 
54  //! Returns the (approximate in some cases) radial size of the particle of this
55  //! form factor's shape. This is used for SSCA calculations
56  virtual double radialExtension() const = 0;
57 
58  virtual complex_t formfactor_at_bottom(C3 q) const = 0;
59 
60  //! Creates the Python constructor of this class (or derived classes)
61  virtual std::string pythonConstructor() const;
62 
63  //! Default implementation only allows rotations along z-axis
64  virtual bool canSliceAnalytically(const IRotation* rot) const;
65 
66  //! Returns scattering amplitude for complex scattering wavevector q=k_i-k_f in case
67  //! of matrix interactions. Default implementation calls formfactor_at_bottom(q) and
68  //! multiplies with the unit matrix.
69  virtual SpinMatrix formfactor_pol(C3 q) const;
70 
71 protected:
72  //! IShape3D object, used to retrieve vertices (which may be approximate in the case
73  //! of round shapes). For soft particles, this will be a hard mean shape.
74  std::unique_ptr<IShape3D> m_shape3D;
75 };
76 
77 #endif // BORNAGAIN_SAMPLE_PARTICLE_IFORMFACTOR_H
78 #endif // USER_API
Defines interface class ISampleNode.
Defines class SpinMatrix.
Defines class Spinor.
Abstract base class for Born form factors.
Definition: IFormFactor.h:36
virtual SpinMatrix thePolFF(const WavevectorInfo &wavevectors) const
Definition: IFormFactor.cpp:45
virtual std::string pythonConstructor() const
Creates the Python constructor of this class (or derived classes)
Definition: IFormFactor.cpp:69
virtual double radialExtension() const =0
Returns the (approximate in some cases) radial size of the particle of this form factor's shape....
virtual double bottomZ(const IRotation *rotation) const
Definition: IFormFactor.cpp:50
virtual complex_t formfactor_at_bottom(C3 q) const =0
virtual SpinMatrix formfactor_pol(C3 q) const
Returns scattering amplitude for complex scattering wavevector q=k_i-k_f in case of matrix interactio...
Definition: IFormFactor.cpp:78
~IFormFactor() override
IFormFactor * clone() const override=0
Returns a clone of this ISampleNode object.
virtual double topZ(const IRotation *rotation) const
Definition: IFormFactor.cpp:57
virtual bool canSliceAnalytically(const IRotation *rot) const
Default implementation only allows rotations along z-axis.
Definition: IFormFactor.cpp:64
virtual complex_t theFF(const WavevectorInfo &wavevectors) const
Definition: IFormFactor.cpp:40
std::unique_ptr< IShape3D > m_shape3D
IShape3D object, used to retrieve vertices (which may be approximate in the case of round shapes)....
Definition: IFormFactor.h:74
std::string shapeName() const
Definition: IFormFactor.cpp:33
virtual double volume() const
Definition: IFormFactor.cpp:83
Abstract base class for rotations.
Definition: Rotations.h:29
Abstract base class for sample components and properties related to scattering.
Definition: ISampleNode.h:27
Abstract base class for different shapes.
Definition: IShape3D.h:31
Holds all wavevector information relevant for calculating form factors.