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