BornAgain  1.19.0
Simulate and fit neutron and x-ray scattering at grazing incidence
FormFactorCrystal.h
Go to the documentation of this file.
1 // ************************************************************************************************
2 //
3 // BornAgain: simulate and fit reflection and scattering
4 //
5 //! @file Sample/Particle/FormFactorCrystal.h
6 //! @brief Defines class FormFactorCrystal.
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_FORMFACTORCRYSTAL_H
16 #define BORNAGAIN_SAMPLE_PARTICLE_FORMFACTORCRYSTAL_H
17 
20 
21 //! The form factor of a MesoCrystal.
22 //! @ingroup formfactors
23 
25 public:
26  FormFactorCrystal(const Lattice3D& lattice, const IFormFactor& basis_form_factor,
27  const IFormFactor& meso_form_factor, double position_variance = 0.0);
28  ~FormFactorCrystal() override;
29 
30  FormFactorCrystal* clone() const override
31  {
34  }
35 
36  void accept(INodeVisitor* visitor) const override { visitor->visit(this); }
37 
38  void setAmbientMaterial(const Material& material) override
39  {
41  }
42 
43  double volume() const override { return m_meso_form_factor->volume(); }
44  double radialExtension() const override { return m_meso_form_factor->radialExtension(); }
45 
46  double bottomZ(const IRotation& rotation) const override;
47 
48  double topZ(const IRotation& rotation) const override;
49 
50  complex_t evaluate(const WavevectorInfo& wavevectors) const override;
51 #ifndef SWIG
52  Eigen::Matrix2cd evaluatePol(const WavevectorInfo& wavevectors) const override;
53 #endif
54 
55 private:
57  complex_t debyeWallerFactor(const kvector_t& q_i) const;
58 
61  IFormFactor* m_meso_form_factor; //!< The outer shape of this mesocrystal
64 };
65 
66 #endif // BORNAGAIN_SAMPLE_PARTICLE_FORMFACTORCRYSTAL_H
std::complex< double > complex_t
Definition: Complex.h:20
Defines interface IBornFF.
Defines class Lattice.
The form factor of a MesoCrystal.
IFormFactor * m_meso_form_factor
The outer shape of this mesocrystal.
void accept(INodeVisitor *visitor) const override
Calls the INodeVisitor's visit method.
double bottomZ(const IRotation &rotation) const override
Returns the z-coordinate of the lowest point in this shape after a given rotation.
Eigen::Matrix2cd evaluatePol(const WavevectorInfo &wavevectors) const override
Returns scattering amplitude for matrix interactions.
double volume() const override
Returns the total volume of the particle of this form factor's shape.
double topZ(const IRotation &rotation) const override
Returns the z-coordinate of the lowest point in this shape after a given rotation.
FormFactorCrystal(const Lattice3D &lattice, const IFormFactor &basis_form_factor, const IFormFactor &meso_form_factor, double position_variance=0.0)
~FormFactorCrystal() override
IFormFactor * m_basis_form_factor
complex_t evaluate(const WavevectorInfo &wavevectors) const override
Returns scattering amplitude for complex wavevectors ki, kf.
void setAmbientMaterial(const Material &material) override
Passes the material in which this particle is embedded.
void calculateLargestReciprocalDistance()
double radialExtension() const override
Returns the (approximate in some cases) radial size of the particle of this form factor's shape.
FormFactorCrystal * clone() const override
Returns a clone of this ISampleNode object.
complex_t debyeWallerFactor(const kvector_t &q_i) const
Abstract base class for all form factors.
Definition: IFormFactor.h:36
virtual double radialExtension() const =0
Returns the (approximate in some cases) radial size of the particle of this form factor's shape.
virtual void setAmbientMaterial(const Material &)=0
Passes the material in which this particle is embedded.
virtual double volume() const
Returns the total volume of the particle of this form factor's shape.
Definition: IFormFactor.cpp:78
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 Bravais lattice, characterized by three basis vectors, and optionally an ISelectionRule.
Definition: Lattice3D.h:29
A wrapper for underlying material implementation.
Definition: Material.h:29
Holds all wavevector information relevant for calculating form factors.