BornAgain  1.19.79
Open-source research software to simulate and fit neutron and x-ray reflectometry and grazing-incidence small-angle scattering
HemiEllipsoid Class Reference

Description

An hemi ellipsoid, obtained by truncating a full ellipsoid in the middle plane spanned by two principal axes.

Definition at line 24 of file HemiEllipsoid.h.

Inheritance diagram for HemiEllipsoid:
[legend]
Collaboration diagram for HemiEllipsoid:
[legend]

Public Member Functions

 HemiEllipsoid (double radius_x, double radius_y, double height)
 
 HemiEllipsoid (std::vector< double > P)
 
 ~HemiEllipsoid () override=default
 
virtual double bottomZ (const IRotation *rotation) const
 
virtual bool canSliceAnalytically (const IRotation *rot) const
 Default implementation only allows rotations along z-axis. More...
 
void checkNodeArgs () const
 Raises exception if a parameter value is invalid. More...
 
std::string className () const final
 Returns the class name, to be hard-coded in each leaf class that inherits from INode. More...
 
HemiEllipsoidclone () const override
 Returns a clone of this ISampleNode object. More...
 
std::vector< const Material * > containedMaterials () const
 Returns set of unique materials contained in this ISampleNode. More...
 
complex_t formfactor_at_bottom (C3 q) const override
 
virtual SpinMatrix formfactor_pol (C3 q) const
 Returns scattering amplitude for complex scattering wavevector q=k_i-k_f in case of matrix interactions. Default implementation calls formfactor_at_bottom(q) and multiplies with the unit matrix. More...
 
double height () const
 
bool isMagnetic () const
 Returns true if there is any magnetic material in this ISampleNode. More...
 
virtual const Materialmaterial () const
 Returns nullptr, unless overwritten to return a specific material. More...
 
virtual std::vector< const INode * > nodeChildren () const
 Returns all children. More...
 
std::vector< const INode * > nodeOffspring () const
 Returns all descendants. More...
 
std::vector< ParaMetaparDefs () const final
 Returns the parameter definitions, to be hard-coded in each leaf class. More...
 
virtual std::string pythonConstructor () const
 Creates the Python constructor of this class (or derived classes) More...
 
double radialExtension () const override
 Returns the (approximate in some cases) radial size of the particle of this form factor's shape. This is used for SSCA calculations. More...
 
double radiusX () const
 
double radiusY () const
 
std::string shapeName () const
 
virtual complex_t theFF (const WavevectorInfo &wavevectors) const
 
virtual SpinMatrix thePolFF (const WavevectorInfo &wavevectors) const
 
virtual double topZ (const IRotation *rotation) const
 
virtual void transferToCPP ()
 Used for Python overriding of clone (see swig/tweaks.py) More...
 
virtual double volume () const
 

Protected Attributes

std::vector< double > m_P
 
std::unique_ptr< IShape3Dm_shape3D
 IShape3D object, used to retrieve vertices (which may be approximate in the case of round shapes). For soft particles, this will be a hard mean shape. More...
 

Private Attributes

const double & m_height
 
const double & m_radius_x
 
const double & m_radius_y
 

Constructor & Destructor Documentation

◆ HemiEllipsoid() [1/2]

HemiEllipsoid::HemiEllipsoid ( double  radius_x,
double  radius_y,
double  height 
)

Definition at line 33 of file HemiEllipsoid.cpp.

34  : HemiEllipsoid(std::vector<double>{radius_x, radius_y, height})
35 {
36 }
double height() const
Definition: HemiEllipsoid.h:44
HemiEllipsoid(double radius_x, double radius_y, double height)

References height().

Referenced by clone().

Here is the call graph for this function:

◆ HemiEllipsoid() [2/2]

HemiEllipsoid::HemiEllipsoid ( std::vector< double >  P)

Definition at line 22 of file HemiEllipsoid.cpp.

