BornAgain  1.19.0
Simulate and fit neutron and x-ray scattering at grazing incidence
IComputeFF.h
Go to the documentation of this file.
1 // ************************************************************************************************
2 //
3 // BornAgain: simulate and fit reflection and scattering
4 //
5 //! @file Sample/FFCompute/IComputeFF.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 #ifdef SWIG
16 #error no need to expose this header to Swig
17 #endif
18 
19 #ifndef USER_API
20 #ifndef BORNAGAIN_SAMPLE_FFCOMPUTE_ICOMPUTEFF_H
21 #define BORNAGAIN_SAMPLE_FFCOMPUTE_ICOMPUTEFF_H
22 
23 #include "Base/Types/Complex.h"
24 #include <Eigen/Core>
25 #include <memory>
26 
27 class IFormFactor;
29 class IRotation;
30 class Material;
31 class WavevectorInfo;
32 
33 //! Abstract base class for form factor evaluations.
34 //!
35 //! Wraps an IFormFactor, and provides functions evaluate or evaluatePol.
36 
37 //! @ingroup formfactors_internal
38 
39 class IComputeFF {
40 
41 public:
42  IComputeFF() = delete;
43  virtual ~IComputeFF();
44  virtual IComputeFF* clone() const = 0;
45 
46  virtual void setAmbientMaterial(const Material& material);
47 
48  virtual double volume() const;
49  virtual double radialExtension() const;
50  virtual double bottomZ(const IRotation& rotation) const;
51  virtual double topZ(const IRotation& rotation) const;
52  virtual complex_t evaluate(const WavevectorInfo& wavevectors) const = 0;
53 #ifndef SWIG
54  //! Returns scattering amplitude for matrix interactions
55  virtual Eigen::Matrix2cd evaluatePol(const WavevectorInfo& wavevectors) const;
56  //! Sets reflection/transmission info
57  virtual void setSpecularInfo(std::unique_ptr<const ILayerRTCoefficients>,
58  std::unique_ptr<const ILayerRTCoefficients>);
59 #endif
60 
61 protected:
62  IComputeFF(const IFormFactor&);
63 
64  std::unique_ptr<IFormFactor> m_ff;
65 };
66 
67 #endif // BORNAGAIN_SAMPLE_FFCOMPUTE_ICOMPUTEFF_H
68 #endif // USER_API
Defines complex_t, and a few elementary functions.
std::complex< double > complex_t
Definition: Complex.h:20
Abstract base class for form factor evaluations.
Definition: IComputeFF.h:39
virtual double topZ(const IRotation &rotation) const
Definition: IComputeFF.cpp:44
virtual complex_t evaluate(const WavevectorInfo &wavevectors) const =0
virtual void setAmbientMaterial(const Material &material)
Definition: IComputeFF.cpp:24
IComputeFF()=delete
virtual IComputeFF * clone() const =0
virtual void setSpecularInfo(std::unique_ptr< const ILayerRTCoefficients >, std::unique_ptr< const ILayerRTCoefficients >)
Sets reflection/transmission info.
Definition: IComputeFF.cpp:54
std::unique_ptr< IFormFactor > m_ff
Definition: IComputeFF.h:64
virtual double bottomZ(const IRotation &rotation) const
Definition: IComputeFF.cpp:39
virtual double radialExtension() const
Definition: IComputeFF.cpp:34
virtual ~IComputeFF()
virtual Eigen::Matrix2cd evaluatePol(const WavevectorInfo &wavevectors) const
Returns scattering amplitude for matrix interactions.
Definition: IComputeFF.cpp:49
virtual double volume() const
Definition: IComputeFF.cpp:29
Abstract base class for all form factors.
Definition: IFormFactor.h:36
Interface to access reflection/transmission coefficients.
Abstract base class for rotations.
Definition: Rotations.h:28
A wrapper for underlying material implementation.
Definition: Material.h:29
Holds all wavevector information relevant for calculating form factors.