BornAgain  1.19.0
Simulate and fit neutron and x-ray scattering at grazing incidence
FormFactorWeighted.h
Go to the documentation of this file.
1 // ************************************************************************************************
2 //
3 // BornAgain: simulate and fit reflection and scattering
4 //
5 //! @file Sample/Particle/FormFactorWeighted.h
6 //! @brief Defines 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 
15 #ifndef BORNAGAIN_SAMPLE_PARTICLE_FORMFACTORWEIGHTED_H
16 #define BORNAGAIN_SAMPLE_PARTICLE_FORMFACTORWEIGHTED_H
17 
19 
20 //! Coherent sum of different scalar IFormFactor%s with different weights.
21 //!
22 //! Used by ParticleComposition.
23 //! If same particles are at different positions, then consider
24 //! FormFactorDecoratorMultiPositionFactor (restore from commit 0500a26de76).
25 
26 //! @ingroup formfactors_internal
27 
29 public:
31  ~FormFactorWeighted() override;
32 
33  FormFactorWeighted* clone() const override;
34 
35  void accept(INodeVisitor* visitor) const override { visitor->visit(this); }
36 
37  double radialExtension() const override;
38 
39  double bottomZ(const IRotation& rotation) const override;
40 
41  double topZ(const IRotation& rotation) const override;
42 
43  void addFormFactor(const IFormFactor& form_factor, double weight = 1.0);
44 
45  void setAmbientMaterial(const Material& material) override;
46 
47  complex_t evaluate(const WavevectorInfo& wavevectors) const override;
48 
49 #ifndef SWIG
50  //! Calculates and returns a polarized form factor calculation in DWBA
51  Eigen::Matrix2cd evaluatePol(const WavevectorInfo& wavevectors) const override;
52 #endif
53 
54 protected:
55  std::vector<IFormFactor*> m_form_factors;
56  std::vector<double> m_weights;
57 };
58 
59 #endif // BORNAGAIN_SAMPLE_PARTICLE_FORMFACTORWEIGHTED_H
std::complex< double > complex_t
Definition: Complex.h:20
Defines and implements interface IFormFactor.
Coherent sum of different scalar IFormFactors with different weights.
void accept(INodeVisitor *visitor) const override
Calls the INodeVisitor's visit method.
void setAmbientMaterial(const Material &material) override
Passes the material in which this particle is embedded.
double radialExtension() const override
Returns the (approximate in some cases) radial size of the particle of this form factor's shape.
complex_t evaluate(const WavevectorInfo &wavevectors) const override
Returns scattering amplitude for complex wavevectors ki, kf.
FormFactorWeighted * clone() const override
Returns a clone of this ISampleNode object.
Eigen::Matrix2cd evaluatePol(const WavevectorInfo &wavevectors) const override
Calculates and returns a polarized form factor calculation in DWBA.
void addFormFactor(const IFormFactor &form_factor, double weight=1.0)
std::vector< IFormFactor * > m_form_factors
double topZ(const IRotation &rotation) const override
Returns the z-coordinate of the lowest point in this shape after a given rotation.
std::vector< double > m_weights
double bottomZ(const IRotation &rotation) const override
Returns the z-coordinate of the lowest point in this shape after a given rotation.
Abstract base class for all form factors.
Definition: IFormFactor.h:36
Visitor interface to visit ISampleNode objects.
Definition: INodeVisitor.h:146
virtual void visit(const BasicLattice2D *)
Definition: INodeVisitor.h:151
Abstract base class for rotations.
Definition: Rotations.h:28
virtual const Material * material() const
Returns nullptr, unless overwritten to return a specific material.
Definition: ISampleNode.h:37
A wrapper for underlying material implementation.
Definition: Material.h:29
Holds all wavevector information relevant for calculating form factors.