Skip to content
Snippets Groups Projects
Commit dccbdc3b authored by Nicolas Pope's avatar Nicolas Pope
Browse files

Implements #224 to use discontinuity mask

parent 4025ac5d
No related branches found
No related tags found
No related merge requests found
Showing
with 268 additions and 22 deletions
......@@ -215,11 +215,11 @@ bool ILW::_phase0(ftl::rgbd::FrameSet &fs, cudaStream_t stream) {
f.createTexture<float>(Channel::Depth2, Format<float>(f.get<GpuMat>(Channel::Colour).size()));
f.createTexture<float>(Channel::Confidence, Format<float>(f.get<GpuMat>(Channel::Colour).size()));
f.createTexture<int>(Channel::Mask, Format<int>(f.get<GpuMat>(Channel::Colour).size()));
//f.createTexture<int>(Channel::Mask, Format<int>(f.get<GpuMat>(Channel::Colour).size()));
f.createTexture<uchar4>(Channel::Colour);
f.createTexture<float>(Channel::Depth);
f.get<GpuMat>(Channel::Mask).setTo(cv::Scalar(0));
//f.get<GpuMat>(Channel::Mask).setTo(cv::Scalar(0));
}
//cudaSafeCall(cudaStreamSynchronize(stream_));
......@@ -232,7 +232,7 @@ bool ILW::_phase1(ftl::rgbd::FrameSet &fs, int win, cudaStream_t stream) {
cv::cuda::Stream cvstream = cv::cuda::StreamAccessor::wrapStream(stream);
// Find discontinuity mask
for (size_t i=0; i<fs.frames.size(); ++i) {
/*for (size_t i=0; i<fs.frames.size(); ++i) {
auto &f = fs.frames[i];
auto s = fs.sources[i];
......@@ -243,7 +243,7 @@ bool ILW::_phase1(ftl::rgbd::FrameSet &fs, int win, cudaStream_t stream) {
discon_mask_,
stream
);
}
}*/
// First do any preprocessing
if (fill_depth_) {
......
......@@ -35,6 +35,7 @@
#include <ftl/operators/normals.hpp>
#include <ftl/operators/filling.hpp>
#include <ftl/operators/segmentation.hpp>
#include <ftl/operators/mask.hpp>
#include <ftl/cuda/normals.hpp>
#include <ftl/registration.hpp>
......@@ -305,13 +306,15 @@ static void run(ftl::Configurable *root) {
// Create the source depth map pipeline
auto *pipeline1 = ftl::config::create<ftl::operators::Graph>(root, "pre_filters");
pipeline1->append<ftl::operators::ColourChannels>("colour"); // Convert BGR to BGRA
pipeline1->append<ftl::operators::HFSmoother>("hfnoise"); // Remove high-frequency noise
//pipeline1->append<ftl::operators::HFSmoother>("hfnoise"); // Remove high-frequency noise
pipeline1->append<ftl::operators::Normals>("normals"); // Estimate surface normals
//pipeline1->append<ftl::operators::SmoothChannel>("smoothing"); // Generate a smoothing channel
//pipeline1->append<ftl::operators::ScanFieldFill>("filling"); // Generate a smoothing channel
pipeline1->append<ftl::operators::CrossSupport>("cross");
pipeline1->append<ftl::operators::DiscontinuityMask>("discontinuity");
pipeline1->append<ftl::operators::CullDiscontinuity>("remove_discontinuity");
pipeline1->append<ftl::operators::ColourMLS>("mls"); // Perform MLS (using smoothing channel)
pipeline1->append<ftl::operators::VisCrossSupport>("viscross");
pipeline1->append<ftl::operators::VisCrossSupport>("viscross")->set("enabled", false);
// Alignment
......
......@@ -14,6 +14,8 @@ set(OPERSRC
src/disparity/bilateral_filter.cpp
src/segmentation.cu
src/segmentation.cpp
src/mask.cu
src/mask.cpp
)
if (HAVE_OPTFLOW)
......
#ifndef _FTL_OPERATORS_MASK_HPP_
#define _FTL_OPERATORS_MASK_HPP_
#include <ftl/operators/operator.hpp>
#include <ftl/cuda_common.hpp>
namespace ftl {
namespace operators {
/**
* Generate a masking channel that indicates depth discontinuities within a
* specified radius from a pixel. This is useful for culling potentially bad
* depth and colour values when merging and smoothing.
*/
class DiscontinuityMask : public ftl::operators::Operator {
public:
explicit DiscontinuityMask(ftl::Configurable*);
~DiscontinuityMask();
inline Operator::Type type() const override { return Operator::Type::OneToOne; }
bool apply(ftl::rgbd::Frame &in, ftl::rgbd::Frame &out, ftl::rgbd::Source *src, cudaStream_t stream) override;
};
/**
* Remove depth values marked with the discontinuity mask.
*/
class CullDiscontinuity : public ftl::operators::Operator {
public:
explicit CullDiscontinuity(ftl::Configurable*);
~CullDiscontinuity();
inline Operator::Type type() const override { return Operator::Type::OneToOne; }
bool apply(ftl::rgbd::Frame &in, ftl::rgbd::Frame &out, ftl::rgbd::Source *src, cudaStream_t stream) override;
};
}
}
#endif // _FTL_OPERATORS_MASK_HPP_
......@@ -99,6 +99,7 @@ class Graph : public ftl::Configurable {
private:
std::list<ftl::operators::detail::OperatorNode> operators_;
std::map<std::string, ftl::Configurable*> configs_;
cudaStream_t stream_;
ftl::Configurable *_append(ftl::operators::detail::ConstructionHelperBase*);
};
......
......@@ -12,6 +12,8 @@ ColourChannels::~ColourChannels() {
}
bool ColourChannels::apply(ftl::rgbd::Frame &in, ftl::rgbd::Frame &out, ftl::rgbd::Source *s, cudaStream_t stream) {
auto cvstream = cv::cuda::StreamAccessor::wrapStream(stream);
// Convert colour from BGR to BGRA if needed
if (in.get<cv::cuda::GpuMat>(Channel::Colour).type() == CV_8UC3) {
//cv::cuda::Stream cvstream = cv::cuda::StreamAccessor::wrapStream(stream);
......@@ -19,7 +21,7 @@ bool ColourChannels::apply(ftl::rgbd::Frame &in, ftl::rgbd::Frame &out, ftl::rgb
auto &col = in.get<cv::cuda::GpuMat>(Channel::Colour);
temp_.create(col.size(), CV_8UC4);
cv::cuda::swap(col, temp_);
cv::cuda::cvtColor(temp_,col, cv::COLOR_BGR2BGRA, 0);
cv::cuda::cvtColor(temp_,col, cv::COLOR_BGR2BGRA, 0, cvstream);
}
return true;
......
......@@ -22,7 +22,7 @@ bool ScanFieldFill::apply(ftl::rgbd::Frame &in, ftl::rgbd::Frame &out, ftl::rgbd
in.createTexture<float>(Channel::Depth),
in.createTexture<float>(Channel::Smoothing),
thresh,
s->parameters(), 0
s->parameters(), stream
);
return true;
......
#include <ftl/operators/mask.hpp>
#include "mask_cuda.hpp"
using ftl::operators::DiscontinuityMask;
using ftl::operators::CullDiscontinuity;
using ftl::codecs::Channel;
using ftl::rgbd::Format;
DiscontinuityMask::DiscontinuityMask(ftl::Configurable *cfg) : ftl::operators::Operator(cfg) {
}
DiscontinuityMask::~DiscontinuityMask() {
}
bool DiscontinuityMask::apply(ftl::rgbd::Frame &in, ftl::rgbd::Frame &out, ftl::rgbd::Source *s, cudaStream_t stream) {
int radius = config()->value("radius", 2);
float threshold = config()->value("threshold", 0.1f);
ftl::cuda::discontinuity(
out.createTexture<int>(Channel::Mask, ftl::rgbd::Format<int>(in.get<cv::cuda::GpuMat>(Channel::Depth).size())),
in.createTexture<uchar4>(Channel::Support1),
in.createTexture<float>(Channel::Depth),
s->parameters(), radius, threshold, stream
);
return true;
}
CullDiscontinuity::CullDiscontinuity(ftl::Configurable *cfg) : ftl::operators::Operator(cfg) {
}
CullDiscontinuity::~CullDiscontinuity() {
}
bool CullDiscontinuity::apply(ftl::rgbd::Frame &in, ftl::rgbd::Frame &out, ftl::rgbd::Source *s, cudaStream_t stream) {
ftl::cuda::cull_discontinuity(
in.createTexture<int>(Channel::Mask),
out.createTexture<float>(Channel::Depth),
stream
);
return true;
}
\ No newline at end of file
#include "mask_cuda.hpp"
#define T_PER_BLOCK 8
using ftl::cuda::Mask;
__global__ void discontinuity_kernel(ftl::cuda::TextureObject<int> mask_out, ftl::cuda::TextureObject<uchar4> support, ftl::cuda::TextureObject<float> depth, ftl::rgbd::Camera params, float threshold, int radius) {
const unsigned int x = blockIdx.x*blockDim.x + threadIdx.x;
const unsigned int y = blockIdx.y*blockDim.y + threadIdx.y;
if (x < params.width && y < params.height) {
Mask mask(0);
const float d = depth.tex2D((int)x, (int)y);
if (d >= params.minDepth && d <= params.maxDepth) {
/* Orts-Escolano S. et al. 2016. Holoportation: Virtual 3D teleportation in real-time. */
// If colour cross support region terminates within the requested
// radius, and the absolute depth difference on the other side is
// greater than threshold, then is is a discontinuity.
// Repeat for left, right, up and down.
const uchar4 sup = support.tex2D((int)x, (int)y);
if (sup.x <= radius) {
float dS = depth.tex2D((int)x - sup.x - radius, (int)y);
if (fabs(dS - d) > threshold) mask.isDiscontinuity(true);
}
if (sup.y <= radius) {
float dS = depth.tex2D((int)x + sup.y + radius, (int)y);
if (fabs(dS - d) > threshold) mask.isDiscontinuity(true);
}
if (sup.z <= radius) {
float dS = depth.tex2D((int)x, (int)y - sup.z - radius);
if (fabs(dS - d) > threshold) mask.isDiscontinuity(true);
}
if (sup.w <= radius) {
float dS = depth.tex2D((int)x, (int)y + sup.w + radius);
if (fabs(dS - d) > threshold) mask.isDiscontinuity(true);
}
}
mask_out(x,y) = (int)mask;
}
}
void ftl::cuda::discontinuity(ftl::cuda::TextureObject<int> &mask_out, ftl::cuda::TextureObject<uchar4> &support, ftl::cuda::TextureObject<float> &depth, const ftl::rgbd::Camera &params, int discon, float thresh, cudaStream_t stream) {
const dim3 gridSize((params.width + T_PER_BLOCK - 1)/T_PER_BLOCK, (params.height + T_PER_BLOCK - 1)/T_PER_BLOCK);
const dim3 blockSize(T_PER_BLOCK, T_PER_BLOCK);
discontinuity_kernel<<<gridSize, blockSize, 0, stream>>>(mask_out, support, depth, params, thresh, discon);
cudaSafeCall( cudaGetLastError() );
#ifdef _DEBUG
cudaSafeCall(cudaDeviceSynchronize());
#endif
}
__global__ void cull_discontinuity_kernel(ftl::cuda::TextureObject<int> mask, ftl::cuda::TextureObject<float> depth) {
const unsigned int x = blockIdx.x*blockDim.x + threadIdx.x;
const unsigned int y = blockIdx.y*blockDim.y + threadIdx.y;
if (x < depth.width() && y < depth.height()) {
Mask m(mask.tex2D((int)x,(int)y));
if (m.isDiscontinuity()) depth(x,y) = 1000.0f;
}
}
void ftl::cuda::cull_discontinuity(ftl::cuda::TextureObject<int> &mask, ftl::cuda::TextureObject<float> &depth, cudaStream_t stream) {
const dim3 gridSize((depth.width() + T_PER_BLOCK - 1)/T_PER_BLOCK, (depth.height() + T_PER_BLOCK - 1)/T_PER_BLOCK);
const dim3 blockSize(T_PER_BLOCK, T_PER_BLOCK);
cull_discontinuity_kernel<<<gridSize, blockSize, 0, stream>>>(mask, depth);
cudaSafeCall( cudaGetLastError() );
#ifdef _DEBUG
cudaSafeCall(cudaDeviceSynchronize());
#endif
}
#ifndef _FTL_CUDA_MASK_HPP_
#define _FTL_CUDA_MASK_HPP_
#include <ftl/cuda_common.hpp>
#include <ftl/rgbd/camera.hpp>
namespace ftl {
namespace cuda {
/**
* Wrap an int mask value used to flag individual depth pixels.
*/
class Mask {
public:
__device__ inline Mask() : v_(0) {}
__device__ explicit inline Mask(int v) : v_(v) {}
#ifdef __CUDACC__
__device__ inline Mask(const ftl::cuda::TextureObject<int> &m, int x, int y) : v_(m.tex2D(x,y)) {}
#endif
__device__ inline operator int() const { return v_; }
__device__ inline bool is(int m) const { return v_ & m; }
__device__ inline bool isFilled() const { return v_ & kMask_Filled; }
__device__ inline bool isDiscontinuity() const { return v_ & kMask_Discontinuity; }
__device__ inline bool hasCorrespondence() const { return v_ & kMask_Correspondence; }
__device__ inline bool isBad() const { return v_ & kMask_Bad; }
__device__ inline void isFilled(bool v) { v_ = (v) ? v_ | kMask_Filled : v_ & (~kMask_Filled); }
__device__ inline void isDiscontinuity(bool v) { v_ = (v) ? v_ | kMask_Discontinuity : v_ & (~kMask_Discontinuity); }
__device__ inline void hasCorrespondence(bool v) { v_ = (v) ? v_ | kMask_Correspondence : v_ & (~kMask_Correspondence); }
__device__ inline void isBad(bool v) { v_ = (v) ? v_ | kMask_Bad : v_ & (~kMask_Bad); }
static constexpr int kMask_Filled = 0x0001;
static constexpr int kMask_Discontinuity = 0x0002;
static constexpr int kMask_Correspondence = 0x0004;
static constexpr int kMask_Bad = 0x0008;
private:
int v_;
};
void discontinuity(
ftl::cuda::TextureObject<int> &mask,
ftl::cuda::TextureObject<uchar4> &support,
ftl::cuda::TextureObject<float> &depth,
const ftl::rgbd::Camera &params,
int radius, float threshold,
cudaStream_t stream);
void cull_discontinuity(
ftl::cuda::TextureObject<int> &mask,
ftl::cuda::TextureObject<float> &depth,
cudaStream_t stream);
}
}
#endif
......@@ -28,7 +28,7 @@ bool Normals::apply(ftl::rgbd::Frame &in, ftl::rgbd::Frame &out, ftl::rgbd::Sour
ftl::cuda::normals(
out.createTexture<float4>(Channel::Normals, ftl::rgbd::Format<float4>(in.get<cv::cuda::GpuMat>(Channel::Depth).size())),
in.createTexture<float>(Channel::Depth),
s->parameters(), 0
s->parameters(), stream
);
return true;
......
......@@ -31,16 +31,18 @@ bool Operator::apply(FrameSet &in, Frame &out, Source *os, cudaStream_t stream)
Graph::Graph(nlohmann::json &config) : ftl::Configurable(config) {
cudaSafeCall( cudaStreamCreate(&stream_) );
}
Graph::~Graph() {
cudaStreamDestroy(stream_);
}
bool Graph::apply(FrameSet &in, FrameSet &out, cudaStream_t stream) {
if (!value("enabled", true)) return false;
auto stream_actual = (stream == 0) ? stream_ : stream;
if (in.frames.size() != out.frames.size()) return false;
for (auto &i : operators_) {
......@@ -53,11 +55,16 @@ bool Graph::apply(FrameSet &in, FrameSet &out, cudaStream_t stream) {
auto *instance = i.instances[j];
if (instance->enabled()) {
instance->apply(in.frames[j], out.frames[j], in.sources[j], stream);
instance->apply(in.frames[j], out.frames[j], in.sources[j], stream_actual);
}
}
}
if (stream == 0) {
cudaStreamSynchronize(stream_actual);
cudaSafeCall( cudaGetLastError() );
}
return true;
}
......
......@@ -22,7 +22,7 @@ bool CrossSupport::apply(ftl::rgbd::Frame &in, ftl::rgbd::Frame &out, ftl::rgbd:
out.createTexture<uchar4>(Channel::Support2, ftl::rgbd::Format<uchar4>(in.get<cv::cuda::GpuMat>(Channel::Colour).size())),
config()->value("depth_tau", 0.04f),
config()->value("v_max", 5),
config()->value("h_max", 5), 0
config()->value("h_max", 5), stream
);
} //else {
ftl::cuda::support_region(
......@@ -30,7 +30,7 @@ bool CrossSupport::apply(ftl::rgbd::Frame &in, ftl::rgbd::Frame &out, ftl::rgbd:
out.createTexture<uchar4>(Channel::Support1, ftl::rgbd::Format<uchar4>(in.get<cv::cuda::GpuMat>(Channel::Colour).size())),
config()->value("tau", 5.0f),
config()->value("v_max", 5),
config()->value("h_max", 5), 0
config()->value("h_max", 5), stream
);
//}
......@@ -62,7 +62,7 @@ bool VisCrossSupport::apply(ftl::rgbd::Frame &in, ftl::rgbd::Frame &out, ftl::rg
in.createTexture<float>(Channel::Depth),
in.createTexture<uchar4>(Channel::Support1),
in.createTexture<uchar4>(Channel::Support2),
0
stream
);
} else {
ftl::cuda::vis_support_region(
......@@ -74,7 +74,7 @@ bool VisCrossSupport::apply(ftl::rgbd::Frame &in, ftl::rgbd::Frame &out, ftl::rg
config()->value("offset_y", 0),
config()->value("spacing_x", 50),
config()->value("spacing_y", 50),
0
stream
);
if (show_depth) {
......@@ -87,7 +87,7 @@ bool VisCrossSupport::apply(ftl::rgbd::Frame &in, ftl::rgbd::Frame &out, ftl::rg
config()->value("offset_y", 0),
config()->value("spacing_x", 50),
config()->value("spacing_y", 50),
0
stream
);
}
}
......
......@@ -33,7 +33,7 @@ bool HFSmoother::apply(ftl::rgbd::Frame &in, ftl::rgbd::Frame &out, ftl::rgbd::S
in.createTexture<float>(Channel::Energy, ftl::rgbd::Format<float>(in.get<cv::cuda::GpuMat>(Channel::Depth).size())),
in.createTexture<float>(Channel::Smoothing, ftl::rgbd::Format<float>(in.get<cv::cuda::GpuMat>(Channel::Depth).size())),
var_thresh,
s->parameters(), 0
s->parameters(), stream
);
}
......@@ -158,7 +158,7 @@ bool SimpleMLS::apply(ftl::rgbd::Frame &in, ftl::rgbd::Frame &out, ftl::rgbd::So
thresh,
radius,
s->parameters(),
0
stream
);
in.swapChannels(Channel::Depth, Channel::Depth2);
......@@ -204,7 +204,7 @@ bool ColourMLS::apply(ftl::rgbd::Frame &in, ftl::rgbd::Frame &out, ftl::rgbd::So
col_smooth,
radius,
s->parameters(),
0
stream
);
} else {
ftl::cuda::colour_mls_smooth_csr(
......@@ -218,7 +218,7 @@ bool ColourMLS::apply(ftl::rgbd::Frame &in, ftl::rgbd::Frame &out, ftl::rgbd::So
col_smooth,
filling,
s->parameters(),
0
stream
);
}
......@@ -259,7 +259,7 @@ bool AdaptiveMLS::apply(ftl::rgbd::Frame &in, ftl::rgbd::Frame &out, ftl::rgbd::
in.createTexture<float>(Channel::Smoothing),
radius,
s->parameters(),
0
stream
);
in.swapChannels(Channel::Depth, Channel::Depth2);
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment