From 6a1697b9c721cc5567ed1a17605d1185a9aa2e70 Mon Sep 17 00:00:00 2001 From: Walter Van Herck <w.van.herck@fz-juelich.de> Date: Mon, 29 Oct 2012 17:36:26 +0100 Subject: [PATCH] Added parameter registration, xi integration flag usage, isgisaxs8 first sample initialization, debug symbols in release --- App/src/FunctionalTestFactory.cpp | 5 +- App/src/StandardSamples.cpp | 7 ++- App/src/TestIsGISAXS8.cpp | 2 + Core/Algorithms/inc/FTDistributions.h | 6 ++- Core/Algorithms/src/FTDistributions.cpp | 9 ++++ .../inc/InterferenceFunction2DParaCrystal.h | 4 +- .../src/InterferenceFunction2DParaCrystal.cpp | 46 ++++++++++++++++--- shared.pri | 2 +- 8 files changed, 70 insertions(+), 11 deletions(-) diff --git a/App/src/FunctionalTestFactory.cpp b/App/src/FunctionalTestFactory.cpp index 3f4c07c111a..2f28ae4ea1f 100644 --- a/App/src/FunctionalTestFactory.cpp +++ b/App/src/FunctionalTestFactory.cpp @@ -7,6 +7,7 @@ #include "TestIsGISAXS2.h" #include "TestIsGISAXS3.h" #include "TestIsGISAXS4.h" +#include "TestIsGISAXS8.h" #include "TestIsGISAXS9.h" #include "TestIsGISAXS10.h" #include "TestIsGISAXS11.h" @@ -44,7 +45,9 @@ FunctionalTestFactory::FunctionalTestFactory() : m_benchmark(0) registerItem("isgisaxs3", IFactoryCreateFunction<TestIsGISAXS3, IFunctionalTest>, "functional test: isgisaxs ex-3 (cylinder FF)"); registerItem("isgisaxs4", IFactoryCreateFunction<TestIsGISAXS4, IFunctionalTest>, - "functional test: isgisaxs ex-4 (paracrystal structure factors)"); + "functional test: isgisaxs ex-4 (paracrystal 1d structure factors)"); + registerItem("isgisaxs8", IFactoryCreateFunction<TestIsGISAXS8, IFunctionalTest>, + "functional test: isgisaxs ex-8 (paracrystal lattice structure factors)"); registerItem("isgisaxs9", IFactoryCreateFunction<TestIsGISAXS9, IFunctionalTest>, "functional test: isgisaxs ex-9 (rotated pyramid FF)"); registerItem("isgisaxs10", IFactoryCreateFunction<TestIsGISAXS10, IFunctionalTest>, diff --git a/App/src/StandardSamples.cpp b/App/src/StandardSamples.cpp index 2bf149fe5e8..77a595d74b9 100644 --- a/App/src/StandardSamples.cpp +++ b/App/src/StandardSamples.cpp @@ -456,7 +456,12 @@ ISample *StandardSamples::IsGISAXS8_2DDL_lattice() air_layer.setMaterial(p_air_material); Layer substrate_layer; substrate_layer.setMaterial(p_substrate_material); - IInterferenceFunction *p_interference_function = new InterferenceFunction1DParaCrystal(20.0*Units::nanometer,7*Units::nanometer, 1e3*Units::nanometer); + InterferenceFunction2DParaCrystal *p_interference_function = new InterferenceFunction2DParaCrystal(10.0*Units::nanometer, 10.0*Units::nanometer, M_PI/2.0, + 0.0, 1e7*Units::nanometer); + p_interference_function->setDomainSizes(20.0*Units::micrometer, 20.0*Units::micrometer); + FTDistribution2DCauchy pdf1(0.5*Units::nanometer, 2.0*Units::nanometer); + FTDistribution2DCauchy pdf2(2.0*Units::nanometer, 0.5*Units::nanometer); + p_interference_function->setProbabilityDistributions(pdf1, pdf2); ParticleDecoration particle_decoration( new Particle(n_particle, new FormFactorCylinder(5*Units::nanometer, 5*Units::nanometer))); particle_decoration.addInterferenceFunction(p_interference_function); LayerDecorator air_layer_decorator(air_layer, particle_decoration); diff --git a/App/src/TestIsGISAXS8.cpp b/App/src/TestIsGISAXS8.cpp index 0d0707664df..23c8861474b 100644 --- a/App/src/TestIsGISAXS8.cpp +++ b/App/src/TestIsGISAXS8.cpp @@ -27,6 +27,8 @@ void TestIsGISAXS8::execute() // 2DDL_lattice p_sample = dynamic_cast<MultiLayer *>(SampleFactory::instance().createItem("IsGISAXS8_2DDL_lattice")); +// ParameterPool *p_pool = p_sample->createParameterTree(); +// std::cout << (*p_pool); experiment.setSample(*p_sample); experiment.runSimulation(); IsGISAXSTools::writeOutputDataToFile(*experiment.getOutputData(), m_data_path+"this_2DDL_lattice.ima"); diff --git a/Core/Algorithms/inc/FTDistributions.h b/Core/Algorithms/inc/FTDistributions.h index 5dbc2a4cb91..d15a793cb41 100644 --- a/Core/Algorithms/inc/FTDistributions.h +++ b/Core/Algorithms/inc/FTDistributions.h @@ -14,6 +14,8 @@ //! @author Scientific Computing Group at FRM II //! @date Oct 26, 2012 +#include "IParameterized.h" +#include "ParameterPool.h" #include <cmath> class IFTDistribution1D @@ -27,7 +29,7 @@ protected: double m_omega; }; -class IFTDistribution2D +class IFTDistribution2D : public IParameterized { public: IFTDistribution2D(double omega_x, double omega_y) : m_omega_x(omega_x), m_omega_y(omega_y), @@ -55,6 +57,8 @@ public: virtual FTDistribution2DCauchy *clone() const; virtual double evaluate(double qx, double qy) const; +protected: + virtual void init_parameters(); }; #endif /* FTDISTRIBUTIONS_H_ */ diff --git a/Core/Algorithms/src/FTDistributions.cpp b/Core/Algorithms/src/FTDistributions.cpp index 33de6450830..3a2420455df 100644 --- a/Core/Algorithms/src/FTDistributions.cpp +++ b/Core/Algorithms/src/FTDistributions.cpp @@ -3,6 +3,8 @@ FTDistribution2DCauchy::FTDistribution2DCauchy(double omega_x, double omega_y) : IFTDistribution2D(omega_x, omega_y) { + setName("2DDistributionCauchy"); + init_parameters(); } FTDistribution2DCauchy* FTDistribution2DCauchy::clone() const @@ -15,3 +17,10 @@ double FTDistribution2DCauchy::evaluate(double qx, double qy) const double sum_sq = qx*qx*m_omega_x*m_omega_x + qy*qy*m_omega_y*m_omega_y; return std::pow(1.0 + sum_sq, -1.5); } + +void FTDistribution2DCauchy::init_parameters() +{ + getParameterPool()->clear(); + getParameterPool()->registerParameter("omega_x", &m_omega_x); + getParameterPool()->registerParameter("omega_y", &m_omega_y); +} diff --git a/Core/Samples/inc/InterferenceFunction2DParaCrystal.h b/Core/Samples/inc/InterferenceFunction2DParaCrystal.h index bb7e144e286..809629a30da 100644 --- a/Core/Samples/inc/InterferenceFunction2DParaCrystal.h +++ b/Core/Samples/inc/InterferenceFunction2DParaCrystal.h @@ -21,7 +21,7 @@ class InterferenceFunction2DParaCrystal : public IInterferenceFunction { public: InterferenceFunction2DParaCrystal(double length_1, double length_2, double alpha_lattice, double xi=0.0, double corr_length=0.0); - virtual ~InterferenceFunction2DParaCrystal() {} + virtual ~InterferenceFunction2DParaCrystal(); virtual InterferenceFunction2DParaCrystal *clone() const { InterferenceFunction2DParaCrystal *p_new = new InterferenceFunction2DParaCrystal(m_lattice_lengths[0], m_lattice_lengths[1], m_alpha_lattice, m_xi, m_corr_length); p_new->setDomainSizes(m_domain_sizes[0], m_domain_sizes[1]); @@ -45,6 +45,8 @@ public: virtual double evaluate(const cvector_t &q) const; + //! add parameters from local pool to external pool and call recursion over direct children + virtual std::string addParametersToExternalPool(std::string path, ParameterPool *external_pool, int copy_number=-1) const; protected: void transformToPrincipalAxes(double qx, double qy, double gamma, double delta, double &q_pa_1, double &q_pa_2) const; double m_lattice_lengths[2]; diff --git a/Core/Samples/src/InterferenceFunction2DParaCrystal.cpp b/Core/Samples/src/InterferenceFunction2DParaCrystal.cpp index b8696de3e1b..ad8469994ba 100644 --- a/Core/Samples/src/InterferenceFunction2DParaCrystal.cpp +++ b/Core/Samples/src/InterferenceFunction2DParaCrystal.cpp @@ -31,6 +31,13 @@ InterferenceFunction2DParaCrystal::InterferenceFunction2DParaCrystal(double leng init_parameters(); } +InterferenceFunction2DParaCrystal::~InterferenceFunction2DParaCrystal() +{ + for (size_t i=0; i<2; ++i) { + if (m_pdfs[i]) delete m_pdfs[i]; + } +} + void InterferenceFunction2DParaCrystal::setProbabilityDistributions( const IFTDistribution2D& pdf_1, const IFTDistribution2D& pdf_2) { @@ -45,14 +52,36 @@ double InterferenceFunction2DParaCrystal::evaluate(const cvector_t &q) const { m_qx = q.x().real(); m_qy = q.y().real(); - std::binder1st<std::const_mem_fun1_t<double, InterferenceFunction2DParaCrystal, double> > f_base = std::bind1st(std::mem_fun(&InterferenceFunction2DParaCrystal::interferenceForXi), this); - gsl_function f; - f.function = &wrapFunctionForGSL; - f.params = &f_base; - double result = MathFunctions::Integrate1D(&f, 0.0, M_PI)/M_PI; + double result; + if (m_integrate_xi) { + std::binder1st<std::const_mem_fun1_t<double, InterferenceFunction2DParaCrystal, double> > f_base = std::bind1st(std::mem_fun(&InterferenceFunction2DParaCrystal::interferenceForXi), this); + gsl_function f; + f.function = &wrapFunctionForGSL; + f.params = &f_base; + result = MathFunctions::Integrate1D(&f, 0.0, M_PI)/M_PI; + } + else { + result = interferenceForXi(m_xi); + } return result; } +std::string InterferenceFunction2DParaCrystal::addParametersToExternalPool(std::string path, + ParameterPool* external_pool, int copy_number) const +{ + // add own parameters + std::string new_path = IParameterized::addParametersToExternalPool(path, external_pool, copy_number); + + // add parameters of the probability density functions + if (m_pdfs[0]) { + m_pdfs[0]->addParametersToExternalPool(new_path, external_pool, 0); + } + if (m_pdfs[1]) { + m_pdfs[1]->addParametersToExternalPool(new_path, external_pool, 1); + } + return new_path; +} + void InterferenceFunction2DParaCrystal::transformToPrincipalAxes(double qx, double qy, double gamma, double delta, double& q_pa_1, double& q_pa_2) const { @@ -63,8 +92,13 @@ void InterferenceFunction2DParaCrystal::transformToPrincipalAxes(double qx, void InterferenceFunction2DParaCrystal::init_parameters() { getParameterPool()->clear(); -// getParameterPool()->registerParameter("peak_distance", &m_peak_distance); + getParameterPool()->registerParameter("lattice_length_1", &m_lattice_lengths[0]); + getParameterPool()->registerParameter("lattice_length_2", &m_lattice_lengths[1]); + getParameterPool()->registerParameter("lattice_angle", &m_alpha_lattice); + getParameterPool()->registerParameter("lattice_orientation", &m_xi); getParameterPool()->registerParameter("corr_length", &m_corr_length); + getParameterPool()->registerParameter("domain_size_1", &m_domain_sizes[0]); + getParameterPool()->registerParameter("domain_size_2", &m_domain_sizes[1]); } InterferenceFunction2DParaCrystal* InterferenceFunction2DParaCrystal::createSquare( diff --git a/shared.pri b/shared.pri index eeef573c34d..f5b8c1263ac 100644 --- a/shared.pri +++ b/shared.pri @@ -90,7 +90,7 @@ NumberOfSuchFiles=$$system(ls $${GENERAL_EXTERNALS_DIR}/lib/libboost_thread-mt* QMAKE_CXXFLAGS_DEBUG += -fdiagnostics-show-option # to find out in gcc which option control warning #QMAKE_CXXFLAGS_RELEASE += -O3 -ffast-math -msse3 QMAKE_CXXFLAGS_RELEASE -= -O2 -QMAKE_CXXFLAGS_RELEASE += -O3 -ffast-math +QMAKE_CXXFLAGS_RELEASE += -g -O3 -ffast-math # uncommenting line below produces non-stripped (very large) libraries #QMAKE_STRIP=: -- GitLab