From d14bb20e01ac5903bf6c5dd9c3df56b5092cbb7f Mon Sep 17 00:00:00 2001 From: pospelov <pospelov@fz-juelich.de> Date: Thu, 19 Apr 2012 17:45:06 +0200 Subject: [PATCH] On the way to add roughnesses to layers --- App/App.pro | 6 ++- App/src/TestFresnelCoeff.cpp | 8 +-- App/src/main.cpp | 12 ++++- Core/Core.pro | 9 +++- Core/inc/Exceptions.h | 6 +++ Core/inc/HomogeneousMaterial.h | 29 +++++++++- Core/inc/IMaterial.h | 35 +++++++++++++ Core/inc/ISample.h | 1 + Core/inc/Layer.h | 50 ++++++++++++++---- Core/inc/LayerRoughness.h | 51 ++++++++++++++++-- Core/inc/MultiLayer.h | 74 +++++++++++++++++++++----- Core/src/HomogeneousMaterial.cpp | 35 ++++++++++++- Core/src/ISample.cpp | 1 + Core/src/Layer.cpp | 28 +++++++++- Core/src/LayerRoughness.cpp | 31 ++++++++++- Core/src/MultiLayer.cpp | 90 ++++++++++++++++++++++++++------ Core/src/OpticalFresnel.cpp | 4 +- 17 files changed, 409 insertions(+), 61 deletions(-) diff --git a/App/App.pro b/App/App.pro index 8b4dfd4b931..419dcc6e08e 100644 --- a/App/App.pro +++ b/App/App.pro @@ -7,13 +7,15 @@ QT -= core gui SOURCES += \ src/main.cpp \ src/TestFresnelCoeff.cpp \ - src/DrawHelper.cpp + src/DrawHelper.cpp \ + src/TestDiffuseScattering.cpp HEADERS += \ inc/DrawHelper.h \ inc/TestFresnelCoeff.h \ inc/App.h \ - inc/AppLinkDef.h + inc/AppLinkDef.h \ + inc/TestDiffuseScattering.h INCLUDEPATH += ./inc diff --git a/App/src/TestFresnelCoeff.cpp b/App/src/TestFresnelCoeff.cpp index 4044ddbcec6..e17657bf03d 100644 --- a/App/src/TestFresnelCoeff.cpp +++ b/App/src/TestFresnelCoeff.cpp @@ -48,14 +48,14 @@ void TestFresnelCoeff::execute() lSubstrate.setMaterial(&mSubstrate, 0); // adding layers - mySample.add(&lAmbience); + mySample.add(lAmbience); const unsigned nrepetitions = 4; for(unsigned i=0; i<nrepetitions; ++i) { - mySample.add(&lAg1); - mySample.add(&lCr1); + mySample.add(lAg1); + mySample.add(lCr1); } - mySample.add(&lSubstrate); + mySample.add(lSubstrate); // ------------------- // calculation of Fresnel coefficients for multi-layers system diff --git a/App/src/main.cpp b/App/src/main.cpp index cff32ae7626..2d19131f80c 100644 --- a/App/src/main.cpp +++ b/App/src/main.cpp @@ -1,13 +1,14 @@ #include <iostream> #include <string> #include "TestFresnelCoeff.h" +#include "TestDiffuseScattering.h" #include "TApplication.h" -using namespace std; + int main(int argc, char **argv) { - cout << "Hello Brave New World!" << endl; + std::cout << "Hello Brave New World!" << std::endl; TApplication theApp("theApp",0,0); if(argc>=2) { @@ -18,6 +19,13 @@ int main(int argc, char **argv) TestFresnelCoeff test; test.execute(); } + + // user algorithm to test reflection/refraction coefficients for multilayer system + if(spar.find("diffuse") != std::string::npos) { + TestDiffuseScattering test; + test.execute(); + } + } if(gApplication) { diff --git a/Core/Core.pro b/Core/Core.pro index b8ba5f082ce..15db156a1be 100644 --- a/Core/Core.pro +++ b/Core/Core.pro @@ -20,7 +20,10 @@ SOURCES += \ src/CalculatorOptical.cpp \ src/ISimulation.cpp \ src/OutputData.cpp \ - src/OpticalFresnel.cpp + src/OpticalFresnel.cpp \ + src/IRoughness.cpp \ + src/LayerInterface.cpp \ + src/IMaterial.cpp HEADERS += \ inc/ISample.h \ @@ -37,7 +40,9 @@ HEADERS += \ inc/ISimulation.h \ inc/OutputData.h \ inc/NamedVector.h \ - inc/OpticalFresnel.h + inc/OpticalFresnel.h \ + inc/IRoughness.h \ + inc/LayerInterface.h INCLUDEPATH += ./inc diff --git a/Core/inc/Exceptions.h b/Core/inc/Exceptions.h index 6c444bc4a73..9eb74fd7283 100644 --- a/Core/inc/Exceptions.h +++ b/Core/inc/Exceptions.h @@ -27,4 +27,10 @@ public: ClassInitializationException(const std::string& message); }; +class SelfReferenceException : public std::logic_error +{ +public: + SelfReferenceException(const std::string &message); +}; + #endif // EXCEPTIONS_H diff --git a/Core/inc/HomogeneousMaterial.h b/Core/inc/HomogeneousMaterial.h index 8ef3549bbac..83165710986 100644 --- a/Core/inc/HomogeneousMaterial.h +++ b/Core/inc/HomogeneousMaterial.h @@ -1,20 +1,45 @@ #ifndef HOMOGENEOUSMATERIAL_H #define HOMOGENEOUSMATERIAL_H +// ******************************************************************** +// * The BornAgain project * +// * Simulation of neutron and x-ray scattering at grazing incidence * +// * * +// * LICENSE AND DISCLAIMER * +// * Lorem ipsum dolor sit amet, consectetur adipiscing elit. Mauris * +// * eget quam orci. Quisque porta varius dui, quis posuere nibh * +// * mollis quis. Mauris commodo rhoncus porttitor. * +// ******************************************************************** +//! @file HomogeneousMaterial.h +//! @brief Definition of HomogeneousMaterial class +//! @author James Bond <j.bond@fz-juelich.de>, Chuck Norris <c.norris@fz-juelich.de> +//! @date 01.04.2012 #include "IMaterial.h" #include "Types.h" + +//- ------------------------------------------------------------------- +/// @class HomogeneousMaterial +/// @brief Definition of homogeneous material with refraction index +//- ------------------------------------------------------------------- class HomogeneousMaterial : public IMaterial { public: + HomogeneousMaterial(); HomogeneousMaterial(complex_t refractive_index); + HomogeneousMaterial(const std::string &name, complex_t refractive_index); + HomogeneousMaterial(const HomogeneousMaterial &other); + HomogeneousMaterial &operator=(const HomogeneousMaterial &other); virtual ~HomogeneousMaterial() {} - complex_t getRefractiveIndex() { return m_refractive_index; } + /// return refractive index of the material + complex_t getRefractiveIndex() const { return m_refractive_index; } + /// set refractive index of he material + void setRefractiveIndex(complex_t refractive_index) { m_refractive_index = refractive_index; } private: - complex_t m_refractive_index; + complex_t m_refractive_index; ///< complex index of refraction }; #endif // HOMOGENEOUSMATERIAL_H diff --git a/Core/inc/IMaterial.h b/Core/inc/IMaterial.h index 1c40c864cc5..45e970efac6 100644 --- a/Core/inc/IMaterial.h +++ b/Core/inc/IMaterial.h @@ -1,10 +1,45 @@ #ifndef IMATERIAL_H #define IMATERIAL_H +// ******************************************************************** +// * The BornAgain project * +// * Simulation of neutron and x-ray scattering at grazing incidence * +// * * +// * LICENSE AND DISCLAIMER * +// * Lorem ipsum dolor sit amet, consectetur adipiscing elit. Mauris * +// * eget quam orci. Quisque porta varius dui, quis posuere nibh * +// * mollis quis. Mauris commodo rhoncus porttitor. * +// ******************************************************************** +//! @file IMaterial.h +//! @brief Definition of IMaterial class +//! @author James Bond <j.bond@fz-juelich.de>, Chuck Norris <c.norris@fz-juelich.de> +//! @date 01.04.2012 +#include <string> + + +//- ------------------------------------------------------------------- +/// @class IMaterial +/// @brief Material definition +//- ------------------------------------------------------------------- class IMaterial { public: + IMaterial() {} + IMaterial(const std::string &name) : m_name(name) {} + /// copy constructor + IMaterial(const IMaterial &other); + /// assignment operator + IMaterial &operator=(const IMaterial &other); virtual ~IMaterial() {} + + /// set name of the material + void setName(const std::string &name) { m_name = name; } + + /// return name of the material + std::string getName() const { return m_name; } + +private: + std::string m_name; }; #endif // IMATERIAL_H diff --git a/Core/inc/ISample.h b/Core/inc/ISample.h index a2f0ed88816..7bdc412fdfb 100644 --- a/Core/inc/ISample.h +++ b/Core/inc/ISample.h @@ -6,6 +6,7 @@ class ISample { public: + ISample(){} virtual ~ISample() {} virtual void add(ISample* p_child); diff --git a/Core/inc/Layer.h b/Core/inc/Layer.h index a35a7678127..a69c58db50b 100644 --- a/Core/inc/Layer.h +++ b/Core/inc/Layer.h @@ -1,30 +1,62 @@ #ifndef LAYER_H #define LAYER_H +// ******************************************************************** +// * The BornAgain project * +// * Simulation of neutron and x-ray scattering at grazing incidence * +// * * +// * LICENSE AND DISCLAIMER * +// * Lorem ipsum dolor sit amet, consectetur adipiscing elit. Mauris * +// * eget quam orci. Quisque porta varius dui, quis posuere nibh * +// * mollis quis. Mauris commodo rhoncus porttitor. * +// ******************************************************************** +//! @file Layer.h +//! @brief Definition of Layer class +//! @author James Bond <j.bond@fz-juelich.de>, Chuck Norris <c.norris@fz-juelich.de> +//! @date 01.04.2012 #include "ISample.h" #include "IMaterial.h" #include "Types.h" #include "HomogeneousMaterial.h" + +//- ------------------------------------------------------------------- +/// @class Layer +/// @brief Definition of Layer with thickness and pointer to the material +//- ------------------------------------------------------------------- class Layer : public ISample { public: Layer(); + Layer(const Layer &other); + Layer &operator=(const Layer &other); virtual ~Layer(); + /// set layer thickness in _angstrom_ virtual void setThickness(double thickness); + + /// return layer thickness in _angstrom_ virtual double getThickness() const { return m_thickness; } - virtual void setMaterial(IMaterial* p_material) { mp_bulk_material = p_material; } - virtual void setMaterial(IMaterial* p_material, double thickness); - virtual IMaterial* getMaterial() { return mp_bulk_material; } - virtual complex_t getRefractiveIndex() const { return (dynamic_cast<HomogeneousMaterial *>(mp_bulk_material))->getRefractiveIndex(); } -private: - IMaterial* mp_bulk_material; - double m_thickness; -// LayerRoughness* mp_top_roughness; -// LayerRoughness* mp_bottom_roughness; + /// set material to the layer + /// @param p_material pointer to the material + virtual void setMaterial(const IMaterial* p_material) { p_material ? mp_material = p_material : throw NullPointerException("The material doesn't exist"); } + + /// set material of given thickness to the layer + /// @param p_material pointer to the material of layer + /// @param thickness thickness of the material in angstrom + virtual void setMaterial(const IMaterial* p_material, double thickness); + /// return layer's material + virtual const IMaterial* getMaterial() const { return mp_material; } + + /// return refractive index of the layer's material + virtual complex_t getRefractiveIndex() const { return (dynamic_cast<const HomogeneousMaterial *>(mp_material))->getRefractiveIndex(); } + +private: + const IMaterial* mp_material; ///< pointer to the material + double m_thickness; ///< layer thickness in _angstrom_ }; + #endif // LAYER_H diff --git a/Core/inc/LayerRoughness.h b/Core/inc/LayerRoughness.h index 8a4a9462391..2e5828ac756 100644 --- a/Core/inc/LayerRoughness.h +++ b/Core/inc/LayerRoughness.h @@ -1,15 +1,58 @@ #ifndef LAYERROUGHNESS_H #define LAYERROUGHNESS_H +// ******************************************************************** +// * The BornAgain project * +// * Simulation of neutron and x-ray scattering at grazing incidence * +// * * +// * LICENSE AND DISCLAIMER * +// * Lorem ipsum dolor sit amet, consectetur adipiscing elit. Mauris * +// * eget quam orci. Quisque porta varius dui, quis posuere nibh * +// * mollis quis. Mauris commodo rhoncus porttitor. * +// ******************************************************************** +//! @file LayerRoughness.h +//! @brief Definition of LayerRoughness class +//! @author James Bond <j.bond@fz-juelich.de>, Chuck Norris <c.norris@fz-juelich.de> +//! @date 01.04.2012 -#include "Layer.h" +#include "Types.h" +#include "IRoughness.h" -class LayerRoughness +//- ------------------------------------------------------------------- +/// @class LayerRoughness +/// @brief Roughness of interface between two layers. +/// +/// Based on the article +/// D.K.G. de Boer, Physical review B, Volume 51, Number 8, 15 February 1995 +/// "X-ray reflection and transmission by rough surfaces" +//- ------------------------------------------------------------------- +class LayerRoughness : public IRoughness { public: - static LayerRoughness* createSmoothLayerInterface(Layer* p_top_layer, Layer* p_bottom_layer); + LayerRoughness(); + LayerRoughness(double sigma, double hurstParameter, double latteralCorrLength); + + /// return power spectral density of the surface roughness + double getPowerSpectralDensity(const kvector_t &kvec); + + /// set rms value of roughness + void setSigma(double sigma) { m_sigma = sigma; } + /// return rms value of roughness + double getSigma() const { return m_sigma; } + + /// Set hurst parameter. It describes how jagged the surface is. + inline void setHurstParameter(double hurstParameter) { m_hurstParameter = hurstParameter; } + /// return hurst parameter + inline double getHurstParameter() const { return m_hurstParameter; } + + /// set lateral correlation length + inline void setLatteralCorrLength(double latteralCorrLength) { m_latteralCorrLength = latteralCorrLength; } + /// return latteral correlation length + inline double getLatteralCorrLength() const { return m_latteralCorrLength; } protected: - LayerRoughness(); + double m_sigma; ///< rms of roughness + double m_hurstParameter; ///< Hurst parameter which describes how jagged the interface, 0<H<=1 + double m_latteralCorrLength; ///< latteral correlation length of the roughness }; #endif // LAYERROUGHNESS_H diff --git a/Core/inc/MultiLayer.h b/Core/inc/MultiLayer.h index a8bc99f72cd..833757e0e65 100644 --- a/Core/inc/MultiLayer.h +++ b/Core/inc/MultiLayer.h @@ -1,37 +1,85 @@ #ifndef MULTILAYER_H #define MULTILAYER_H +// ******************************************************************** +// * The BornAgain project * +// * Simulation of neutron and x-ray scattering at grazing incidence * +// * * +// * LICENSE AND DISCLAIMER * +// * Lorem ipsum dolor sit amet, consectetur adipiscing elit. Mauris * +// * eget quam orci. Quisque porta varius dui, quis posuere nibh * +// * mollis quis. Mauris commodo rhoncus porttitor. * +// ******************************************************************** +//! @file MultiLayer.h +//! @brief Definition of MultiLayer class +//! @author James Bond <j.bond@fz-juelich.de>, Chuck Norris <c.norris@fz-juelich.de> +//! @date 01.04.2012 #include <vector> #include "Layer.h" +#include "LayerInterface.h" #include "LayerRoughness.h" + +//- ------------------------------------------------------------------- +/// @class MultiLayer +/// @brief Stack of layers one on top of the other +/// +/// Example of system of 4 layers: +/// +/// ambience layer #0 z=getLayerBottomZ(0)=0.0 +/// --------- interface #0 +/// Fe, 20A layer #1 z=getLayerBottomZ(1)=-20.0 +/// --------- interface #1 +/// Cr, 40A layer #2 z=getLayerBottomZ(2)=-60.0 +/// --------- interface #2 +/// substrate layer #3 z=getLayerBottomZ(3)=-60.0 +/// +//- ------------------------------------------------------------------- class MultiLayer : public ISample { public: MultiLayer(); ~MultiLayer(); - // adds layer with top roughness - void addLayerWithTopRoughness(Layer* p_layer, LayerRoughness* p_roughness); + /// return number of layers in multilayer + inline size_t getNumberOfLayers() const { return m_layers.size(); } + + /// add layer with top roughness + void addLayerWithTopRoughness(const Layer &layer, const LayerRoughness &roughness); + + /// add object to multilayer, overrides from ISample + void add(const ISample &p_child); + + /// return layer with given index + const Layer *getLayer(size_t i_layer) const { return &m_layers[ check_layer_index(i_layer) ]; } - // returns number of layers in multilayer - size_t getNumberOfLayers() const { return m_layers.size(); } + /// return z-coordinate of the layer's bottom + double getLayerBottomZ(size_t i_layer) const { return m_layers_z[ check_layer_index(i_layer) ]; } - // Overrides from ISample: - void add(ISample* p_child); + /// return bottom interface of layer + //const LayerInterface *getLayerBottomInterface(size_t i_layer); - // returns layer with given index - const Layer *getLayer(size_t i_layer) const { return m_layers[ check_layer_index(i_layer) ]; } + /// destruct allocated objects + void clear(); - // returns layer z-coordinate - double getLayerZ(size_t i_layer) const { return m_layers_z[ check_layer_index(i_layer) ]; } + /// return clone of multilayer with clones of all layers and new interfaces between layers + MultiLayer *clone(); private: + /// hiding copy constructor + MultiLayer(const MultiLayer &other); + + /// hiding assignment operator + MultiLayer &operator=(const MultiLayer &other); + + /// check index of layer w.r.t. vector length inline size_t check_layer_index(size_t i_layer) const { return i_layer < m_layers.size() ? i_layer : throw OutOfBoundsException("Layer index is out of bounds"); } - std::vector<Layer*> m_layers; - std::vector<LayerRoughness*> m_roughnesses; - std::vector<double > m_layers_z; + + std::vector<Layer > m_layers; ///< stack of layers [nlayers] + std::vector<double > m_layers_z; ///< coordinate of layer's bottoms [nlayers] + std::vector<LayerInterface *> m_interfaces; ///< stack of layer interfaces [nlayers-1] }; #endif // MULTILAYER_H + diff --git a/Core/src/HomogeneousMaterial.cpp b/Core/src/HomogeneousMaterial.cpp index 2bf7e38e613..e45832c838a 100644 --- a/Core/src/HomogeneousMaterial.cpp +++ b/Core/src/HomogeneousMaterial.cpp @@ -1,6 +1,39 @@ #include "HomogeneousMaterial.h" + +HomogeneousMaterial::HomogeneousMaterial() +{ + +} + + HomogeneousMaterial::HomogeneousMaterial(complex_t refractive_index) - : m_refractive_index(refractive_index) + : IMaterial("noname"), m_refractive_index(refractive_index) { + } + + +HomogeneousMaterial::HomogeneousMaterial(const std::string &name, complex_t refractive_index) + : IMaterial(name), m_refractive_index(refractive_index) +{ + +} + + +HomogeneousMaterial::HomogeneousMaterial(const HomogeneousMaterial &other) : IMaterial(other) +{ + m_refractive_index = m_refractive_index; +} + + +HomogeneousMaterial &HomogeneousMaterial::operator=(const HomogeneousMaterial &other) +{ + if(this != &other) + { + IMaterial::operator=(other); + m_refractive_index = other.m_refractive_index; + } + return *this; +} + diff --git a/Core/src/ISample.cpp b/Core/src/ISample.cpp index 30b7c3f1115..ec8946be2af 100644 --- a/Core/src/ISample.cpp +++ b/Core/src/ISample.cpp @@ -1,5 +1,6 @@ #include "ISample.h" + void ISample::add(ISample* p_child) { throw NotImplementedException("This sample class is not allowed to have subsamples."); diff --git a/Core/src/Layer.cpp b/Core/src/Layer.cpp index 9fa21670fa2..67f4bfe2513 100644 --- a/Core/src/Layer.cpp +++ b/Core/src/Layer.cpp @@ -1,14 +1,38 @@ #include "Layer.h" #include <stdexcept> -Layer::Layer() : mp_bulk_material(0), m_thickness(0) +/* ************************************************************************* */ +Layer::Layer() : mp_material(0), m_thickness(0) { } + +Layer::Layer(const Layer &other) : ISample(other) +{ + mp_material = other.mp_material; + m_thickness = other.m_thickness; +} + + +Layer &Layer::operator=(const Layer &other) +{ + if( this != &other) + { + ISample::operator=(other); + mp_material = other.mp_material; + m_thickness = other.m_thickness; + } + return *this; +} + + Layer::~Layer() { + } + +/* ************************************************************************* */ void Layer::setThickness(double thickness) { if (thickness>=0.0) @@ -20,7 +44,7 @@ void Layer::setThickness(double thickness) } -void Layer::setMaterial(IMaterial* p_material, double thickness) +void Layer::setMaterial(const IMaterial* p_material, double thickness) { setMaterial(p_material); setThickness(thickness); diff --git a/Core/src/LayerRoughness.cpp b/Core/src/LayerRoughness.cpp index f0afaf8031f..f15b74a9f5b 100644 --- a/Core/src/LayerRoughness.cpp +++ b/Core/src/LayerRoughness.cpp @@ -1,6 +1,33 @@ #include "LayerRoughness.h" +#include <math.h> -LayerRoughness *LayerRoughness::createSmoothLayerInterface(Layer *p_top_layer, Layer *p_bottom_layer) + +LayerRoughness::LayerRoughness() : m_sigma(0), m_hurstParameter(0), m_latteralCorrLength(0) +{ + +} + + +LayerRoughness::LayerRoughness(double sigma, double hurstParameter, double latteralCorrLength) : + m_sigma(sigma), m_hurstParameter(hurstParameter), m_latteralCorrLength(latteralCorrLength) { - return 0; + +} + + +/* ************************************************************************* */ +//! Power spectral density of the surface roughness is a result of two-dimensional +//! Fourier transform of the correlation function of the roughness profile. +//! +//! Based on the article +//! D.K.G. de Boer, Physical review B, Volume 51, Number 8, 15 February 1995 +//! "X-ray reflection and transmission by rough surfaces" +/* ************************************************************************* */ +double LayerRoughness::getPowerSpectralDensity(const kvector_t &kvec) +{ + double H = m_hurstParameter; + double xsi = m_latteralCorrLength; + double Qpar2 = kvec.x()*kvec.x() + kvec.y()*kvec.y(); + return 4.0*M_PI*H * m_sigma*m_sigma * xsi*xsi * std::pow( (1.0 + Qpar2*xsi*xsi), (-1-H) ); } + diff --git a/Core/src/MultiLayer.cpp b/Core/src/MultiLayer.cpp index 5eb6c36d092..ba94f414f6a 100644 --- a/Core/src/MultiLayer.cpp +++ b/Core/src/MultiLayer.cpp @@ -7,38 +7,96 @@ MultiLayer::MultiLayer() { } + +//MultiLayer::MultiLayer(const MultiLayer &other) : ISample(other) +//{ +// m_layers = other.m_layers; +// m_layers_z = other.m_layers_z; +// m_interfaces = other.m_interfaces; +//} + + +//MultiLayer &MultiLayer::operator=(const MultiLayer &other) +//{ +// if( this != &other) +// { + +// clean(); +// m_layers = other.m_layers; +// m_layers_z = other.m_layers_z; +// //m_interfaces = other.m_interfaces; +// } +// return *this; +//} + + MultiLayer::~MultiLayer() { + } -void MultiLayer::addLayerWithTopRoughness(Layer *p_layer, LayerRoughness *p_roughness) + +/* ************************************************************************* */ +// clear MultiLayer contents including interfaces +/* ************************************************************************* */ +void MultiLayer::clear() { - if (!p_layer) + m_layers.clear(); + m_layers_z.clear(); + for(size_t i=0; i<m_interfaces.size(); i++) { - throw NullPointerException("The layer to add does not exist."); + if( m_interfaces[i] ) { + delete m_interfaces[i]; + m_interfaces[i] = 0; + } + m_interfaces.clear(); } +} + + +/* ************************************************************************* */ +// clear MultiLayer contents including interfaces +/* ************************************************************************* */ +MultiLayer *MultiLayer::clone() +{ +// MultiLayer *newML = new MultiLayer; +// newML->m_layers = m_layers; +// newML->m_layers_z = m_layers_z; +// for(size_t i=0; i<m_interfaces.size(); i++) { + +// } +// newMl->m_interfaces.resize(m_layers.) + return 0; +} + + + +/* ************************************************************************* */ +// clear MultiLayer contents including interfaces +/* ************************************************************************* */ +void MultiLayer::addLayerWithTopRoughness(const Layer &layer, const LayerRoughness &roughness) +{ if (getNumberOfLayers()) { - Layer* p_last_layer = dynamic_cast<Layer*>(m_layers[getNumberOfLayers()-1]); - LayerRoughness* p_top_roughness = p_roughness; - if (!p_top_roughness) - { - p_top_roughness = LayerRoughness::createSmoothLayerInterface(p_last_layer, p_layer); - } - m_layers.push_back(p_layer); - m_roughnesses.push_back(p_top_roughness); - m_layers_z.push_back(m_layers_z.back() - p_layer->getThickness() ); + //const Layer *p_last_layer = getLayer( getNumberOfLayers()-1 ); + const Layer *p_previous_layer = &m_layers.back(); + m_layers.push_back(layer); + const Layer *p_new_layer = &m_layers.back(); + LayerInterface *interface = LayerInterface::createRoughInterface( p_previous_layer, p_new_layer, roughness); + m_interfaces.push_back(interface); + m_layers_z.push_back(m_layers_z.back() - layer.getThickness() ); return; } - m_layers.push_back(p_layer); + m_layers.push_back(layer); m_layers_z.push_back(0.0); } -void MultiLayer::add(ISample* p_child) +void MultiLayer::add(const ISample &p_child) { - Layer* p_layer = dynamic_cast<Layer*>(p_child); - addLayerWithTopRoughness(p_layer, 0); + LayerRoughness defRough; + const Layer &layer = dynamic_cast<const Layer&>(p_child); + addLayerWithTopRoughness(layer, defRough); } diff --git a/Core/src/OpticalFresnel.cpp b/Core/src/OpticalFresnel.cpp index 27ad1f1973f..449e027b115 100644 --- a/Core/src/OpticalFresnel.cpp +++ b/Core/src/OpticalFresnel.cpp @@ -29,7 +29,7 @@ int OpticalFresnel::execute(const MultiLayer &sample, const kvector_t &kvec, Mul // ratio of amplitudes of outgoing and incoming waves for(int i=fC.size()-2; i>=0; --i) { - double z = sample.getLayerZ(i); + double z = sample.getLayerBottomZ(i); fC[i].X = std::exp(complex_t(0,-2)*fC[i].kz*z) * (fC[i].r + fC[i+1].X*std::exp(complex_t(0,2)*fC[i+1].kz*z) ) / ( complex_t(1,0)+fC[i].r*fC[i+1].X*std::exp(complex_t(0,2)*fC[i+1].kz*z) ); @@ -38,7 +38,7 @@ int OpticalFresnel::execute(const MultiLayer &sample, const kvector_t &kvec, Mul fC[0].R = fC[0].X; fC[0].T = 1; for(size_t i=0; i<fC.size()-1; ++i) { - double z = sample.getLayerZ(i); + double z = sample.getLayerBottomZ(i); fC[i+1].R = complex_t(1,0)/fC[i].tb * (fC[i].T*fC[i].rb*std::exp( complex_t(0,-1)*(fC[i+1].kz+fC[i].kz)*z ) + fC[i].R*std::exp(complex_t(0,-1)*(fC[i+1].kz-fC[i].kz)*z) ); -- GitLab