23  : IFormFactor(P)
24  , m_radius_x(m_P[0])
25  , m_radius_y(m_P[1])
26  , m_height(m_P[2])
27 {
28  checkNodeArgs();
29  m_shape3D =
30  std::make_unique<TruncatedEllipsoidNet>(m_radius_x, m_radius_x, m_height, m_height, 0.0);
31 }
const double & m_height
Definition: HemiEllipsoid.h:55
const double & m_radius_x
Definition: HemiEllipsoid.h:53
const double & m_radius_y
Definition: HemiEllipsoid.h:54
std::unique_ptr< IShape3D > m_shape3D
IShape3D object, used to retrieve vertices (which may be approximate in the case of round shapes)....
Definition: IFormFactor.h:74
void checkNodeArgs() const
Raises exception if a parameter value is invalid.
Definition: INode.cpp:27
std::vector< double > m_P
Definition: INode.h:63

References INode::checkNodeArgs(), m_height, m_radius_x, and IFormFactor::m_shape3D.

Here is the call graph for this function:

◆ ~HemiEllipsoid()

HemiEllipsoid::~HemiEllipsoid ( )
overridedefault

Member Function Documentation

◆ bottomZ()

double IFormFactor::bottomZ ( const IRotation rotation) const
virtualinherited

Reimplemented in Sphere, IFormFactorPrism, and IFormFactorPolyhedron.

Definition at line 50 of file IFormFactor.cpp.

51 {
52  if (!m_shape3D)
53  throw std::runtime_error("Bug: Form factor has no m_shape3D, cannot compute bottom z");
54  return PolyhedralUtil::BottomZ(m_shape3D->vertices(), rotation);
55 }
double BottomZ(const std::vector< R3 > &vertices, const IRotation *rotation)
Calculates the z-coordinate of the lowest vertex after rotation.

References PolyhedralUtil::BottomZ(), and IFormFactor::m_shape3D.

Here is the call graph for this function:

◆ canSliceAnalytically()

bool IFormFactor::canSliceAnalytically ( const IRotation rot) const
virtualinherited

Default implementation only allows rotations along z-axis.

Reimplemented in Sphere.

Definition at line 64 of file IFormFactor.cpp.

65 {
66  return !rotation || rotation->zInvariant();
67 }

References IRotation::zInvariant().

Here is the call graph for this function:

◆ checkNodeArgs()

void INode::checkNodeArgs ( ) const
inherited

Raises exception if a parameter value is invalid.

Definition at line 27 of file INode.cpp.

28 {
29  size_t nP = m_P.size();
30  if (parDefs().size() != nP) {
31  std::cerr << "BUG in class " << className() << std::endl;
32  std::cerr << "#m_P = " << nP << std::endl;
33  std::cerr << "#PDf = " << parDefs().size() << std::endl;
34  for (const ParaMeta& pm : parDefs())
35  std::cerr << " PDf: " << pm.name << std::endl;
36  ASSERT(0);
37  }
38  ASSERT(parDefs().size() == nP);
39  for (size_t i = 0; i < nP; ++i) {
40  const ParaMeta pm = parDefs()[i];
41 
43  if (pm.vMin == -INF) {
44  ASSERT(pm.vMax == +INF);
45  // nothing to do
46  } else if (pm.vMax == +INF) {
47  ASSERT(pm.vMin == 0);
48  limits = RealLimits::nonnegative();
49  } else {
50  limits = RealLimits::limited(pm.vMin, pm.vMax);
51  }
52  limits.check(pm.name, m_P[i]);
53  }
54 }
#define ASSERT(condition)
Definition: Assert.h:45
const double INF
Definition: INode.h:26
virtual std::vector< ParaMeta > parDefs() const
Returns the parameter definitions, to be hard-coded in each leaf class.
Definition: INode.h:51
virtual std::string className() const =0
Returns the class name, to be hard-coded in each leaf class that inherits from INode.
Limits for a real fit parameter.
Definition: RealLimits.h:24
static RealLimits limitless()
Creates an object without bounds (default)
Definition: RealLimits.cpp:139
void check(const std::string &name, double value) const
Throws if value is outside limits. Parameter 'name' is for exception message.
Definition: RealLimits.cpp:170
static RealLimits nonnegative()
Creates an object which can have only positive values with 0. included.
Definition: RealLimits.cpp:124
static RealLimits limited(double left_bound_value, double right_bound_value)
Creates an object bounded from the left and right.
Definition: RealLimits.cpp:134
Metadata of one model parameter.
Definition: INode.h:29
double vMin
Definition: INode.h:33
double vMax
Definition: INode.h:34
std::string name
Definition: INode.h:30

