BornAgain  1.19.0
Simulate and fit neutron and x-ray scattering at grazing incidence
ComputeDWBA Class Reference

Provides scalar DWBA computation for given IFormFactor. More...

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

Public Member Functions

 ComputeDWBA (const IFormFactor &ff)
 
 ~ComputeDWBA () override
 
virtual double bottomZ (const IRotation &rotation) const
 
ComputeDWBAclone () const override
 
complex_t evaluate (const WavevectorInfo &wavevectors) const override
 Returns the coherent sum of the four DWBA terms for scalar scattering. More...
 
virtual Eigen::Matrix2cd evaluatePol (const WavevectorInfo &wavevectors) const
 Returns scattering amplitude for matrix interactions. More...
 
virtual double radialExtension () const
 
virtual void setAmbientMaterial (const Material &material)
 
void setSpecularInfo (std::unique_ptr< const ILayerRTCoefficients > p_in_coeffs, std::unique_ptr< const ILayerRTCoefficients > p_out_coeffs) override
 Sets reflection/transmission info. More...
 
virtual double topZ (const IRotation &rotation) const
 
virtual double volume () const
 

Protected Attributes

std::unique_ptr< IFormFactorm_ff
 

Private Attributes

std::unique_ptr< const ILayerRTCoefficientsm_in_coeffs
 
std::unique_ptr< const ILayerRTCoefficientsm_out_coeffs
 

Friends

class TestPolarizedDWBATerms
 

Detailed Description

Provides scalar DWBA computation for given IFormFactor.

Definition at line 32 of file ComputeDWBA.h.

Constructor & Destructor Documentation

◆ ComputeDWBA()

ComputeDWBA::ComputeDWBA ( const IFormFactor ff)

Definition at line 20 of file ComputeDWBA.cpp.

20 : IComputeFF(ff) {}
IComputeFF()=delete

Referenced by clone().

◆ ~ComputeDWBA()

ComputeDWBA::~ComputeDWBA ( )
overridedefault

Member Function Documentation

◆ bottomZ()

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

Definition at line 39 of file IComputeFF.cpp.

40 {
41  return m_ff->bottomZ(rotation);
42 }
std::unique_ptr< IFormFactor > m_ff
Definition: IComputeFF.h:64

References IComputeFF::m_ff.

◆ clone()

ComputeDWBA * ComputeDWBA::clone ( ) const
overridevirtual

Implements IComputeFF.

Definition at line 24 of file ComputeDWBA.cpp.

25 {
26  ComputeDWBA* result = new ComputeDWBA(*m_ff);
27  std::unique_ptr<const ILayerRTCoefficients> p_in_coefs =
28  m_in_coeffs ? std::unique_ptr<const ILayerRTCoefficients>(m_in_coeffs->clone()) : nullptr;
29  std::unique_ptr<const ILayerRTCoefficients> p_out_coefs =
30  m_out_coeffs ? std::unique_ptr<const ILayerRTCoefficients>(m_out_coeffs->clone()) : nullptr;
31  result->setSpecularInfo(std::move(p_in_coefs), std::move(p_out_coefs));
32  return result;
33 }
Provides scalar DWBA computation for given IFormFactor.
Definition: ComputeDWBA.h:32
void setSpecularInfo(std::unique_ptr< const ILayerRTCoefficients > p_in_coeffs, std::unique_ptr< const ILayerRTCoefficients > p_out_coeffs) override
Sets reflection/transmission info.
Definition: ComputeDWBA.cpp:72
std::unique_ptr< const ILayerRTCoefficients > m_out_coeffs
Definition: ComputeDWBA.h:49
std::unique_ptr< const ILayerRTCoefficients > m_in_coeffs
Definition: ComputeDWBA.h:48
ComputeDWBA(const IFormFactor &ff)
Definition: ComputeDWBA.cpp:20

References ComputeDWBA(), IComputeFF::m_ff, m_in_coeffs, m_out_coeffs, and setSpecularInfo().

Here is the call graph for this function:

◆ evaluate()

complex_t ComputeDWBA::evaluate ( const WavevectorInfo wavevectors) const
overridevirtual

Returns the coherent sum of the four DWBA terms for scalar scattering.

Implements IComputeFF.

Definition at line 35 of file ComputeDWBA.cpp.

