BornAgain  1.18.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 scattering at grazing incidence
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_CORE_PARTICLE_FORMFACTORCRYSTAL_H
16 #define BORNAGAIN_CORE_PARTICLE_FORMFACTORCRYSTAL_H
17 
18 #include "Sample/Lattice/Lattice.h"
20 
21 //! The form factor of a MesoCrystal.
22 //! @ingroup formfactors
23 
25 {
26 public:
27  FormFactorCrystal(const Lattice& lattice, const IFormFactor& basis_form_factor,
28  const IFormFactor& meso_form_factor, double position_variance = 0.0);
29  ~FormFactorCrystal() override final;
30 
31  FormFactorCrystal* clone() const override final
32  {
33  return new FormFactorCrystal(m_lattice, *mp_basis_form_factor, *mp_meso_form_factor,
34  m_position_variance);
35  }
36 
37  void accept(INodeVisitor* visitor) const override final { visitor->visit(this); }
38 
39  void setAmbientMaterial(const Material& material) override
40  {
41  mp_basis_form_factor->setAmbientMaterial(material);
42  }
43 
44  double volume() const override final { return mp_meso_form_factor->volume(); }
45  double radialExtension() const override final { return mp_meso_form_factor->radialExtension(); }
46 
47  double bottomZ(const IRotation& rotation) const override;
48 
49  double topZ(const IRotation& rotation) const override final;
50 
51  complex_t evaluate(const WavevectorInfo& wavevectors) const override final;
52 #ifndef SWIG
53  Eigen::Matrix2cd evaluatePol(const WavevectorInfo& wavevectors) const override final;
54 #endif
55 
56 private:
57  void calculateLargestReciprocalDistance();
58  complex_t debyeWallerFactor(const kvector_t& q_i) const;
59 
60  Lattice m_lattice;
61  IFormFactor* mp_basis_form_factor;
62  IFormFactor* mp_meso_form_factor; //!< The outer shape of this mesocrystal
63  double m_position_variance;
64  double m_max_rec_length;
65 };
66 
67 #endif // BORNAGAIN_CORE_PARTICLE_FORMFACTORCRYSTAL_H
Defines pure virtual interface class IFormFactorBorn.
Defines class Lattice.
The form factor of a MesoCrystal.
Eigen::Matrix2cd evaluatePol(const WavevectorInfo &wavevectors) const override final
Returns scattering amplitude for matrix interactions.
FormFactorCrystal * clone() const override final
Returns a clone of this ISample object.
double bottomZ(const IRotation &rotation) const override
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.
double radialExtension() const override final
Returns the (approximate in some cases) radial size of the particle of this form factor's shape.
double volume() const override final
Returns the total volume 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.
void accept(INodeVisitor *visitor) const override final
Calls the INodeVisitor's visit method.
void setAmbientMaterial(const Material &material) override
Passes the material in which this particle is embedded.
Pure virtual base class for all form factors.
Definition: IFormFactor.h:40
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:56
Visitor interface to visit ISample objects.
Definition: INodeVisitor.h:149
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 lattice with three basis vectors.
Definition: Lattice.h:28
A wrapper for underlying material implementation.
Definition: Material.h:29
Holds all wavevector information relevant for calculating form factors.