References ASSERT, RealLimits::check(), INode::className(), INF, RealLimits::limited(), RealLimits::limitless(), INode::m_P, ParaMeta::name, RealLimits::nonnegative(), INode::parDefs(), ParaMeta::vMax, and ParaMeta::vMin.

Referenced by BarGauss::BarGauss(), BarLorentz::BarLorentz(), Bipyramid4::Bipyramid4(), Box::Box(), CantellatedCube::CantellatedCube(), Cone::Cone(), ConstantBackground::ConstantBackground(), CosineRippleBox::CosineRippleBox(), CosineRippleGauss::CosineRippleGauss(), CosineRippleLorentz::CosineRippleLorentz(), Cylinder::Cylinder(), DistributionCosine::DistributionCosine(), DistributionGate::DistributionGate(), DistributionGaussian::DistributionGaussian(), DistributionLogNormal::DistributionLogNormal(), DistributionLorentz::DistributionLorentz(), DistributionTrapezoid::DistributionTrapezoid(), Dodecahedron::Dodecahedron(), EllipsoidalCylinder::EllipsoidalCylinder(), FootprintGauss::FootprintGauss(), FootprintSquare::FootprintSquare(), FuzzySphere::FuzzySphere(), GaussSphere::GaussSphere(), HemiEllipsoid(), HollowSphere::HollowSphere(), HorizontalCylinder::HorizontalCylinder(), Icosahedron::Icosahedron(), LongBoxGauss::LongBoxGauss(), LongBoxLorentz::LongBoxLorentz(), PlatonicOctahedron::PlatonicOctahedron(), PlatonicTetrahedron::PlatonicTetrahedron(), Prism3::Prism3(), Prism6::Prism6(), Profile1DCauchy::Profile1DCauchy(), Profile1DCosine::Profile1DCosine(), Profile1DGate::Profile1DGate(), Profile1DGauss::Profile1DGauss(), Profile1DTriangle::Profile1DTriangle(), Profile1DVoigt::Profile1DVoigt(), Profile2DCauchy::Profile2DCauchy(), Profile2DCone::Profile2DCone(), Profile2DGate::Profile2DGate(), Profile2DGauss::Profile2DGauss(), Profile2DVoigt::Profile2DVoigt(), Pyramid2::Pyramid2(), Pyramid3::Pyramid3(), Pyramid4::Pyramid4(), Pyramid6::Pyramid6(), RotationEuler::RotationEuler(), RotationX::RotationX(), RotationY::RotationY(), RotationZ::RotationZ(), SawtoothRippleBox::SawtoothRippleBox(), SawtoothRippleGauss::SawtoothRippleGauss(), SawtoothRippleLorentz::SawtoothRippleLorentz(), Sphere::Sphere(), Spheroid::Spheroid(), TruncatedCube::TruncatedCube(), TruncatedSphere::TruncatedSphere(), and TruncatedSpheroid::TruncatedSpheroid().

Here is the call graph for this function:

◆ className()

std::string HemiEllipsoid::className ( ) const
inlinefinalvirtual

Returns the class name, to be hard-coded in each leaf class that inherits from INode.

Implements INode.

Definition at line 35 of file HemiEllipsoid.h.