36 {
37  // Retrieve the two different incoming wavevectors in the layer
38  cvector_t k_i_T = wavevectors.getKi();
39  k_i_T.setZ(-m_in_coeffs->getScalarKz());
40  cvector_t k_i_R = k_i_T;
41  k_i_R.setZ(-k_i_T.z());
42 
43  // Retrieve the two different outgoing wavevector bins in the layer
44  cvector_t k_f_T = wavevectors.getKf();
45  k_f_T.setZ(m_out_coeffs->getScalarKz());
46  cvector_t k_f_R = k_f_T;
47  k_f_R.setZ(-k_f_T.z());
48 
49  // Construct the four different scattering contributions wavevector infos
50  double wavelength = wavevectors.wavelength();
51  WavevectorInfo k_TT(k_i_T, k_f_T, wavelength);
52  WavevectorInfo k_RT(k_i_R, k_f_T, wavelength);
53  WavevectorInfo k_TR(k_i_T, k_f_R, wavelength);
54  WavevectorInfo k_RR(k_i_R, k_f_R, wavelength);
55 
56  // Get the four R,T coefficients
57  complex_t T_in = m_in_coeffs->getScalarT();
58  complex_t R_in = m_in_coeffs->getScalarR();
59  complex_t T_out = m_out_coeffs->getScalarT();
60  complex_t R_out = m_out_coeffs->getScalarR();
61 
62  // The four different scattering contributions; S stands for scattering
63  // off the particle, R for reflection off the layer interface
64  complex_t term_S = T_in * m_ff->evaluate(k_TT) * T_out;
65  complex_t term_RS = R_in * m_ff->evaluate(k_RT) * T_out;
66  complex_t term_SR = T_in * m_ff->evaluate(k_TR) * R_out;
67  complex_t term_RSR = R_in * m_ff->evaluate(k_RR) * R_out;
68 
69  return term_S + term_RS + term_SR + term_RSR;
70 }
std::complex< double > complex_t
Definition: Complex.h:20
T z() const
Returns z-component in cartesian coordinate system.
Definition: BasicVector3D.h:67
void setZ(const T &a)
Sets z-component in cartesian coordinate system.
Definition: BasicVector3D.h:74
Holds all wavevector information relevant for calculating form factors.
cvector_t getKf() const
cvector_t getKi() const
double wavelength() const

References WavevectorInfo::getKf(), WavevectorInfo::getKi(), IComputeFF::m_ff, m_in_coeffs, m_out_coeffs, BasicVector3D< T >::setZ(), WavevectorInfo::wavelength(), and BasicVector3D< T >::z().

Here is the call graph for this function:

◆ evaluatePol()

Eigen::Matrix2cd IComputeFF::evaluatePol ( const WavevectorInfo wavevectors) const
virtualinherited

Returns scattering amplitude for matrix interactions.

Reimplemented in ComputeDWBAPol, and ComputeBAPol.

Definition at line 49 of file IComputeFF.cpp.

50 {
51  throw std::runtime_error("Bug: impossible call to FFCompute::evaluatePol");
52 }

◆ radialExtension()

double IComputeFF::radialExtension ( ) const
virtualinherited

Definition at line 34 of file IComputeFF.cpp.

35 {
36  return m_ff->radialExtension();
37 }

References IComputeFF::m_ff.

◆ setAmbientMaterial()

void IComputeFF::setAmbientMaterial ( const Material material)
virtualinherited

Definition at line 24 of file IComputeFF.cpp.

25 {
26  m_ff->setAmbientMaterial(material);
27 }

References IComputeFF::m_ff.

◆ setSpecularInfo()

void ComputeDWBA::setSpecularInfo ( std::unique_ptr< const ILayerRTCoefficients ,
std::unique_ptr< const ILayerRTCoefficients  
)
overridevirtual

Sets reflection/transmission info.

Reimplemented from IComputeFF.

Definition at line 72 of file ComputeDWBA.cpp.

74 {
75  m_in_coeffs = std::move(p_in_coeffs);
76  m_out_coeffs = std::move(p_out_coeffs);
77 }

References m_in_coeffs, and m_out_coeffs.

Referenced by clone().

◆ topZ()

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

Definition at line 44 of file IComputeFF.cpp.

45 {
46  return m_ff->topZ(rotation);
47 }

References IComputeFF::m_ff.

◆ volume()

double IComputeFF::volume ( ) const
virtualinherited

Definition at line 29 of file IComputeFF.cpp.

30 {
31  return m_ff->volume();
32 }

References IComputeFF::m_ff.

Friends And Related Function Documentation

◆ TestPolarizedDWBATerms

friend class TestPolarizedDWBATerms
friend

Definition at line 45 of file ComputeDWBA.h.

Member Data Documentation

◆ m_ff

◆ m_in_coeffs

std::unique_ptr<const ILayerRTCoefficients> ComputeDWBA::m_in_coeffs
private

Definition at line 48 of file ComputeDWBA.h.

Referenced by clone(), evaluate(), and setSpecularInfo().

◆ m_out_coeffs

std::unique_ptr<const ILayerRTCoefficients> ComputeDWBA::m_out_coeffs
private

Definition at line 49 of file ComputeDWBA.h.

Referenced by clone(), evaluate(), and setSpecularInfo().


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