BornAgain  1.18.0
Simulate and fit neutron and x-ray scattering at grazing incidence
IFormFactorBorn.h
Go to the documentation of this file.
1 // ************************************************************************** //
2 //
3 // BornAgain: simulate and fit scattering at grazing incidence
4 //
5 //! @file Sample/Scattering/IFormFactorBorn.h
6 //! @brief Defines pure virtual interface class IFormFactorBorn.
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_IFORMFACTORBORN_H
16 #define BORNAGAIN_CORE_SCATTERING_IFORMFACTORBORN_H
17 
19 
20 class IShape;
21 struct SlicingEffects; // defined below
22 
23 //! Pure virtual base class for Born form factors.
24 //!
25 //! In contrast to the generic IFormFactor, a Born form factor does not depend
26 //! on the incoming and outgoing wave vectors ki and kf, except through their
27 //! difference, the scattering vector q=ki-kf.
28 
29 //! @ingroup formfactors_internal
30 
32 {
33 public:
35  IFormFactorBorn(const NodeMeta& meta, const std::vector<double>& PValues);
37 
38  IFormFactorBorn* clone() const override = 0;
39 
40  void setAmbientMaterial(const Material&) override {}
41 
42  complex_t evaluate(const WavevectorInfo& wavevectors) const override;
43 
44 #ifndef SWIG
45  Eigen::Matrix2cd evaluatePol(const WavevectorInfo& wavevectors) const override;
46 #endif
47 
48  virtual double bottomZ(const IRotation& rotation) const override;
49  virtual double topZ(const IRotation& rotation) const override;
50 
51  //! Returns scattering amplitude for complex scattering wavevector q=k_i-k_f.
52  //! This method is public only for convenience of plotting form factors in Python.
53  virtual complex_t evaluate_for_q(cvector_t q) const = 0;
54 
55 protected:
56  //! Default implementation only allows rotations along z-axis
57  bool canSliceAnalytically(const IRotation& rot) const override;
58 
59 #ifndef SWIG
60  //! Returns scattering amplitude for complex scattering wavevector q=k_i-k_f in case
61  //! of matrix interactions. Default implementation calls evaluate_for_q(q) and
62  //! multiplies with the unit matrix.
63  virtual Eigen::Matrix2cd evaluate_for_q_pol(cvector_t q) const;
64 #endif
65 
66  //! IShape object, used to retrieve vertices (which may be approximate in the case
67  //! of round shapes). For soft particles, this will be a hard mean shape.
68  std::unique_ptr<IShape> mP_shape;
69 
70  //! Helper method for slicing
71  SlicingEffects computeSlicingEffects(ZLimits limits, const kvector_t& position,
72  double height) const;
73 
74  //! Calculates the z-coordinate of the lowest vertex after rotation
75  static double BottomZ(const std::vector<kvector_t>& vertices, const IRotation& rotation);
76 
77  //! Calculates the z-coordinate of the highest vertex after rotation
78  static double TopZ(const std::vector<kvector_t>& vertices, const IRotation& rotation);
79 };
80 
81 //! Nested structure that holds slicing effects on position and removed parts.
82 
83 //! @ingroup formfactors_internal
84 
87  double dz_bottom;
88  double dz_top;
89 };
90 
91 #ifdef POLYHEDRAL_DIAGNOSTIC
92 //! Information about the latest form factor evaluation. Not thread-safe.
93 //! Used only by external test program.
94 class Diagnosis
95 {
96 public:
97  int maxOrder;
98  int nExpandedFaces;
99  int debmsg;
100  bool request_convergence;
101  bool operator!=(const Diagnosis& other) const
102  {
103  return maxOrder != other.maxOrder || nExpandedFaces != other.nExpandedFaces;
104  }
105  friend std::ostream& operator<<(std::ostream& stream, const Diagnosis& diag)
106  {
107  return stream << " [" << diag.nExpandedFaces << ":" << diag.maxOrder << "]";
108  }
109 };
110 #endif
111 
112 #endif // BORNAGAIN_CORE_SCATTERING_IFORMFACTORBORN_H
std::complex< double > complex_t
Definition: Complex.h:20
Defines and implements pure virtual interface IFormFactor.
bool operator!=(const Material &left, const Material &right)
Comparison operator for material wrapper (inequality check)
Definition: Material.cpp:129
std::ostream & operator<<(std::ostream &os, const RangedDistribution &distribution)
Pure virtual base class for Born form factors.
SlicingEffects computeSlicingEffects(ZLimits limits, const kvector_t &position, double height) const
Helper method for slicing.
bool canSliceAnalytically(const IRotation &rot) const override
Default implementation only allows rotations along z-axis.
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...
IFormFactorBorn * clone() const override=0
Returns a clone of this ISample object.
void setAmbientMaterial(const Material &) override
Passes the material in which this particle is embedded.
virtual complex_t evaluate_for_q(cvector_t q) const =0
Returns scattering amplitude for complex scattering wavevector q=k_i-k_f.
static double TopZ(const std::vector< kvector_t > &vertices, const IRotation &rotation)
Calculates the z-coordinate of the highest vertex after rotation.
virtual double bottomZ(const IRotation &rotation) const override
Returns the z-coordinate of the lowest point in this shape after a given rotation.
complex_t evaluate(const WavevectorInfo &wavevectors) const override
Returns scattering amplitude for complex wavevectors ki, kf.
Eigen::Matrix2cd evaluatePol(const WavevectorInfo &wavevectors) const override
Returns scattering amplitude for matrix interactions.
static double BottomZ(const std::vector< kvector_t > &vertices, const IRotation &rotation)
Calculates the z-coordinate of the lowest vertex after rotation.
std::unique_ptr< IShape > mP_shape
IShape object, used to retrieve vertices (which may be approximate in the case of round shapes).
virtual double topZ(const IRotation &rotation) const override
Returns the z-coordinate of the lowest point in this shape after a given rotation.
Pure virtual base class for all form factors.
Definition: IFormFactor.h:40
Pure virtual interface for rotations.
Definition: Rotations.h:27
Pure virtual base class for different shapes.
Definition: IShape.h:29
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
Nested structure that holds slicing effects on position and removed parts.
kvector_t position