35 { return "HemiEllipsoid"; }

◆ clone()

HemiEllipsoid* HemiEllipsoid::clone ( ) const
inlineoverridevirtual

Returns a clone of this ISampleNode object.

Implements IFormFactor.

Definition at line 31 of file HemiEllipsoid.h.

32  {
34  }

References HemiEllipsoid(), m_height, m_radius_x, and m_radius_y.

Here is the call graph for this function:

◆ containedMaterials()

std::vector< const Material * > ISampleNode::containedMaterials ( ) const
inherited

Returns set of unique materials contained in this ISampleNode.

Definition at line 25 of file ISampleNode.cpp.

26 {
27  std::vector<const Material*> result;
28  if (const Material* p_material = material())
29  result.push_back(p_material);
30  for (const auto* child : nodeChildren()) {
31  if (const auto* sample = dynamic_cast<const ISampleNode*>(child)) {
32  for (const Material* p_material : sample->containedMaterials())
33  result.push_back(p_material);
34  }
35  }
36  return result;
37 }
virtual std::vector< const INode * > nodeChildren() const
Returns all children.
Definition: INode.cpp:56
Abstract base class for sample components and properties related to scattering.
Definition: ISampleNode.h:27
virtual const Material * material() const
Returns nullptr, unless overwritten to return a specific material.
Definition: ISampleNode.h:36
A wrapper for underlying material implementation.
Definition: Material.h:35

References ISampleNode::material(), and INode::nodeChildren().

Referenced by SampleUtils::Multilayer::ContainsCompatibleMaterials(), SampleToPython::initLabels(), and ISampleNode::isMagnetic().

Here is the call graph for this function:

◆ formfactor_at_bottom()

complex_t HemiEllipsoid::formfactor_at_bottom ( C3  q) const
overridevirtual

Implements IFormFactor.

Definition at line 43 of file HemiEllipsoid.cpp.

44 {
45  double R = m_radius_x;
46  double W = m_radius_y;
47  double H = m_height;
48 
49  if (std::abs(q.mag()) <= std::numeric_limits<double>::epsilon())
50  return M_TWOPI * R * W * H / 3.;
51 
52  return M_TWOPI
54  [=](double Z) {
55  double R = m_radius_x;
56  double W = m_radius_y;
57  double H = m_height;
58 
59  double Rz = R * std::sqrt(1.0 - Z * Z / (H * H));
60  double Wz = W * std::sqrt(1.0 - Z * Z / (H * H));
61 
62  complex_t qxRz = q.x() * Rz;
63  complex_t qyWz = q.y() * Wz;
64 
65  complex_t gamma = std::sqrt(qxRz * qxRz + qyWz * qyWz);
66  complex_t J1_gamma_div_gamma = Math::Bessel::J1c(gamma);
67 
68  return Rz * Wz * J1_gamma_div_gamma * exp_I(q.z() * Z);
69  },
70  0., H);
71 }
#define M_TWOPI
Definition: Constants.h:54
To integrate a complex function of a real variable.
Definition: IntegratorGK.h:44
complex_t integrate(const std::function< complex_t(double)> &f, double lmin, double lmax)
double J1c(double x)
Bessel function J1(x)/x.
Definition: Bessel.cpp:172
double gamma(double x)
constexpr Double_t H()
Definition: TMath.h:140
constexpr Double_t R()
Definition: TMath.h:213

References ROOT::Math::Cephes::gamma(), TMath::H(), ComplexIntegrator::integrate(), Math::Bessel::J1c(), m_height, m_radius_x, m_radius_y, M_TWOPI, and TMath::R().

Here is the call graph for this function:

◆ formfactor_pol()

SpinMatrix IFormFactor::formfactor_pol ( C3  q) const
virtualinherited

Returns scattering amplitude for complex scattering wavevector q=k_i-k_f in case of matrix interactions. Default implementation calls formfactor_at_bottom(q) and multiplies with the unit matrix.

