BornAgain  1.19.79
Open-source research software to simulate and fit neutron and x-ray reflectometry and grazing-incidence small-angle scattering
IReParticle.h
Go to the documentation of this file.
1 // ************************************************************************************************
2 //
3 // BornAgain: simulate and fit reflection and scattering
4 //
5 //! @file Resample/Particle/IReParticle.h
6 //! @brief Defines and implements interface IReParticle.
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 #ifdef SWIG
16 #error no need to expose this header to Swig
17 #endif
18 
19 #ifndef USER_API
20 #ifndef BORNAGAIN_RESAMPLE_PARTICLE_IREPARTICLE_H
21 #define BORNAGAIN_RESAMPLE_PARTICLE_IREPARTICLE_H
22 
23 #include "Base/Types/ICloneable.h"
24 #include <heinz/Complex.h>
25 #include <heinz/Vectors3D.h>
26 
27 class IRotation;
28 class Material;
29 class SpinMatrix;
30 class WavevectorInfo;
31 
32 //! Abstract base class for reprocessed particles.
33 //!
34 //! Reprocessing is necessary to handle particles that cross material layers.
35 //! These particles are divided into several.
36 
37 class IReParticle : public ICloneable {
38 protected:
39  IReParticle() = default;
40 
41 public:
42  ~IReParticle() override = default;
43 
44  IReParticle* clone() const override = 0;
45 
46  //! Passes the material in which this particle is embedded.
47  virtual void setAmbientMaterial(const Material&) {}
48 
49  //! Returns scattering amplitude for complex wavevectors ki, kf.
50  virtual complex_t theFF(const WavevectorInfo& wavevectors) const = 0;
51 
52  //! Returns the (approximate in some cases) radial size of the particle of this
53  //! form factor's shape. This is used for SSCA calculations
54  virtual double radialExtension() const = 0;
55 
56  //! Returns the z-coordinate of the lowest point in this shape after a given rotation
57  virtual double bottomZ(const IRotation* rotation) const = 0;
58 
59  //! Returns the z-coordinate of the lowest point in this shape after a given rotation
60  virtual double topZ(const IRotation* rotation) const = 0;
61 
62  //! Returns scattering amplitude for matrix interactions
63  virtual SpinMatrix thePolFF(const WavevectorInfo& wavevectors) const;
64 
65  //! Returns the total volume of the particle of this form factor's shape
66  virtual double volume() const;
67 };
68 
69 #endif // BORNAGAIN_RESAMPLE_PARTICLE_IREPARTICLE_H
70 #endif // USER_API
Defines and implements the standard mix-in ICloneable.
Interface for polymorphic classes that should not be copied, except by explicit cloning.
Definition: ICloneable.h:23
Abstract base class for reprocessed particles.
Definition: IReParticle.h:37
virtual double topZ(const IRotation *rotation) const =0
Returns the z-coordinate of the lowest point in this shape after a given rotation.
virtual SpinMatrix thePolFF(const WavevectorInfo &wavevectors) const
Returns scattering amplitude for matrix interactions.
Definition: IReParticle.cpp:20
IReParticle * clone() const override=0
IReParticle()=default
virtual complex_t theFF(const WavevectorInfo &wavevectors) const =0
Returns scattering amplitude for complex wavevectors ki, kf.
~IReParticle() override=default
virtual double bottomZ(const IRotation *rotation) const =0
Returns the z-coordinate of the lowest point in this shape after a given rotation.
virtual double volume() const
Returns the total volume of the particle of this form factor's shape.
Definition: IReParticle.cpp:26
virtual void setAmbientMaterial(const Material &)
Passes the material in which this particle is embedded.
Definition: IReParticle.h:47
virtual double radialExtension() const =0
Returns the (approximate in some cases) radial size of the particle of this form factor's shape....
Abstract base class for rotations.
Definition: Rotations.h:29
A wrapper for underlying material implementation.
Definition: Material.h:35
Holds all wavevector information relevant for calculating form factors.