further filtering/resampling code improvements

metadata
Wenzel Jakob 2014-06-06 01:55:02 +02:00
parent 6b3246d15a
commit 8be3692ca6
8 changed files with 335 additions and 175 deletions

View File

@ -1023,8 +1023,8 @@ public:
* horizontal and vertical boundary conditions when looking up data
* outside of the input domain.
*
* A maximum and maximum image value can be specified to prevent to prevent
* out-of-range values that are created by the resampling process.
* A maximum and maximum image value can be specified to prevent to
* prevent out-of-range values that are created by the resampling process.
*
* The optional 'temp' parameter can be used to pass an image of
* resolution <tt>Vector2i(target->getWidth(), this->getHeight())</tt>
@ -1036,25 +1036,6 @@ public:
Bitmap *target, Bitmap *temp = NULL,
Float minValue = 0.0f, Float maxValue = 1.0f) const;
/**
* \brief Apply a separable filter to the image
*
* Applies the provided filter while observing the specified
* horizontal and vertical boundary conditions when looking up data
* outside of the input domain.
*
* A maximum and maximum image value can be specified to prevent to prevent
* out-of-range values that are created by the filtering process.
*
* The optional 'temp' parameter can be used to pass an image of
* the same dimensions as the source and target image.
*/
void filter(const ReconstructionFilter *rfilter,
ReconstructionFilter::EBoundaryCondition bch,
ReconstructionFilter::EBoundaryCondition bcv,
Bitmap *target, Bitmap *temp = NULL,
Float minValue = 0.0f, Float maxValue = 1.0f) const;
/**
* \brief Up- or down-sample this image to a different resolution
*
@ -1064,6 +1045,8 @@ public:
*
* A maximum and maximum image value can be specified to prevent to prevent
* out-of-range values that are created by the resampling process.
*
* This function allocates a new output image and returns it.
*/
ref<Bitmap> resample(const ReconstructionFilter *rfilter,
ReconstructionFilter::EBoundaryCondition bch,
@ -1071,6 +1054,48 @@ public:
const Vector2i &size, Float minValue = 0.0f,
Float maxValue = 1.0f) const;
/**
* \brief Apply a separable convolution filter to the image
*
* Applies the provided filter while observing the specified
* horizontal and vertical boundary conditions when looking up data
* outside of the input domain.
*
* In comparison to \ref convolve(), this function operates
* in the primal domain.
*
* A maximum and maximum image value can be specified to prevent to
* prevent out-of-range values that are created by the filtering process.
*
* The optional 'temp' parameter can be used to pass an image of
* the same dimensions as the source and target image.
*/
void filter(const ReconstructionFilter *rfilter,
ReconstructionFilter::EBoundaryCondition bch,
ReconstructionFilter::EBoundaryCondition bcv,
Bitmap *target, Bitmap *temp = NULL,
Float minValue = 0.0f, Float maxValue = 1.0f) const;
/**
* \brief Apply a separable convolution filter to the image
*
* Applies the provided filter while observing the specified
* horizontal and vertical boundary conditions when looking up data
* outside of the input domain.
*
* In comparison to \ref convolve(), this function operates
* in the primal domain.
*
* A maximum and maximum image value can be specified to prevent to
* prevent out-of-range values that are created by the filtering process.
*
* This function allocates a new output image and returns it.
*/
ref<Bitmap> filter(const ReconstructionFilter *rfilter,
ReconstructionFilter::EBoundaryCondition bch,
ReconstructionFilter::EBoundaryCondition bcv,
Float minValue = 0.0f, Float maxValue = 1.0f) const;
//! @}
// ======================================================================

View File

@ -84,8 +84,8 @@ public:
MTS_DECLARE_CLASS()
protected:
/// Create a new reconstruction filter
ReconstructionFilter(const Properties &props);
/// Create a new reconstruction filter
ReconstructionFilter(const Properties &props);
/// Unserialize a filter
ReconstructionFilter(Stream *stream, InstanceManager *manager);
@ -121,126 +121,88 @@ template <typename Scalar> struct Resampler {
* outside of the defined input domain.
*/
Resampler(const ReconstructionFilter *rfilter, ReconstructionFilter::EBoundaryCondition bc,
int sourceRes, int targetRes) : m_bc(bc), m_sourceRes(sourceRes), m_targetRes(targetRes) {
int sourceRes, int targetRes) : m_bc(bc), m_sourceRes(sourceRes), m_targetRes(targetRes),
m_start(NULL), m_weights(NULL) {
SAssert(sourceRes > 0 && targetRes > 0);
Float filterRadius = rfilter->getRadius(), scale = 1.0f, invScale = 1.0f;
/* Low-pass filter: scale reconstruction filters when downsampling */
if (targetRes < sourceRes) {
scale = (Float) sourceRes / (Float) targetRes;
invScale = 1 / scale;
invScale = 1 / scale;
filterRadius *= scale;
}
m_taps = floorToInt(filterRadius * 2);
m_start = new int[targetRes];
m_weights = new Scalar[m_taps * targetRes];
m_fastStart = 0;
m_fastEnd = m_targetRes;
m_taps = ceilToInt(filterRadius * 2);
if (sourceRes == targetRes && (m_taps % 2) != 1)
--m_taps;
m_halfTaps = m_taps / 2;
for (int i=0; i<targetRes; i++) {
/* Compute the fractional coordinates of the new sample i in the original coordinates */
Float center = (i + (Float) 0.5f) / targetRes * sourceRes;
if (sourceRes != targetRes) { /* Resampling mode */
m_start = new int[targetRes];
m_weights = new Scalar[m_taps * targetRes];
m_fastStart = 0;
m_fastEnd = m_targetRes;
/* Determine the index of the first original sample that might contribute */
m_start[i] = floorToInt(center - filterRadius + (Float) 0.5f);
for (int i=0; i<targetRes; i++) {
/* Compute the fractional coordinates of the new sample i in the original coordinates */
Float center = (i + (Float) 0.5f) / targetRes * sourceRes;
/* Determine the size of center region, on which to run fast non condition-aware code */
if (m_start[i] < 0)
m_fastStart = std::max(m_fastStart, i + 1);
else if (m_start[i] + m_taps - 1 >= m_sourceRes)
m_fastEnd = std::min(m_fastEnd, i - 1);
/* Determine the index of the first original sample that might contribute */
m_start[i] = floorToInt(center - filterRadius + (Float) 0.5f);
/* Determine the size of center region, on which to run fast non condition-aware code */
if (m_start[i] < 0)
m_fastStart = std::max(m_fastStart, i + 1);
else if (m_start[i] + m_taps - 1 >= m_sourceRes)
m_fastEnd = std::min(m_fastEnd, i - 1);
Float sum = 0;
for (int j=0; j<m_taps; j++) {
/* Compute the the position where the filter should be evaluated */
Float pos = m_start[i] + j + (Float) 0.5f - center;
/* Perform the evaluation and record the weight */
Float weight = rfilter->eval(pos * invScale);
m_weights[i * m_taps + j] = (Scalar) weight;
sum += weight;
}
/* Normalize the contribution of each sample */
Float normalization = 1.0f / sum;
for (int j=0; j<m_taps; j++) {
Scalar &value = m_weights[i * m_taps + j];
value = (Scalar) ((Float) value * normalization);
}
}
} else { /* Filtering mode */
m_weights = new Scalar[m_taps];
Float sum = 0;
for (int j=0; j<m_taps; j++) {
/* Compute the the position where the filter should be evaluated */
Float pos = m_start[i] + j + (Float) 0.5f - center;
/* Perform the evaluation and record the weight */
Float weight = rfilter->eval(pos * invScale);
m_weights[i * m_taps + j] = (Scalar) weight;
for (int i=0; i<m_taps; i++) {
Float weight = (Scalar) rfilter->eval((Float) (i-m_halfTaps));
m_weights[i] = weight;
sum += weight;
}
/* Normalize the contribution of each sample */
Float normalization = 1.0f / sum;
for (int j=0; j<m_taps; j++) {
Scalar &value = m_weights[i * m_taps + j];
for (int i=0; i<m_taps; i++) {
Scalar &value = m_weights[i];
value = (Scalar) ((Float) value * normalization);
}
m_fastStart = std::min(m_halfTaps, m_targetRes-1);
m_fastEnd = std::max(m_targetRes-m_halfTaps-1, 0);
}
/* Don't have overlapping fast start/end intervals when
the target image is very small compared to the source image */
m_fastStart = std::min(m_fastStart, m_fastEnd);
}
/// Release all memory
~Resampler() {
delete[] m_start;
delete[] m_weights;
}
/**
* \brief Resample a multi-channel array
*
* \param source
* Source array of samples
* \param target
* Target array of samples
* \param sourceStride
* Stride of samples in the source array. A value
* of '1' implies that they are densely packed.
* \param targetStride
* Stride of samples in the source array. A value
* of '1' implies that they are densely packed.
* \param channels
* Number of channels to be resampled
*/
void resample(const Scalar *source, size_t sourceStride,
Scalar *target, size_t targetStride, int channels) {
const int taps = m_taps;
targetStride = channels * (targetStride - 1);
sourceStride *= channels;
/* Resample the left border region, while accounting for the boundary conditions */
for (int i=0; i<m_fastStart; ++i) {
int start = m_start[i];
for (int ch=0; ch<channels; ++ch) {
Scalar result = 0;
for (int j=0; j<taps; ++j)
result += lookup(source, start + j, sourceStride, ch) * m_weights[i * taps + j];
*target++ = result;
}
target += targetStride;
}
/* Use a faster, vectorizable loop for resampling the main portion */
for (int i=m_fastStart; i<m_fastEnd; ++i) {
int start = m_start[i];
for (int ch=0; ch<channels; ++ch) {
Scalar result = 0;
for (int j=0; j<taps; ++j)
result += source[sourceStride * (start + j) + ch] * m_weights[i * taps + j];
*target++ = result;
}
target += targetStride;
}
/* Resample the right border region, while accounting for the boundary conditions */
for (int i=m_fastEnd; i<m_targetRes; ++i) {
int start = m_start[i];
for (int ch=0; ch<channels; ++ch) {
Scalar result = 0;
for (int j=0; j<taps; ++j)
result += lookup(source, start + j, sourceStride, ch) * m_weights[i * taps + j];
*target++ = result;
}
target += targetStride;
}
if (m_start)
delete[] m_start;
if (m_weights)
delete[] m_weights;
}
/**
@ -270,53 +232,206 @@ template <typename Scalar> struct Resampler {
void resampleAndClamp(const Scalar *source, size_t sourceStride,
Scalar *target, size_t targetStride, int channels,
Scalar min = (Scalar) 0, Scalar max = (Scalar) 1) {
const int taps = m_taps;
const int taps = m_taps, halfTaps = m_halfTaps;
targetStride = channels * (targetStride - 1);
sourceStride *= channels;
/* Resample the left border region, while accounting for the boundary conditions */
for (int i=0; i<m_fastStart; ++i) {
int start = m_start[i];
if (m_start) {
/* Resample the left border region, while accounting for the boundary conditions */
for (int i=0; i<m_fastStart; ++i) {
int start = m_start[i];
for (int ch=0; ch<channels; ++ch) {
Scalar result = 0;
for (int j=0; j<taps; ++j)
result += lookup(source, start + j, sourceStride, ch) * m_weights[i * taps + j];
*target++ = std::min(max, std::max(min, result));
for (int ch=0; ch<channels; ++ch) {
Scalar result = 0;
for (int j=0; j<taps; ++j)
result += lookup(source, start + j, sourceStride, ch) * m_weights[i * taps + j];
*target++ = std::min(max, std::max(min, result));
}
target += targetStride;
}
target += targetStride;
}
/* Use a faster, vectorizable loop for resampling the main portion */
for (int i=m_fastStart; i<m_fastEnd; ++i) {
int start = m_start[i];
/* Use a faster, vectorizable loop for resampling the main portion */
for (int i=m_fastStart; i<m_fastEnd; ++i) {
int start = m_start[i];
for (int ch=0; ch<channels; ++ch) {
Scalar result = 0;
for (int j=0; j<taps; ++j)
result += source[sourceStride * (start + j) + ch] * m_weights[i * taps + j];
*target++ = std::min(max, std::max(min, result));
}
for (int ch=0; ch<channels; ++ch) {
Scalar result = 0;
for (int j=0; j<taps; ++j)
result += source[sourceStride * (start + j) + ch] * m_weights[i * taps + j];
*target++ = std::min(max, std::max(min, result));
target += targetStride;
}
target += targetStride;
}
/* Resample the right border region, while accounting for the boundary conditions */
for (int i=m_fastEnd; i<m_targetRes; ++i) {
int start = m_start[i];
/* Resample the right border region, while accounting for the boundary conditions */
for (int i=m_fastEnd; i<m_targetRes; ++i) {
int start = m_start[i];
for (int ch=0; ch<channels; ++ch) {
Scalar result = 0;
for (int j=0; j<taps; ++j)
result += lookup(source, start + j, sourceStride, ch) * m_weights[i * taps + j];
*target++ = std::min(max, std::max(min, result));
}
for (int ch=0; ch<channels; ++ch) {
Scalar result = 0;
for (int j=0; j<taps; ++j)
result += lookup(source, start + j, sourceStride, ch) * m_weights[i * taps + j];
*target++ = std::min(max, std::max(min, result));
target += targetStride;
}
} else {
/* Filter the left border region, while accounting for the boundary conditions */
for (int i=0; i<m_fastStart; ++i) {
int start = i - halfTaps;
for (int ch=0; ch<channels; ++ch) {
Scalar result = 0;
for (int j=0; j<taps; ++j)
result += lookup(source, start + j, sourceStride, ch) * m_weights[j];
*target++ = std::min(max, std::max(min, result));
}
target += targetStride;
}
target += targetStride;
/* Use a faster, vectorizable loop for filtering the main portion */
for (int i=m_fastStart; i<m_fastEnd; ++i) {
int start = i - halfTaps;
for (int ch=0; ch<channels; ++ch) {
Scalar result = 0;
for (int j=0; j<taps; ++j)
result += source[sourceStride * (start + j) + ch] * m_weights[j];
*target++ = std::min(max, std::max(min, result));
}
target += targetStride;
}
/* Filter the right border region, while accounting for the boundary conditions */
for (int i=m_fastEnd; i<m_targetRes; ++i) {
int start = i - halfTaps;
for (int ch=0; ch<channels; ++ch) {
Scalar result = 0;
for (int j=0; j<taps; ++j)
result += lookup(source, start + j, sourceStride, ch) * m_weights[j];
*target++ = std::min(max, std::max(min, result));
}
target += targetStride;
}
}
}
/**
* \brief Resample a multi-channel array
*
* \param source
* Source array of samples
* \param target
* Target array of samples
* \param sourceStride
* Stride of samples in the source array. A value
* of '1' implies that they are densely packed.
* \param targetStride
* Stride of samples in the source array. A value
* of '1' implies that they are densely packed.
* \param channels
* Number of channels to be resampled
*/
void resample(const Scalar *source, size_t sourceStride,
Scalar *target, size_t targetStride, int channels) {
const int taps = m_taps, halfTaps = m_halfTaps;
targetStride = channels * (targetStride - 1);
sourceStride *= channels;
if (m_start) {
/* Resample the left border region, while accounting for the boundary conditions */
for (int i=0; i<m_fastStart; ++i) {
int start = m_start[i];
for (int ch=0; ch<channels; ++ch) {
Scalar result = 0;
for (int j=0; j<taps; ++j)
result += lookup(source, start + j, sourceStride, ch) * m_weights[i * taps + j];
*target++ = result;
}
target += targetStride;
}
/* Use a faster, vectorizable loop for resampling the main portion */
for (int i=m_fastStart; i<m_fastEnd; ++i) {
int start = m_start[i];
for (int ch=0; ch<channels; ++ch) {
Scalar result = 0;
for (int j=0; j<taps; ++j)
result += source[sourceStride * (start + j) + ch] * m_weights[i * taps + j];
*target++ = result;
}
target += targetStride;
}
/* Resample the right border region, while accounting for the boundary conditions */
for (int i=m_fastEnd; i<m_targetRes; ++i) {
int start = m_start[i];
for (int ch=0; ch<channels; ++ch) {
Scalar result = 0;
for (int j=0; j<taps; ++j)
result += lookup(source, start + j, sourceStride, ch) * m_weights[i * taps + j];
*target++ = result;
}
target += targetStride;
}
} else {
/* Filter the left border region, while accounting for the boundary conditions */
for (int i=0; i<m_fastStart; ++i) {
int start = i - halfTaps;
for (int ch=0; ch<channels; ++ch) {
Scalar result = 0;
for (int j=0; j<taps; ++j)
result += lookup(source, start + j, sourceStride, ch) * m_weights[j];
*target++ = result;
}
target += targetStride;
}
/* Use a faster, vectorizable loop for filtering the main portion */
for (int i=m_fastStart; i<m_fastEnd; ++i) {
int start = i - halfTaps;
for (int ch=0; ch<channels; ++ch) {
Scalar result = 0;
for (int j=0; j<taps; ++j)
result += source[sourceStride * (start + j) + ch] * m_weights[j];
*target++ = result;
}
target += targetStride;
}
/* Filter the right border region, while accounting for the boundary conditions */
for (int i=m_fastEnd; i<m_targetRes; ++i) {
int start = i - halfTaps;
for (int ch=0; ch<channels; ++ch) {
Scalar result = 0;
for (int j=0; j<taps; ++j)
result += lookup(source, start + j, sourceStride, ch) * m_weights[j];
*target++ = result;
}
target += targetStride;
}
}
}
private:
FINLINE Scalar lookup(const Scalar *source, int pos, size_t stride, int offset) const {
@ -349,7 +464,7 @@ private:
int *m_start;
Scalar *m_weights;
int m_fastStart, m_fastEnd;
int m_taps;
int m_taps, m_halfTaps;
};
extern MTS_EXPORT_CORE std::ostream &operator<<(std::ostream &os, const ReconstructionFilter::EBoundaryCondition &value);

View File

@ -51,8 +51,8 @@ extern MTS_EXPORT_CORE std::string indent(const std::string &string, int amount=
extern MTS_EXPORT_CORE std::string formatString(const char *pFmt, ...);
/**
* \brief Convert a time difference (in ms) to a string representation
* \param time Time value in milliseconds
* \brief Convert a time difference (in seconds) to a string representation
* \param time Time difference in (fractional) sections
* \param precise When set to true, a higher-precision string representation
* is generated.
*/

View File

@ -121,7 +121,7 @@ public:
/// Move the intersection forward or backward through time
inline void adjustTime(Float time);
/// Calls the suitable implementaiton of \ref Shape::getNormalDerivative()
/// Calls the suitable implementation of \ref Shape::getNormalDerivative()
inline void getNormalDerivative(Vector &dndu, Vector &dndv,
bool shadingFrame = true) const;

View File

@ -23,6 +23,7 @@
#include <boost/algorithm/string.hpp>
#include <boost/scoped_array.hpp>
#include <boost/thread/mutex.hpp>
#include <set>
#if defined(__WINDOWS__)
#undef _CRT_SECURE_NO_WARNINGS
@ -2216,12 +2217,12 @@ void Bitmap::applyMatrix(Float matrix_[3][3]) {
}
}
/// Bitmap resampling utility function
/// Bitmap filtering / resampling utility function
template <typename Scalar> static void resample(ref<const ReconstructionFilter> rfilter,
ReconstructionFilter::EBoundaryCondition bch,
ReconstructionFilter::EBoundaryCondition bcv,
const Bitmap *source, Bitmap *target, ref<Bitmap> temp, Float minValue,
Float maxValue, bool forceFiltering) {
Float maxValue, bool filter) {
if (!rfilter) {
/* Resample using a 2-lobed Lanczos reconstruction filter */
@ -2235,7 +2236,7 @@ template <typename Scalar> static void resample(ref<const ReconstructionFilter>
}
if (source->getHeight() == target->getHeight() &&
source->getWidth() == target->getWidth() && !forceFiltering) {
source->getWidth() == target->getWidth() && !filter) {
memcpy(target->getData(), source->getData(), source->getBufferSize());
return;
}
@ -2245,13 +2246,13 @@ template <typename Scalar> static void resample(ref<const ReconstructionFilter>
minValue != -std::numeric_limits<Float>::infinity() ||
maxValue != std::numeric_limits<Float>::infinity();
if (source->getWidth() != target->getWidth() || forceFiltering) {
if (source->getWidth() != target->getWidth() || filter) {
/* Re-sample along the X direction */
Resampler<Scalar> r(rfilter, bch, source->getWidth(), target->getWidth());
/* Create a bitmap for intermediate storage */
if (!temp) {
if (source->getHeight() == target->getHeight() && !forceFiltering)
if (source->getHeight() == target->getHeight() && !filter)
temp = target; // write directly to the output bitmap
else // otherwise: write to a temporary bitmap
temp = new Bitmap(source->getPixelFormat(), source->getComponentFormat(),
@ -2289,7 +2290,7 @@ template <typename Scalar> static void resample(ref<const ReconstructionFilter>
source = temp;
}
if (source->getHeight() != target->getHeight() || forceFiltering) {
if (source->getHeight() != target->getHeight() || filter) {
/* Re-sample along the Y direction */
Resampler<Scalar> r(rfilter, bcv, source->getHeight(), target->getHeight());
@ -2370,7 +2371,6 @@ void Bitmap::filter(const ReconstructionFilter *rfilter,
}
}
ref<Bitmap> Bitmap::resample(const ReconstructionFilter *rfilter,
ReconstructionFilter::EBoundaryCondition bch,
ReconstructionFilter::EBoundaryCondition bcv,
@ -2378,10 +2378,23 @@ ref<Bitmap> Bitmap::resample(const ReconstructionFilter *rfilter,
ref<Bitmap> result = new Bitmap(m_pixelFormat, m_componentFormat, size);
result->m_metadata = m_metadata;
result->m_gamma = m_gamma;
result->m_channelNames = m_channelNames;
resample(rfilter, bch, bcv, result, NULL, minValue, maxValue);
return result;
}
ref<Bitmap> Bitmap::filter(const ReconstructionFilter *rfilter,
ReconstructionFilter::EBoundaryCondition bch,
ReconstructionFilter::EBoundaryCondition bcv,
Float minValue, Float maxValue) const {
ref<Bitmap> result = new Bitmap(m_pixelFormat, m_componentFormat, getSize());
result->m_metadata = m_metadata;
result->m_gamma = m_gamma;
result->m_channelNames = m_channelNames;
filter(rfilter, bch, bcv, result, NULL, minValue, maxValue);
return result;
}
bool Bitmap::operator==(const Bitmap &bitmap) const {
return m_pixelFormat == bitmap.m_pixelFormat &&
m_componentFormat == bitmap.m_componentFormat &&

View File

@ -26,6 +26,11 @@
#pragma warning(disable : 4267) // 'return' : conversion from 'size_t' to 'long', possible loss of data
#endif
#define BP_RETURN_CONSTREF bp::return_value_policy<bp::copy_const_reference>()
#define BP_RETURN_NONCONSTREF bp::return_value_policy<bp::copy_non_const_reference>()
#define BP_RETURN_VALUE bp::return_value_policy<bp::return_by_value>()
#define BP_RETURN_INTREF bp::return_internal_reference<>()
#define BP_STRUCT_DECL(Name, Init) \
bp::class_<Name> Name ##_struct(#Name, Init); \
bp::register_ptr_to_python<Name*>();
@ -149,11 +154,6 @@
bp::detail::current_scope = value.ptr(); \
} while (0);
#define BP_RETURN_CONSTREF bp::return_value_policy<bp::copy_const_reference>()
#define BP_RETURN_NONCONSTREF bp::return_value_policy<bp::copy_non_const_reference>()
#define BP_RETURN_VALUE bp::return_value_policy<bp::return_by_value>()
#define BP_RETURN_INTREF bp::return_internal_reference<>()
namespace boost {
namespace python {
template <typename T> T* get_pointer(mitsuba::ref<T> & p) {

View File

@ -604,13 +604,13 @@ static ref<ConfigurableObject> pluginmgr_create_helper(PluginManager *manager, b
return object;
}
static ref<ConfigurableObject> pluginmgr_create(PluginManager *manager, bp::dict dict) {
static bp::object pluginmgr_create(PluginManager *manager, bp::dict dict) {
std::map<std::string, ConfigurableObject *> objs;
ref<ConfigurableObject> result = pluginmgr_create_helper(manager, dict, objs);
for (std::map<std::string, ConfigurableObject *>::iterator it = objs.begin();
it != objs.end(); ++it)
it->second->decRef();
return result;
return cast(result);
}
static bp::tuple mkCoordinateSystem(const Vector &n) {
@ -1099,8 +1099,10 @@ BOOST_PYTHON_MEMBER_FUNCTION_OVERLOADS(fromLinearRGB_overloads, fromLinearRGB, 3
BOOST_PYTHON_MEMBER_FUNCTION_OVERLOADS(fromXYZ_overloads, fromXYZ, 3, 4)
BOOST_PYTHON_MEMBER_FUNCTION_OVERLOADS(fromIPT_overloads, fromIPT, 3, 4)
BOOST_PYTHON_MEMBER_FUNCTION_OVERLOADS(reset_overloads, reset, 0, 1)
BOOST_PYTHON_MEMBER_FUNCTION_OVERLOADS(filter_overloads, filter, 4, 7)
BOOST_PYTHON_MEMBER_FUNCTION_OVERLOADS(resample_overloads, resample, 4, 7)
BOOST_PYTHON_MEMBER_FUNCTION_OVERLOADS(resample1_overloads, resample, 4, 7)
BOOST_PYTHON_MEMBER_FUNCTION_OVERLOADS(resample2_overloads, resample, 6, 4)
BOOST_PYTHON_MEMBER_FUNCTION_OVERLOADS(filter1_overloads, filter, 4, 7)
BOOST_PYTHON_MEMBER_FUNCTION_OVERLOADS(filter2_overloads, filter, 3, 5)
#define IMPLEMENT_ANIMATION_TRACK(Name) \
BP_CLASS(Name, AbstractAnimationTrack, (bp::init<AbstractAnimationTrack::EType, size_t>())) \
@ -1475,13 +1477,17 @@ void export_core() {
ReconstructionFilter::EBoundaryCondition, ReconstructionFilter::EBoundaryCondition,
Bitmap *, Bitmap *, Float, Float) const = &Bitmap::resample;
ref<Bitmap> (Bitmap::*resample_2)(const ReconstructionFilter *,
ReconstructionFilter::EBoundaryCondition, ReconstructionFilter::EBoundaryCondition,
const Vector2i &, Float, Float) const = &Bitmap::resample;
void (Bitmap::*filter_1)(const ReconstructionFilter *,
ReconstructionFilter::EBoundaryCondition, ReconstructionFilter::EBoundaryCondition,
Bitmap *, Bitmap *, Float, Float) const = &Bitmap::filter;
ref<Bitmap> (Bitmap::*resample_2)(const ReconstructionFilter *,
ref<Bitmap> (Bitmap::*filter_2)(const ReconstructionFilter *,
ReconstructionFilter::EBoundaryCondition, ReconstructionFilter::EBoundaryCondition,
const Vector2i &, Float, Float) const = &Bitmap::resample;
Float, Float) const = &Bitmap::filter;
const std::vector<std::string> & (Bitmap::*getChannelNames_1)() const = &Bitmap::getChannelNames;
@ -1541,9 +1547,10 @@ void export_core() {
.def("copyFrom", copyFrom_3)
.def("convolve", &Bitmap::convolve)
.def("arithmeticOperation", &Bitmap::arithmeticOperation, BP_RETURN_VALUE)
.def("filter", filter_1, filter_overloads())
.def("resample", resample_1, resample_overloads())
.def("resample", resample_2, BP_RETURN_VALUE)
.def("resample", resample_1, resample1_overloads())
.def("resample", resample_2, resample2_overloads()[BP_RETURN_VALUE])
.def("filter", filter_1, filter1_overloads())
.def("filter", filter_2, filter2_overloads()[BP_RETURN_VALUE])
.def("setGamma", &Bitmap::setGamma)
.def("getGamma", &Bitmap::getGamma)
.def("setMetadataString", &Bitmap::setMetadataString)

View File

@ -48,7 +48,7 @@ public:
}
std::string toString() const {
return "BoxFilter[]";
return formatString("GaussianFilter[radius=%f]", m_radius);
}
MTS_DECLARE_CLASS()