Definition at line 78 of file IFormFactor.cpp.

79 {
81 }
virtual complex_t formfactor_at_bottom(C3 q) const =0
static SpinMatrix One()
Definition: SpinMatrix.cpp:36

References IFormFactor::formfactor_at_bottom(), and SpinMatrix::One().

Referenced by IFormFactor::thePolFF().

Here is the call graph for this function:

◆ height()

double HemiEllipsoid::height ( ) const
inline

Definition at line 44 of file HemiEllipsoid.h.

44 { return m_height; }

References m_height.

Referenced by HemiEllipsoid().

◆ isMagnetic()

bool ISampleNode::isMagnetic ( ) const
inherited

Returns true if there is any magnetic material in this ISampleNode.

Definition at line 39 of file ISampleNode.cpp.

40 {
41  const auto materials = containedMaterials();
42  return std::any_of(materials.cbegin(), materials.cend(),
43  [](const Material* mat) { return mat->isMagneticMaterial(); });
44 }
std::vector< const Material * > containedMaterials() const
Returns set of unique materials contained in this ISampleNode.
Definition: ISampleNode.cpp:25

References ISampleNode::containedMaterials().

Referenced by reSample::make().

Here is the call graph for this function:

◆ material()

virtual const Material* ISampleNode::material ( ) const
inlinevirtualinherited

Returns nullptr, unless overwritten to return a specific material.

Reimplemented in Particle, and Layer.

Definition at line 36 of file ISampleNode.h.

36 { return nullptr; }

Referenced by ISampleNode::containedMaterials().

◆ nodeChildren()

◆ nodeOffspring()

std::vector< const INode * > INode::nodeOffspring ( ) const
inherited

Returns all descendants.

Definition at line 61 of file INode.cpp.

62 {
63  std::vector<const INode*> result;
64  result.push_back(this);
65  for (const auto* child : nodeChildren()) {
66  for (const auto* p : child->nodeOffspring())
67  result.push_back(p);
68  }
69  return result;
70 }

References INode::nodeChildren().

Here is the call graph for this function:

◆ parDefs()

std::vector<ParaMeta> HemiEllipsoid::parDefs ( ) const
inlinefinalvirtual

Returns the parameter definitions, to be hard-coded in each leaf class.

Reimplemented from INode.

Definition at line 37 of file HemiEllipsoid.h.

38  {
39  return {{"RadiusX", "nm", "radius in x direction", 0, +INF, 0},
40  {"RadiusY", "nm", "radius in y direction", 0, +INF, 0},
41  {"Height", "nm", "height = radius in z direction", 0, +INF, 0}};
42  }

References INF.

◆ pythonConstructor()

std::string IFormFactor::pythonConstructor ( ) const
virtualinherited

Creates the Python constructor of this class (or derived classes)

Definition at line 69 of file IFormFactor.cpp.

70 {
71  std::vector<std::pair<double, std::string>> arguments;
72  for (size_t i = 0; i < parDefs().size(); i++)
73  arguments.emplace_back(m_P[i], parDefs()[i].unit);
74 
75  return Py::Fmt::printFunction(className(), arguments);
76 }
std::string printFunction(const std::string &name, const std::vector< std::pair< double, std::string >> &arguments)
Print a function in the form "<name>(<arguments>)". arguments will be processed by printArguments(),...
Definition: PyFmt.cpp:168

References INode::className(), INode::m_P, INode::parDefs(), and Py::Fmt::printFunction().

Here is the call graph for this function:

◆ radialExtension()

double HemiEllipsoid::radialExtension ( ) const
overridevirtual

Returns the (approximate in some cases) radial size of the particle of this form factor's shape. This is used for SSCA calculations.

Implements IFormFactor.

Definition at line 38 of file HemiEllipsoid.cpp.

39 {
40  return (m_radius_x + m_radius_y) / 2.0;
41 }

