BornAgain  1.18.0
Simulate and fit neutron and x-ray scattering at grazing incidence
FormFactorWeighted.cpp
Go to the documentation of this file.
1 // ************************************************************************** //
2 //
3 // BornAgain: simulate and fit scattering at grazing incidence
4 //
5 //! @file Sample/Particle/FormFactorWeighted.cpp
6 //! @brief Implements class FormFactorWeighted.
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 
16 #include "Base/Utils/Algorithms.h"
17 
19 {
20  setName("FormFactorWeighted");
21 }
22 
24 {
25  for (size_t index = 0; index < m_form_factors.size(); ++index)
26  delete m_form_factors[index];
27 }
28 
30 {
31  FormFactorWeighted* result = new FormFactorWeighted();
32  for (size_t index = 0; index < m_form_factors.size(); ++index)
33  result->addFormFactor(*m_form_factors[index], m_weights[index]);
34  return result;
35 }
36 
38 {
39  double result{0.0};
40  for (size_t index = 0; index < m_form_factors.size(); ++index)
41  result += m_weights[index] * m_form_factors[index]->radialExtension();
42  return result;
43 }
44 
45 double FormFactorWeighted::bottomZ(const IRotation& rotation) const
46 {
47  if (m_form_factors.empty())
48  throw std::runtime_error("FormFactorWeighted::bottomZ() -> Error: "
49  "'this' contains no form factors.");
50  return algo::min_value(m_form_factors.begin(), m_form_factors.end(),
51  [&rotation](IFormFactor* ff) { return ff->bottomZ(rotation); });
52 }
53 
54 double FormFactorWeighted::topZ(const IRotation& rotation) const
55 {
56  if (m_form_factors.empty())
57  throw std::runtime_error("FormFactorWeighted::topZ() -> Error: "
58  "'this' contains no form factors.");
59  return algo::max_value(m_form_factors.begin(), m_form_factors.end(),
60  [&rotation](IFormFactor* ff) { return ff->topZ(rotation); });
61 }
62 
63 void FormFactorWeighted::addFormFactor(const IFormFactor& form_factor, double weight)
64 {
65  m_form_factors.push_back(form_factor.clone());
66  m_weights.push_back(weight);
67 }
68 
70 {
71  for (size_t index = 0; index < m_form_factors.size(); ++index)
73 }
74 
76 {
77  complex_t result(0.0, 0.0);
78  for (size_t index = 0; index < m_form_factors.size(); ++index)
79  result += m_weights[index] * m_form_factors[index]->evaluate(wavevectors);
80  return result;
81 }
82 
83 Eigen::Matrix2cd FormFactorWeighted::evaluatePol(const WavevectorInfo& wavevectors) const
84 {
85  Eigen::Matrix2cd result = Eigen::Matrix2cd::Zero();
86  for (size_t index = 0; index < m_form_factors.size(); ++index)
87  result += m_weights[index] * m_form_factors[index]->evaluatePol(wavevectors);
88  return result;
89 }
Defines and implements namespace algo with some algorithms.
std::complex< double > complex_t
Definition: Complex.h:20
Defines class FormFactorWeighted.
Coherent sum of different scalar IFormFactor's with different weights.
void setAmbientMaterial(const Material &material) override final
Passes the material in which this particle is embedded.
double bottomZ(const IRotation &rotation) const override final
Returns the z-coordinate of the lowest point in this shape after a given rotation.
~FormFactorWeighted() override final
void addFormFactor(const IFormFactor &form_factor, double weight=1.0)
Eigen::Matrix2cd evaluatePol(const WavevectorInfo &wavevectors) const override final
Calculates and returns a polarized form factor calculation in DWBA.
double topZ(const IRotation &rotation) const override final
Returns the z-coordinate of the lowest point in this shape after a given rotation.
complex_t evaluate(const WavevectorInfo &wavevectors) const override final
Returns scattering amplitude for complex wavevectors ki, kf.
std::vector< IFormFactor * > m_form_factors
std::vector< double > m_weights
FormFactorWeighted * clone() const override final
Returns a clone of this ISample object.
double radialExtension() const override final
Returns the (approximate in some cases) radial size of the particle of this form factor's shape.
Pure virtual base class for all form factors.
Definition: IFormFactor.h:40
IFormFactor * clone() const override=0
Returns a clone of this ISample object.
void setName(const std::string &name)
Pure virtual interface for rotations.
Definition: Rotations.h:27
virtual const Material * material() const
Returns nullptr, unless overwritten to return a specific material.
Definition: ISample.h:37
A wrapper for underlying material implementation.
Definition: Material.h:29
Holds all wavevector information relevant for calculating form factors.
double max_value(const Iterator &begin, const Iterator &end, const Evaluator &evaluate)
Returns the maximum value of function evaluate as applied to the elements of an iterator range.
Definition: Algorithms.h:65
double min_value(const Iterator &begin, const Iterator &end, const Evaluator &evaluate)
Returns the minimum value of function evaluate as applied to the elements of an iterator range.
Definition: Algorithms.h:54