BornAgain  1.18.0
Simulate and fit neutron and x-ray scattering at grazing incidence
FormFactorCoreShell.cpp
Go to the documentation of this file.
1 // ************************************************************************** //
2 //
3 // BornAgain: simulate and fit scattering at grazing incidence
4 //
5 //! @file Sample/Particle/FormFactorCoreShell.cpp
6 //! @brief Implements class FormFactorCoreShell.
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 
18  : mP_core(core), mP_shell(shell)
19 {
20  setName("FormFactorCoreShell");
21 }
22 
24 
26 {
27  return new FormFactorCoreShell(mP_core->clone(), mP_shell->clone());
28 }
29 
31 {
32  return mP_shell->radialExtension();
33 }
34 
35 double FormFactorCoreShell::bottomZ(const IRotation& rotation) const
36 {
37  return mP_shell->bottomZ(rotation);
38 }
39 
40 double FormFactorCoreShell::topZ(const IRotation& rotation) const
41 {
42  return mP_shell->topZ(rotation);
43 }
44 
46 {
47  mP_shell->setAmbientMaterial(material);
48 }
49 
51 {
52  return mP_shell->evaluate(wavevectors) + mP_core->evaluate(wavevectors);
53 }
54 
55 Eigen::Matrix2cd FormFactorCoreShell::evaluatePol(const WavevectorInfo& wavevectors) const
56 {
57  return mP_shell->evaluatePol(wavevectors) + mP_core->evaluatePol(wavevectors);
58 }
std::complex< double > complex_t
Definition: Complex.h:20
Defines class FormFactorCoreShell.
Form Factor for a core shell particle.
std::unique_ptr< IFormFactor > mP_core
FormFactorCoreShell * clone() const override final
Returns a clone of this ISample object.
void setAmbientMaterial(const Material &material) override final
Passes the material in which this particle is embedded.
Eigen::Matrix2cd evaluatePol(const WavevectorInfo &wavevectors) const override final
Calculates and returns a polarized form factor calculation in DWBA.
~FormFactorCoreShell() override final
double radialExtension() const override final
Returns the (approximate in some cases) radial size of the particle of this form factor's shape.
double topZ(const IRotation &rotation) const override final
Returns the z-coordinate of the lowest point in this shape after a given rotation.
double bottomZ(const IRotation &rotation) const override final
Returns the z-coordinate of the lowest point in this shape after a given rotation.
FormFactorCoreShell(IFormFactor *core, IFormFactor *shell)
complex_t evaluate(const WavevectorInfo &wavevectors) const override final
Returns scattering amplitude for complex wavevectors ki, kf.
std::unique_ptr< IFormFactor > mP_shell
Pure virtual base class for all form factors.
Definition: IFormFactor.h:40
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.