References m_radius_x, and m_radius_y.

◆ radiusX()

double HemiEllipsoid::radiusX ( ) const
inline

Definition at line 45 of file HemiEllipsoid.h.

45 { return m_radius_x; }

References m_radius_x.

◆ radiusY()

double HemiEllipsoid::radiusY ( ) const
inline

Definition at line 46 of file HemiEllipsoid.h.

46 { return m_radius_y; }

References m_radius_y.

◆ shapeName()

std::string IFormFactor::shapeName ( ) const
inherited

Definition at line 33 of file IFormFactor.cpp.

34 {
35  if (className().substr(0, 10) == "FormFactor")
36  return className().substr(10);
37  return className();
38 }

References INode::className().

Here is the call graph for this function:

◆ theFF()

complex_t IFormFactor::theFF ( const WavevectorInfo wavevectors) const
virtualinherited

Definition at line 40 of file IFormFactor.cpp.

41 {
42  return formfactor_at_bottom(wavevectors.getQ());
43 }
C3 getQ() const

References IFormFactor::formfactor_at_bottom(), and WavevectorInfo::getQ().

Here is the call graph for this function:

◆ thePolFF()

SpinMatrix IFormFactor::thePolFF ( const WavevectorInfo wavevectors) const
virtualinherited

Definition at line 45 of file IFormFactor.cpp.

46 {
47  return formfactor_pol(wavevectors.getQ());
48 }
virtual SpinMatrix formfactor_pol(C3 q) const
Returns scattering amplitude for complex scattering wavevector q=k_i-k_f in case of matrix interactio...
Definition: IFormFactor.cpp:78

References IFormFactor::formfactor_pol(), and WavevectorInfo::getQ().

Here is the call graph for this function:

◆ topZ()

double IFormFactor::topZ ( const IRotation rotation) const
virtualinherited

Reimplemented in Sphere, IFormFactorPrism, and IFormFactorPolyhedron.

Definition at line 57 of file IFormFactor.cpp.

58 {
59  if (!m_shape3D)
60  throw std::runtime_error("Bug: Form factor has no m_shape3D, cannot compute top z");
61  return PolyhedralUtil::TopZ(m_shape3D->vertices(), rotation);
62 }
double TopZ(const std::vector< R3 > &vertices, const IRotation *rotation)
Calculates the z-coordinate of the highest vertex after rotation.

References IFormFactor::m_shape3D, and PolyhedralUtil::TopZ().

Here is the call graph for this function:

◆ transferToCPP()

virtual void ICloneable::transferToCPP ( )
inlinevirtualinherited

Used for Python overriding of clone (see swig/tweaks.py)

Definition at line 32 of file ICloneable.h.

◆ volume()

double IFormFactor::volume ( ) const
virtualinherited

Reimplemented in IFormFactorPrism, IFormFactorPolyhedron, and Box.

Definition at line 83 of file IFormFactor.cpp.

84 {
85  return std::abs(formfactor_at_bottom(C3()));
86 }

References IFormFactor::formfactor_at_bottom().

Referenced by Compute::Slicing::createParticleInSlice().

Here is the call graph for this function:

Member Data Documentation

◆ m_height

const double& HemiEllipsoid::m_height
private

Definition at line 55 of file HemiEllipsoid.h.

Referenced by HemiEllipsoid(), clone(), formfactor_at_bottom(), and height().

◆ m_P

◆ m_radius_x

const double& HemiEllipsoid::m_radius_x
private

Definition at line 53 of file HemiEllipsoid.h.

Referenced by HemiEllipsoid(), clone(), formfactor_at_bottom(), radialExtension(), and radiusX().

◆ m_radius_y

const double& HemiEllipsoid::m_radius_y
private

Definition at line 54 of file HemiEllipsoid.h.

Referenced by clone(), formfactor_at_bottom(), radialExtension(), and radiusY().

◆ m_shape3D


The documentation for this class was generated from the following files: