BornAgain  1.19.0
Simulate and fit neutron and x-ray scattering at grazing incidence
IFormFactor.cpp
Go to the documentation of this file.
1 // ************************************************************************************************
2 //
3 // BornAgain: simulate and fit reflection and scattering
4 //
5 //! @file Sample/Scattering/IFormFactor.cpp
6 //! @brief Implements interface class 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 
21 #include <memory>
22 #include <utility>
23 
24 namespace {
25 bool shapeIsContainedInLimits(const IFormFactor& formfactor, ZLimits limits, const IRotation& rot,
26  kvector_t translation)
27 {
28  double zbottom = formfactor.bottomZ(rot) + translation.z();
29  double ztop = formfactor.topZ(rot) + translation.z();
30  OneSidedLimit lower_limit = limits.lowerLimit();
31  OneSidedLimit upper_limit = limits.upperLimit();
32  if (!upper_limit.m_limitless && ztop > upper_limit.m_value)
33  return false;
34  if (!lower_limit.m_limitless && zbottom < lower_limit.m_value)
35  return false;
36  return true;
37 }
38 bool shapeOutsideLimits(const IFormFactor& formfactor, ZLimits limits, const IRotation& rot,
39  kvector_t translation)
40 {
41  double zbottom = formfactor.bottomZ(rot) + translation.z();
42  double ztop = formfactor.topZ(rot) + translation.z();
43  OneSidedLimit lower_limit = limits.lowerLimit();
44  OneSidedLimit upper_limit = limits.upperLimit();
45  if (!upper_limit.m_limitless && zbottom >= upper_limit.m_value)
46  return true;
47  if (!lower_limit.m_limitless && ztop <= lower_limit.m_value)
48  return true;
49  return false;
50 }
51 } // namespace
52 
53 IFormFactor::IFormFactor(const NodeMeta& meta, const std::vector<double>& PValues)
54  : ISampleNode(meta, PValues)
55 {
56 }
57 
59  kvector_t translation) const
60 {
61  if (shapeIsContainedInLimits(*this, limits, rot, translation))
62  return createTransformedFormFactor(*this, rot, translation);
63  if (shapeOutsideLimits(*this, limits, rot, translation))
64  return nullptr;
65  if (canSliceAnalytically(rot))
66  return sliceFormFactor(limits, rot, translation);
67  throw std::runtime_error(getName()
68  + "::createSlicedFormFactor error: not supported for "
69  "the given rotation!");
70 }
71 
72 Eigen::Matrix2cd IFormFactor::evaluatePol(const WavevectorInfo&) const
73 {
74  // Throws to prevent unanticipated behaviour
75  throw std::runtime_error("IFormFactor::evaluatePol: is not implemented by default");
76 }
77 
78 double IFormFactor::volume() const
79 {
80  auto zero_wavevectors = WavevectorInfo::GetZeroQ();
81  return std::abs(evaluate(zero_wavevectors));
82 }
83 
84 void IFormFactor::setSpecularInfo(std::unique_ptr<const ILayerRTCoefficients>,
85  std::unique_ptr<const ILayerRTCoefficients>)
86 {
87 }
88 
90 {
91  return false;
92 }
93 
95 {
96  throw std::runtime_error(getName() + "::sliceFormFactor error: not implemented!");
97 }
98 
100  const IRotation& rot, kvector_t translation)
101 {
102  std::unique_ptr<IFormFactor> P_fftemp, P_result;
103  if (!rot.isIdentity())
104  P_fftemp = std::make_unique<FormFactorDecoratorRotation>(formfactor, rot);
105  else
106  P_fftemp.reset(formfactor.clone());
107  if (translation != kvector_t())
108  P_result = std::make_unique<FormFactorDecoratorPositionFactor>(*P_fftemp, translation);
109  else
110  std::swap(P_fftemp, P_result);
111  return P_result.release();
112 }
Defines class FormFactorDecoratorPositionFactor.
Defines class FormFactorDecoratorRotation.
Defines and implements interface IFormFactor.
Defines and implements class ILayerRTCoefficients.
void swap(OutputDataIterator< TValue, TContainer > &left, OutputDataIterator< TValue, TContainer > &right)
make Swappable
Defines IRotation classes.
BasicVector3D< double > kvector_t
Definition: Vectors3D.h:21
Defines WavevectorInfo.
T z() const
Returns z-component in cartesian coordinate system.
Definition: BasicVector3D.h:67
Abstract base class for all form factors.
Definition: IFormFactor.h:36
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 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
const std::string & getName() const
Abstract base class for rotations.
Definition: Rotations.h:28
virtual bool isIdentity() const
Returns true if rotation matrix is identity matrix (no rotations)
Definition: Rotations.cpp:58
Abstract base class for sample components and properties related to scattering.
Definition: ISampleNode.h:28
Holds all wavevector information relevant for calculating form factors.
static WavevectorInfo GetZeroQ()
Class that contains upper and lower limits of the z-coordinate for the slicing of form factors.
Definition: ZLimits.h:45
OneSidedLimit lowerLimit() const
Definition: ZLimits.cpp:39
OneSidedLimit upperLimit() const
Definition: ZLimits.cpp:44
Metadata of one model node.
Definition: INode.h:38
Helper class that represents a onesided limit.
Definition: ZLimits.h:35
double m_value
Definition: ZLimits.h:37
bool m_limitless
Definition: ZLimits.h:36