mirror of
https://github.com/opencv/opencv.git
synced 2025-01-18 14:13:15 +08:00
Add OpenCL accelerated implementation for Retina.
This commit is contained in:
parent
d8c6e89a54
commit
29eefe52bb
@ -1,2 +1,2 @@
|
||||
set(the_description "Biologically inspired algorithms")
|
||||
ocv_define_module(bioinspired opencv_core OPTIONAL opencv_highgui)
|
||||
ocv_define_module(bioinspired opencv_core OPTIONAL opencv_highgui opencv_ocl)
|
||||
|
@ -304,7 +304,8 @@ public:
|
||||
CV_EXPORTS Ptr<Retina> createRetina(Size inputSize);
|
||||
CV_EXPORTS Ptr<Retina> createRetina(Size inputSize, const bool colorMode, int colorSamplingMethod=RETINA_COLOR_BAYER, const bool useRetinaLogSampling=false, const double reductionFactor=1.0, const double samplingStrenght=10.0);
|
||||
|
||||
|
||||
CV_EXPORTS Ptr<Retina> createRetina_OCL(Size inputSize);
|
||||
CV_EXPORTS Ptr<Retina> createRetina_OCL(Size inputSize, const bool colorMode, int colorSamplingMethod=RETINA_COLOR_BAYER, const bool useRetinaLogSampling=false, const double reductionFactor=1.0, const double samplingStrenght=10.0);
|
||||
}
|
||||
}
|
||||
#endif /* __OPENCV_BIOINSPIRED_RETINA_HPP__ */
|
||||
|
753
modules/bioinspired/src/opencl/retina_kernel.cl
Normal file
753
modules/bioinspired/src/opencl/retina_kernel.cl
Normal file
@ -0,0 +1,753 @@
|
||||
/*M///////////////////////////////////////////////////////////////////////////////////////
|
||||
//
|
||||
// IMPORTANT: READ BEFORE DOWNLOADING, COPYING, INSTALLING OR USING.
|
||||
//
|
||||
// By downloading, copying, installing or using the software you agree to this license.
|
||||
// If you do not agree to this license, do not download, install,
|
||||
// copy or use the software.
|
||||
//
|
||||
//
|
||||
// License Agreement
|
||||
// For Open Source Computer Vision Library
|
||||
//
|
||||
// Copyright (C) 2010-2013, Multicoreware, Inc., all rights reserved.
|
||||
// Copyright (C) 2010-2013, Advanced Micro Devices, Inc., all rights reserved.
|
||||
// Third party copyrights are property of their respective owners.
|
||||
//
|
||||
// @Authors
|
||||
// Peng Xiao, pengxiao@multicorewareinc.com
|
||||
//
|
||||
// Redistribution and use in source and binary forms, with or without modification,
|
||||
// are permitted provided that the following conditions are met:
|
||||
//
|
||||
// * Redistribution's of source code must retain the above copyright notice,
|
||||
// this list of conditions and the following disclaimer.
|
||||
//
|
||||
// * Redistribution's in binary form must reproduce the above copyright notice,
|
||||
// this list of conditions and the following disclaimer in the documentation
|
||||
// and/or other oclMaterials provided with the distribution.
|
||||
//
|
||||
// * The name of the copyright holders may not be used to endorse or promote products
|
||||
// derived from this software without specific prior written permission.
|
||||
//
|
||||
// This software is provided by the copyright holders and contributors as is and
|
||||
// any express or implied warranties, including, but not limited to, the implied
|
||||
// warranties of merchantability and fitness for a particular purpose are disclaimed.
|
||||
// In no event shall the Intel Corporation or contributors be liable for any direct,
|
||||
// indirect, incidental, special, exemplary, or consequential damages
|
||||
// (including, but not limited to, procurement of substitute goods or services;
|
||||
// loss of use, data, or profits; or business interruption) however caused
|
||||
// and on any theory of liability, whether in contract, strict liability,
|
||||
// or tort (including negligence or otherwise) arising in any way out of
|
||||
// the use of this software, even if advised of the possibility of such damage.
|
||||
//
|
||||
//M*/
|
||||
|
||||
/////////////////////////////////////////////////////////
|
||||
//*******************************************************
|
||||
// basicretinafilter
|
||||
//////////////// _spatiotemporalLPfilter ////////////////
|
||||
//_horizontalCausalFilter_addInput
|
||||
kernel void horizontalCausalFilter_addInput(
|
||||
global const float * input,
|
||||
global float * output,
|
||||
const int cols,
|
||||
const int rows,
|
||||
const int elements_per_row,
|
||||
const int in_offset,
|
||||
const int out_offset,
|
||||
const float _tau,
|
||||
const float _a
|
||||
)
|
||||
{
|
||||
int gid = get_global_id(0);
|
||||
if(gid >= rows)
|
||||
{
|
||||
return;
|
||||
}
|
||||
|
||||
global const float * iptr =
|
||||
input + mad24(gid, elements_per_row, in_offset / 4);
|
||||
global float * optr =
|
||||
output + mad24(gid, elements_per_row, out_offset / 4);
|
||||
|
||||
float res;
|
||||
float4 in_v4, out_v4, res_v4 = (float4)(0);
|
||||
//vectorize to increase throughput
|
||||
for(int i = 0; i < cols / 4; ++i, iptr += 4, optr += 4)
|
||||
{
|
||||
in_v4 = vload4(0, iptr);
|
||||
out_v4 = vload4(0, optr);
|
||||
|
||||
res_v4.x = in_v4.x + _tau * out_v4.x + _a * res_v4.w;
|
||||
res_v4.y = in_v4.y + _tau * out_v4.y + _a * res_v4.x;
|
||||
res_v4.z = in_v4.z + _tau * out_v4.z + _a * res_v4.y;
|
||||
res_v4.w = in_v4.w + _tau * out_v4.w + _a * res_v4.z;
|
||||
|
||||
vstore4(res_v4, 0, optr);
|
||||
}
|
||||
res = res_v4.w;
|
||||
// there may be left some
|
||||
for(int i = 0; i < cols % 4; ++i, ++iptr, ++optr)
|
||||
{
|
||||
res = *iptr + _tau * *optr + _a * res;
|
||||
*optr = res;
|
||||
}
|
||||
}
|
||||
|
||||
//_horizontalAnticausalFilter
|
||||
kernel void horizontalAnticausalFilter(
|
||||
global float * output,
|
||||
const int cols,
|
||||
const int rows,
|
||||
const int elements_per_row,
|
||||
const int out_offset,
|
||||
const float _a
|
||||
)
|
||||
{
|
||||
int gid = get_global_id(0);
|
||||
if(gid >= rows)
|
||||
{
|
||||
return;
|
||||
}
|
||||
|
||||
global float * optr = output +
|
||||
mad24(gid + 1, elements_per_row, - 1 + out_offset / 4);
|
||||
|
||||
float4 result = (float4)(0), out_v4;
|
||||
// we assume elements_per_row is multple of 4
|
||||
for(int i = 0; i < elements_per_row / 4; ++i, optr -= 4)
|
||||
{
|
||||
// shift left, `offset` is type `size_t` so it cannot be negative
|
||||
out_v4 = vload4(0, optr - 3);
|
||||
|
||||
result.w = out_v4.w + _a * result.x;
|
||||
result.z = out_v4.z + _a * result.w;
|
||||
result.y = out_v4.y + _a * result.z;
|
||||
result.x = out_v4.x + _a * result.y;
|
||||
|
||||
vstore4(result, 0, optr - 3);
|
||||
}
|
||||
}
|
||||
|
||||
//_verticalCausalFilter
|
||||
kernel void verticalCausalFilter(
|
||||
global float * output,
|
||||
const int cols,
|
||||
const int rows,
|
||||
const int elements_per_row,
|
||||
const int out_offset,
|
||||
const float _a
|
||||
)
|
||||
{
|
||||
int gid = get_global_id(0);
|
||||
if(gid >= cols)
|
||||
{
|
||||
return;
|
||||
}
|
||||
|
||||
global float * optr = output + gid + out_offset / 4;
|
||||
float result = 0;
|
||||
for(int i = 0; i < rows; ++i, optr += elements_per_row)
|
||||
{
|
||||
result = *optr + _a * result;
|
||||
*optr = result;
|
||||
}
|
||||
}
|
||||
|
||||
//_verticalCausalFilter
|
||||
kernel void verticalAnticausalFilter_multGain(
|
||||
global float * output,
|
||||
const int cols,
|
||||
const int rows,
|
||||
const int elements_per_row,
|
||||
const int out_offset,
|
||||
const float _a,
|
||||
const float _gain
|
||||
)
|
||||
{
|
||||
int gid = get_global_id(0);
|
||||
if(gid >= cols)
|
||||
{
|
||||
return;
|
||||
}
|
||||
|
||||
global float * optr = output + (rows - 1) * elements_per_row + gid + out_offset / 4;
|
||||
float result = 0;
|
||||
for(int i = 0; i < rows; ++i, optr -= elements_per_row)
|
||||
{
|
||||
result = *optr + _a * result;
|
||||
*optr = _gain * result;
|
||||
}
|
||||
}
|
||||
//
|
||||
// end of _spatiotemporalLPfilter
|
||||
/////////////////////////////////////////////////////////////////////
|
||||
|
||||
//////////////// horizontalAnticausalFilter_Irregular ////////////////
|
||||
kernel void horizontalAnticausalFilter_Irregular(
|
||||
global float * output,
|
||||
global float * buffer,
|
||||
const int cols,
|
||||
const int rows,
|
||||
const int elements_per_row,
|
||||
const int out_offset,
|
||||
const int buffer_offset
|
||||
)
|
||||
{
|
||||
int gid = get_global_id(0);
|
||||
if(gid >= rows)
|
||||
{
|
||||
return;
|
||||
}
|
||||
|
||||
global float * optr =
|
||||
output + mad24(rows - gid, elements_per_row, -1 + out_offset / 4);
|
||||
global float * bptr =
|
||||
buffer + mad24(rows - gid, elements_per_row, -1 + buffer_offset / 4);
|
||||
|
||||
float4 buf_v4, out_v4, res_v4 = (float4)(0);
|
||||
|
||||
for(int i = 0; i < elements_per_row / 4; ++i, optr -= 4, bptr -= 4)
|
||||
{
|
||||
buf_v4 = vload4(0, bptr - 3);
|
||||
out_v4 = vload4(0, optr - 3);
|
||||
|
||||
res_v4.w = out_v4.w + buf_v4.w * res_v4.x;
|
||||
res_v4.z = out_v4.z + buf_v4.z * res_v4.w;
|
||||
res_v4.y = out_v4.y + buf_v4.y * res_v4.z;
|
||||
res_v4.x = out_v4.x + buf_v4.x * res_v4.y;
|
||||
|
||||
vstore4(res_v4, 0, optr - 3);
|
||||
}
|
||||
}
|
||||
|
||||
//////////////// verticalCausalFilter_Irregular ////////////////
|
||||
kernel void verticalCausalFilter_Irregular(
|
||||
global float * output,
|
||||
global float * buffer,
|
||||
const int cols,
|
||||
const int rows,
|
||||
const int elements_per_row,
|
||||
const int out_offset,
|
||||
const int buffer_offset
|
||||
)
|
||||
{
|
||||
int gid = get_global_id(0);
|
||||
if(gid >= cols)
|
||||
{
|
||||
return;
|
||||
}
|
||||
|
||||
global float * optr = output + gid + out_offset / 4;
|
||||
global float * bptr = buffer + gid + buffer_offset / 4;
|
||||
float result = 0;
|
||||
for(int i = 0; i < rows; ++i, optr += elements_per_row, bptr += elements_per_row)
|
||||
{
|
||||
result = *optr + *bptr * result;
|
||||
*optr = result;
|
||||
}
|
||||
}
|
||||
|
||||
//////////////// _adaptiveHorizontalCausalFilter_addInput ////////////////
|
||||
kernel void adaptiveHorizontalCausalFilter_addInput(
|
||||
global const float * input,
|
||||
global const float * gradient,
|
||||
global float * output,
|
||||
const int cols,
|
||||
const int rows,
|
||||
const int elements_per_row,
|
||||
const int in_offset,
|
||||
const int grad_offset,
|
||||
const int out_offset
|
||||
)
|
||||
{
|
||||
int gid = get_global_id(0);
|
||||
if(gid >= rows)
|
||||
{
|
||||
return;
|
||||
}
|
||||
|
||||
global const float * iptr =
|
||||
input + mad24(gid, elements_per_row, in_offset / 4);
|
||||
global const float * gptr =
|
||||
gradient + mad24(gid, elements_per_row, grad_offset / 4);
|
||||
global float * optr =
|
||||
output + mad24(gid, elements_per_row, out_offset / 4);
|
||||
|
||||
float4 in_v4, grad_v4, out_v4, res_v4 = (float4)(0);
|
||||
for(int i = 0; i < cols / 4; ++i, iptr += 4, gptr += 4, optr += 4)
|
||||
{
|
||||
in_v4 = vload4(0, iptr);
|
||||
grad_v4 = vload4(0, gptr);
|
||||
|
||||
res_v4.x = in_v4.x + grad_v4.x * res_v4.w;
|
||||
res_v4.y = in_v4.y + grad_v4.y * res_v4.x;
|
||||
res_v4.z = in_v4.z + grad_v4.z * res_v4.y;
|
||||
res_v4.w = in_v4.w + grad_v4.w * res_v4.z;
|
||||
|
||||
vstore4(res_v4, 0, optr);
|
||||
}
|
||||
for(int i = 0; i < cols % 4; ++i, ++iptr, ++gptr, ++optr)
|
||||
{
|
||||
res_v4.w = *iptr + *gptr * res_v4.w;
|
||||
*optr = res_v4.w;
|
||||
}
|
||||
}
|
||||
|
||||
//////////////// _adaptiveVerticalAnticausalFilter_multGain ////////////////
|
||||
kernel void adaptiveVerticalAnticausalFilter_multGain(
|
||||
global const float * gradient,
|
||||
global float * output,
|
||||
const int cols,
|
||||
const int rows,
|
||||
const int elements_per_row,
|
||||
const int grad_offset,
|
||||
const int out_offset,
|
||||
const float gain
|
||||
)
|
||||
{
|
||||
int gid = get_global_id(0);
|
||||
if(gid >= cols)
|
||||
{
|
||||
return;
|
||||
}
|
||||
|
||||
int start_idx = mad24(rows - 1, elements_per_row, gid);
|
||||
|
||||
global const float * gptr = gradient + start_idx + grad_offset / 4;
|
||||
global float * optr = output + start_idx + out_offset / 4;
|
||||
|
||||
float result = 0;
|
||||
for(int i = 0; i < rows; ++i, gptr -= elements_per_row, optr -= elements_per_row)
|
||||
{
|
||||
result = *optr + *gptr * result;
|
||||
*optr = gain * result;
|
||||
}
|
||||
}
|
||||
|
||||
//////////////// _localLuminanceAdaptation ////////////////
|
||||
// FIXME:
|
||||
// This kernel seems to have precision problem on GPU
|
||||
kernel void localLuminanceAdaptation(
|
||||
global const float * luma,
|
||||
global const float * input,
|
||||
global float * output,
|
||||
const int cols,
|
||||
const int rows,
|
||||
const int elements_per_row,
|
||||
const float _localLuminanceAddon,
|
||||
const float _localLuminanceFactor,
|
||||
const float _maxInputValue
|
||||
)
|
||||
{
|
||||
int gidx = get_global_id(0), gidy = get_global_id(1);
|
||||
if(gidx >= cols || gidy >= rows)
|
||||
{
|
||||
return;
|
||||
}
|
||||
int offset = mad24(gidy, elements_per_row, gidx);
|
||||
|
||||
float X0 = luma[offset] * _localLuminanceFactor + _localLuminanceAddon;
|
||||
float input_val = input[offset];
|
||||
// output of the following line may be different between GPU and CPU
|
||||
output[offset] = (_maxInputValue + X0) * input_val / (input_val + X0 + 0.00000000001f);
|
||||
}
|
||||
// end of basicretinafilter
|
||||
//*******************************************************
|
||||
/////////////////////////////////////////////////////////
|
||||
|
||||
|
||||
|
||||
/////////////////////////////////////////////////////////
|
||||
//******************************************************
|
||||
// magno
|
||||
// TODO: this kernel has too many buffer accesses, better to make it
|
||||
// vector read/write for fetch efficiency
|
||||
kernel void amacrineCellsComputing(
|
||||
global const float * opl_on,
|
||||
global const float * opl_off,
|
||||
global float * prev_in_on,
|
||||
global float * prev_in_off,
|
||||
global float * out_on,
|
||||
global float * out_off,
|
||||
const int cols,
|
||||
const int rows,
|
||||
const int elements_per_row,
|
||||
const float coeff
|
||||
)
|
||||
{
|
||||
int gidx = get_global_id(0), gidy = get_global_id(1);
|
||||
if(gidx >= cols || gidy >= rows)
|
||||
{
|
||||
return;
|
||||
}
|
||||
|
||||
int offset = mad24(gidy, elements_per_row, gidx);
|
||||
opl_on += offset;
|
||||
opl_off += offset;
|
||||
prev_in_on += offset;
|
||||
prev_in_off += offset;
|
||||
out_on += offset;
|
||||
out_off += offset;
|
||||
|
||||
float magnoXonPixelResult = coeff * (*out_on + *opl_on - *prev_in_on);
|
||||
*out_on = fmax(magnoXonPixelResult, 0);
|
||||
float magnoXoffPixelResult = coeff * (*out_off + *opl_off - *prev_in_off);
|
||||
*out_off = fmax(magnoXoffPixelResult, 0);
|
||||
|
||||
*prev_in_on = *opl_on;
|
||||
*prev_in_off = *opl_off;
|
||||
}
|
||||
|
||||
/////////////////////////////////////////////////////////
|
||||
//******************************************************
|
||||
// parvo
|
||||
// TODO: this kernel has too many buffer accesses, needs optimization
|
||||
kernel void OPL_OnOffWaysComputing(
|
||||
global float4 * photo_out,
|
||||
global float4 * horiz_out,
|
||||
global float4 * bipol_on,
|
||||
global float4 * bipol_off,
|
||||
global float4 * parvo_on,
|
||||
global float4 * parvo_off,
|
||||
const int cols,
|
||||
const int rows,
|
||||
const int elements_per_row
|
||||
)
|
||||
{
|
||||
int gidx = get_global_id(0), gidy = get_global_id(1);
|
||||
if(gidx * 4 >= cols || gidy >= rows)
|
||||
{
|
||||
return;
|
||||
}
|
||||
// we assume elements_per_row must be multiples of 4
|
||||
int offset = mad24(gidy, elements_per_row >> 2, gidx);
|
||||
photo_out += offset;
|
||||
horiz_out += offset;
|
||||
bipol_on += offset;
|
||||
bipol_off += offset;
|
||||
parvo_on += offset;
|
||||
parvo_off += offset;
|
||||
|
||||
float4 diff = *photo_out - *horiz_out;
|
||||
float4 isPositive;// = convert_float4(diff > (float4)(0.0f, 0.0f, 0.0f, 0.0f));
|
||||
isPositive.x = diff.x > 0.0f;
|
||||
isPositive.y = diff.y > 0.0f;
|
||||
isPositive.z = diff.z > 0.0f;
|
||||
isPositive.w = diff.w > 0.0f;
|
||||
float4 res_on = isPositive * diff;
|
||||
float4 res_off = (isPositive - (float4)(1.0f)) * diff;
|
||||
|
||||
*bipol_on = res_on;
|
||||
*parvo_on = res_on;
|
||||
|
||||
*bipol_off = res_off;
|
||||
*parvo_off = res_off;
|
||||
}
|
||||
|
||||
/////////////////////////////////////////////////////////
|
||||
//******************************************************
|
||||
// retinacolor
|
||||
inline int bayerSampleOffset(int step, int rows, int x, int y)
|
||||
{
|
||||
return mad24(y, step, x) +
|
||||
((y % 2) + (x % 2)) * rows * step;
|
||||
}
|
||||
|
||||
|
||||
/////// colorMultiplexing //////
|
||||
kernel void runColorMultiplexingBayer(
|
||||
global const float * input,
|
||||
global float * output,
|
||||
const int cols,
|
||||
const int rows,
|
||||
const int elements_per_row
|
||||
)
|
||||
{
|
||||
int gidx = get_global_id(0), gidy = get_global_id(1);
|
||||
if(gidx >= cols || gidy >= rows)
|
||||
{
|
||||
return;
|
||||
}
|
||||
|
||||
int offset = mad24(gidy, elements_per_row, gidx);
|
||||
output[offset] = input[bayerSampleOffset(elements_per_row, rows, gidx, gidy)];
|
||||
}
|
||||
|
||||
kernel void runColorDemultiplexingBayer(
|
||||
global const float * input,
|
||||
global float * output,
|
||||
const int cols,
|
||||
const int rows,
|
||||
const int elements_per_row
|
||||
)
|
||||
{
|
||||
int gidx = get_global_id(0), gidy = get_global_id(1);
|
||||
if(gidx >= cols || gidy >= rows)
|
||||
{
|
||||
return;
|
||||
}
|
||||
|
||||
int offset = mad24(gidy, elements_per_row, gidx);
|
||||
output[bayerSampleOffset(elements_per_row, rows, gidx, gidy)] = input[offset];
|
||||
}
|
||||
|
||||
kernel void demultiplexAssign(
|
||||
global const float * input,
|
||||
global float * output,
|
||||
const int cols,
|
||||
const int rows,
|
||||
const int elements_per_row
|
||||
)
|
||||
{
|
||||
int gidx = get_global_id(0), gidy = get_global_id(1);
|
||||
if(gidx >= cols || gidy >= rows)
|
||||
{
|
||||
return;
|
||||
}
|
||||
|
||||
int offset = bayerSampleOffset(elements_per_row, rows, gidx, gidy);
|
||||
output[offset] = input[offset];
|
||||
}
|
||||
|
||||
|
||||
//// normalizeGrayOutputCentredSigmoide
|
||||
kernel void normalizeGrayOutputCentredSigmoide(
|
||||
global const float * input,
|
||||
global float * output,
|
||||
const int cols,
|
||||
const int rows,
|
||||
const int elements_per_row,
|
||||
const float meanval,
|
||||
const float X0
|
||||
)
|
||||
|
||||
{
|
||||
int gidx = get_global_id(0), gidy = get_global_id(1);
|
||||
if(gidx >= cols || gidy >= rows)
|
||||
{
|
||||
return;
|
||||
}
|
||||
int offset = mad24(gidy, elements_per_row, gidx);
|
||||
|
||||
float input_val = input[offset];
|
||||
output[offset] = meanval +
|
||||
(meanval + X0) * (input_val - meanval) / (fabs(input_val - meanval) + X0);
|
||||
}
|
||||
|
||||
//// normalize by photoreceptors density
|
||||
kernel void normalizePhotoDensity(
|
||||
global const float * chroma,
|
||||
global const float * colorDensity,
|
||||
global const float * multiplex,
|
||||
global float * luma,
|
||||
global float * demultiplex,
|
||||
const int cols,
|
||||
const int rows,
|
||||
const int elements_per_row,
|
||||
const float pG
|
||||
)
|
||||
{
|
||||
const int gidx = get_global_id(0), gidy = get_global_id(1);
|
||||
if(gidx >= cols || gidy >= rows)
|
||||
{
|
||||
return;
|
||||
}
|
||||
const int offset = mad24(gidy, elements_per_row, gidx);
|
||||
int index = offset;
|
||||
|
||||
float Cr = chroma[index] * colorDensity[index];
|
||||
index += elements_per_row * rows;
|
||||
float Cg = chroma[index] * colorDensity[index];
|
||||
index += elements_per_row * rows;
|
||||
float Cb = chroma[index] * colorDensity[index];
|
||||
|
||||
const float luma_res = (Cr + Cg + Cb) * pG;
|
||||
luma[offset] = luma_res;
|
||||
demultiplex[bayerSampleOffset(elements_per_row, rows, gidx, gidy)] =
|
||||
multiplex[offset] - luma_res;
|
||||
}
|
||||
|
||||
|
||||
|
||||
//////// computeGradient ///////
|
||||
// TODO:
|
||||
// this function maybe accelerated by image2d_t or lds
|
||||
kernel void computeGradient(
|
||||
global const float * luma,
|
||||
global float * gradient,
|
||||
const int cols,
|
||||
const int rows,
|
||||
const int elements_per_row
|
||||
)
|
||||
{
|
||||
int gidx = get_global_id(0) + 2, gidy = get_global_id(1) + 2;
|
||||
if(gidx >= cols - 2 || gidy >= rows - 2)
|
||||
{
|
||||
return;
|
||||
}
|
||||
int offset = mad24(gidy, elements_per_row, gidx);
|
||||
luma += offset;
|
||||
|
||||
// horizontal and vertical local gradients
|
||||
const float v_grad = fabs(luma[elements_per_row] - luma[- elements_per_row]);
|
||||
const float h_grad = fabs(luma[1] - luma[-1]);
|
||||
|
||||
// neighborhood horizontal and vertical gradients
|
||||
const float cur_val = luma[0];
|
||||
const float v_grad_p = fabs(cur_val - luma[- 2 * elements_per_row]);
|
||||
const float h_grad_p = fabs(cur_val - luma[- 2]);
|
||||
const float v_grad_n = fabs(cur_val - luma[2 * elements_per_row]);
|
||||
const float h_grad_n = fabs(cur_val - luma[2]);
|
||||
|
||||
const float horiz_grad = 0.5f * h_grad + 0.25f * (h_grad_p + h_grad_n);
|
||||
const float verti_grad = 0.5f * v_grad + 0.25f * (v_grad_p + v_grad_n);
|
||||
const bool is_vertical_greater = horiz_grad < verti_grad;
|
||||
|
||||
gradient[offset + elements_per_row * rows] = is_vertical_greater ? 0.06f : 0.57f;
|
||||
gradient[offset ] = is_vertical_greater ? 0.57f : 0.06f;
|
||||
}
|
||||
|
||||
|
||||
/////// substractResidual ///////
|
||||
kernel void substractResidual(
|
||||
global float * input,
|
||||
const int cols,
|
||||
const int rows,
|
||||
const int elements_per_row,
|
||||
const float pR,
|
||||
const float pG,
|
||||
const float pB
|
||||
)
|
||||
{
|
||||
const int gidx = get_global_id(0), gidy = get_global_id(1);
|
||||
if(gidx >= cols || gidy >= rows)
|
||||
{
|
||||
return;
|
||||
}
|
||||
int indices [3] =
|
||||
{
|
||||
mad24(gidy, elements_per_row, gidx),
|
||||
mad24(gidy + rows, elements_per_row, gidx),
|
||||
mad24(gidy + 2 * rows, elements_per_row, gidx)
|
||||
};
|
||||
float vals[3] = {input[indices[0]], input[indices[1]], input[indices[2]]};
|
||||
float residu = pR * vals[0] + pG * vals[1] + pB * vals[2];
|
||||
|
||||
input[indices[0]] = vals[0] - residu;
|
||||
input[indices[1]] = vals[1] - residu;
|
||||
input[indices[2]] = vals[2] - residu;
|
||||
}
|
||||
|
||||
///// clipRGBOutput_0_maxInputValue /////
|
||||
kernel void clipRGBOutput_0_maxInputValue(
|
||||
global float * input,
|
||||
const int cols,
|
||||
const int rows,
|
||||
const int elements_per_row,
|
||||
const float maxVal
|
||||
)
|
||||
{
|
||||
const int gidx = get_global_id(0), gidy = get_global_id(1);
|
||||
if(gidx >= cols || gidy >= rows)
|
||||
{
|
||||
return;
|
||||
}
|
||||
const int offset = mad24(gidy, elements_per_row, gidx);
|
||||
float val = input[offset];
|
||||
val = clamp(val, 0.0f, maxVal);
|
||||
input[offset] = val;
|
||||
}
|
||||
|
||||
//// normalizeGrayOutputNearZeroCentreredSigmoide ////
|
||||
kernel void normalizeGrayOutputNearZeroCentreredSigmoide(
|
||||
global float * input,
|
||||
global float * output,
|
||||
const int cols,
|
||||
const int rows,
|
||||
const int elements_per_row,
|
||||
const float maxVal,
|
||||
const float X0cube
|
||||
)
|
||||
{
|
||||
const int gidx = get_global_id(0), gidy = get_global_id(1);
|
||||
if(gidx >= cols || gidy >= rows)
|
||||
{
|
||||
return;
|
||||
}
|
||||
const int offset = mad24(gidy, elements_per_row, gidx);
|
||||
float currentCubeLuminance = input[offset];
|
||||
currentCubeLuminance = currentCubeLuminance * currentCubeLuminance * currentCubeLuminance;
|
||||
output[offset] = currentCubeLuminance * X0cube / (X0cube + currentCubeLuminance);
|
||||
}
|
||||
|
||||
//// centerReductImageLuminance ////
|
||||
kernel void centerReductImageLuminance(
|
||||
global float * input,
|
||||
const int cols,
|
||||
const int rows,
|
||||
const int elements_per_row,
|
||||
const float mean,
|
||||
const float std_dev
|
||||
)
|
||||
{
|
||||
const int gidx = get_global_id(0), gidy = get_global_id(1);
|
||||
if(gidx >= cols || gidy >= rows)
|
||||
{
|
||||
return;
|
||||
}
|
||||
const int offset = mad24(gidy, elements_per_row, gidx);
|
||||
|
||||
float val = input[offset];
|
||||
input[offset] = (val - mean) / std_dev;
|
||||
}
|
||||
|
||||
//// inverseValue ////
|
||||
kernel void inverseValue(
|
||||
global float * input,
|
||||
const int cols,
|
||||
const int rows,
|
||||
const int elements_per_row
|
||||
)
|
||||
{
|
||||
const int gidx = get_global_id(0), gidy = get_global_id(1);
|
||||
if(gidx >= cols || gidy >= rows)
|
||||
{
|
||||
return;
|
||||
}
|
||||
const int offset = mad24(gidy, elements_per_row, gidx);
|
||||
input[offset] = 1.f / input[offset];
|
||||
}
|
||||
|
||||
#define CV_PI 3.1415926535897932384626433832795
|
||||
|
||||
//// _processRetinaParvoMagnoMapping ////
|
||||
kernel void processRetinaParvoMagnoMapping(
|
||||
global float * parvo,
|
||||
global float * magno,
|
||||
global float * output,
|
||||
const int cols,
|
||||
const int rows,
|
||||
const int halfCols,
|
||||
const int halfRows,
|
||||
const int elements_per_row,
|
||||
const float minDistance
|
||||
)
|
||||
{
|
||||
const int gidx = get_global_id(0), gidy = get_global_id(1);
|
||||
if(gidx >= cols || gidy >= rows)
|
||||
{
|
||||
return;
|
||||
}
|
||||
const int offset = mad24(gidy, elements_per_row, gidx);
|
||||
|
||||
float distanceToCenter =
|
||||
sqrt(((float)(gidy - halfRows) * (gidy - halfRows) + (gidx - halfCols) * (gidx - halfCols)));
|
||||
|
||||
float a = distanceToCenter < minDistance ?
|
||||
(0.5f + 0.5f * (float)cos(CV_PI * distanceToCenter / minDistance)) : 0;
|
||||
float b = 1.f - a;
|
||||
|
||||
output[offset] = parvo[offset] * a + magno[offset] * b;
|
||||
}
|
@ -43,11 +43,17 @@
|
||||
#ifndef __OPENCV_PRECOMP_H__
|
||||
#define __OPENCV_PRECOMP_H__
|
||||
|
||||
#include "opencv2/opencv_modules.hpp"
|
||||
#include "opencv2/bioinspired.hpp"
|
||||
#include "opencv2/core/utility.hpp"
|
||||
#include "opencv2/core/private.hpp"
|
||||
|
||||
#include <valarray>
|
||||
|
||||
#ifdef HAVE_OPENCV_OCL
|
||||
#include "opencv2/ocl/private/util.hpp"
|
||||
#endif
|
||||
|
||||
namespace cv
|
||||
{
|
||||
|
||||
|
1651
modules/bioinspired/src/retina_ocl.cpp
Normal file
1651
modules/bioinspired/src/retina_ocl.cpp
Normal file
File diff suppressed because it is too large
Load Diff
633
modules/bioinspired/src/retina_ocl.hpp
Normal file
633
modules/bioinspired/src/retina_ocl.hpp
Normal file
@ -0,0 +1,633 @@
|
||||
/*M///////////////////////////////////////////////////////////////////////////////////////
|
||||
//
|
||||
// IMPORTANT: READ BEFORE DOWNLOADING, COPYING, INSTALLING OR USING.
|
||||
//
|
||||
// By downloading, copying, installing or using the software you agree to this license.
|
||||
// If you do not agree to this license, do not download, install,
|
||||
// copy or use the software.
|
||||
//
|
||||
//
|
||||
// License Agreement
|
||||
// For Open Source Computer Vision Library
|
||||
//
|
||||
// Copyright (C) 2010-2013, Multicoreware, Inc., all rights reserved.
|
||||
// Copyright (C) 2010-2013, Advanced Micro Devices, Inc., all rights reserved.
|
||||
// Third party copyrights are property of their respective owners.
|
||||
//
|
||||
// @Authors
|
||||
// Peng Xiao, pengxiao@multicorewareinc.com
|
||||
//
|
||||
// Redistribution and use in source and binary forms, with or without modification,
|
||||
// are permitted provided that the following conditions are met:
|
||||
//
|
||||
// * Redistribution's of source code must retain the above copyright notice,
|
||||
// this list of conditions and the following disclaimer.
|
||||
//
|
||||
// * Redistribution's in binary form must reproduce the above copyright notice,
|
||||
// this list of conditions and the following disclaimer in the documentation
|
||||
// and/or other oclMaterials provided with the distribution.
|
||||
//
|
||||
// * The name of the copyright holders may not be used to endorse or promote products
|
||||
// derived from this software without specific prior written permission.
|
||||
//
|
||||
// This software is provided by the copyright holders and contributors as is and
|
||||
// any express or implied warranties, including, but not limited to, the implied
|
||||
// warranties of merchantability and fitness for a particular purpose are disclaimed.
|
||||
// In no event shall the Intel Corporation or contributors be liable for any direct,
|
||||
// indirect, incidental, special, exemplary, or consequential damages
|
||||
// (including, but not limited to, procurement of substitute goods or services;
|
||||
// loss of use, data, or profits; or business interruption) however caused
|
||||
// and on any theory of liability, whether in contract, strict liability,
|
||||
// or tort (including negligence or otherwise) arising in any way out of
|
||||
// the use of this software, even if advised of the possibility of such damage.
|
||||
//
|
||||
//M*/
|
||||
|
||||
#ifndef __OCL_RETINA_HPP__
|
||||
#define __OCL_RETINA_HPP__
|
||||
|
||||
#include "precomp.hpp"
|
||||
|
||||
#ifdef HAVE_OPENCV_OCL
|
||||
|
||||
// please refer to c++ headers for API comments
|
||||
namespace cv
|
||||
{
|
||||
namespace bioinspired
|
||||
{
|
||||
namespace ocl
|
||||
{
|
||||
void normalizeGrayOutputCentredSigmoide(const float meanValue, const float sensitivity, cv::ocl::oclMat &in, cv::ocl::oclMat &out, const float maxValue = 255.f);
|
||||
void normalizeGrayOutput_0_maxOutputValue(cv::ocl::oclMat &inputOutputBuffer, const float maxOutputValue = 255.0);
|
||||
void normalizeGrayOutputNearZeroCentreredSigmoide(cv::ocl::oclMat &inputPicture, cv::ocl::oclMat &outputBuffer, const float sensitivity = 40, const float maxOutputValue = 255.0f);
|
||||
void centerReductImageLuminance(cv::ocl::oclMat &inputOutputBuffer);
|
||||
|
||||
class BasicRetinaFilter
|
||||
{
|
||||
public:
|
||||
BasicRetinaFilter(const unsigned int NBrows, const unsigned int NBcolumns, const unsigned int parametersListSize = 1, const bool useProgressiveFilter = false);
|
||||
~BasicRetinaFilter();
|
||||
inline void clearOutputBuffer()
|
||||
{
|
||||
_filterOutput = 0;
|
||||
};
|
||||
inline void clearSecondaryBuffer()
|
||||
{
|
||||
_localBuffer = 0;
|
||||
};
|
||||
inline void clearAllBuffers()
|
||||
{
|
||||
clearOutputBuffer();
|
||||
clearSecondaryBuffer();
|
||||
};
|
||||
void resize(const unsigned int NBrows, const unsigned int NBcolumns);
|
||||
const cv::ocl::oclMat &runFilter_LPfilter(const cv::ocl::oclMat &inputFrame, const unsigned int filterIndex = 0);
|
||||
void runFilter_LPfilter(const cv::ocl::oclMat &inputFrame, cv::ocl::oclMat &outputFrame, const unsigned int filterIndex = 0);
|
||||
void runFilter_LPfilter_Autonomous(cv::ocl::oclMat &inputOutputFrame, const unsigned int filterIndex = 0);
|
||||
const cv::ocl::oclMat &runFilter_LocalAdapdation(const cv::ocl::oclMat &inputOutputFrame, const cv::ocl::oclMat &localLuminance);
|
||||
void runFilter_LocalAdapdation(const cv::ocl::oclMat &inputFrame, const cv::ocl::oclMat &localLuminance, cv::ocl::oclMat &outputFrame);
|
||||
const cv::ocl::oclMat &runFilter_LocalAdapdation_autonomous(const cv::ocl::oclMat &inputFrame);
|
||||
void runFilter_LocalAdapdation_autonomous(const cv::ocl::oclMat &inputFrame, cv::ocl::oclMat &outputFrame);
|
||||
void setLPfilterParameters(const float beta, const float tau, const float k, const unsigned int filterIndex = 0);
|
||||
inline void setV0CompressionParameter(const float v0, const float maxInputValue, const float)
|
||||
{
|
||||
_v0 = v0 * maxInputValue;
|
||||
_localLuminanceFactor = v0;
|
||||
_localLuminanceAddon = maxInputValue * (1.0f - v0);
|
||||
_maxInputValue = maxInputValue;
|
||||
};
|
||||
inline void setV0CompressionParameter(const float v0, const float meanLuminance)
|
||||
{
|
||||
this->setV0CompressionParameter(v0, _maxInputValue, meanLuminance);
|
||||
};
|
||||
inline void setV0CompressionParameter(const float v0)
|
||||
{
|
||||
_v0 = v0 * _maxInputValue;
|
||||
_localLuminanceFactor = v0;
|
||||
_localLuminanceAddon = _maxInputValue * (1.0f - v0);
|
||||
};
|
||||
inline void setV0CompressionParameterToneMapping(const float v0, const float maxInputValue, const float meanLuminance = 128.0f)
|
||||
{
|
||||
_v0 = v0 * maxInputValue;
|
||||
_localLuminanceFactor = 1.0f;
|
||||
_localLuminanceAddon = meanLuminance * _v0;
|
||||
_maxInputValue = maxInputValue;
|
||||
};
|
||||
inline void updateCompressionParameter(const float meanLuminance)
|
||||
{
|
||||
_localLuminanceFactor = 1;
|
||||
_localLuminanceAddon = meanLuminance * _v0;
|
||||
};
|
||||
inline float getV0CompressionParameter()
|
||||
{
|
||||
return _v0 / _maxInputValue;
|
||||
};
|
||||
inline const cv::ocl::oclMat &getOutput() const
|
||||
{
|
||||
return _filterOutput;
|
||||
};
|
||||
inline unsigned int getNBrows()
|
||||
{
|
||||
return _filterOutput.rows;
|
||||
};
|
||||
inline unsigned int getNBcolumns()
|
||||
{
|
||||
return _filterOutput.cols;
|
||||
};
|
||||
inline unsigned int getNBpixels()
|
||||
{
|
||||
return _filterOutput.size().area();
|
||||
};
|
||||
inline void normalizeGrayOutput_0_maxOutputValue(const float maxValue)
|
||||
{
|
||||
ocl::normalizeGrayOutput_0_maxOutputValue(_filterOutput, maxValue);
|
||||
};
|
||||
inline void normalizeGrayOutputCentredSigmoide()
|
||||
{
|
||||
ocl::normalizeGrayOutputCentredSigmoide(0.0, 2.0, _filterOutput, _filterOutput);
|
||||
};
|
||||
inline void centerReductImageLuminance()
|
||||
{
|
||||
ocl::centerReductImageLuminance(_filterOutput);
|
||||
};
|
||||
inline float getMaxInputValue()
|
||||
{
|
||||
return this->_maxInputValue;
|
||||
};
|
||||
inline void setMaxInputValue(const float newMaxInputValue)
|
||||
{
|
||||
this->_maxInputValue = newMaxInputValue;
|
||||
};
|
||||
|
||||
protected:
|
||||
cv::ocl::oclMat _filterOutput;
|
||||
cv::ocl::oclMat _localBuffer;
|
||||
|
||||
int _NBrows;
|
||||
int _NBcols;
|
||||
unsigned int _halfNBrows;
|
||||
unsigned int _halfNBcolumns;
|
||||
|
||||
std::valarray <float>_filteringCoeficientsTable;
|
||||
float _v0;
|
||||
float _maxInputValue;
|
||||
float _meanInputValue;
|
||||
float _localLuminanceFactor;
|
||||
float _localLuminanceAddon;
|
||||
|
||||
float _a;
|
||||
float _tau;
|
||||
float _gain;
|
||||
|
||||
void _spatiotemporalLPfilter(const cv::ocl::oclMat &inputFrame, cv::ocl::oclMat &LPfilterOutput, const unsigned int coefTableOffset = 0);
|
||||
float _squaringSpatiotemporalLPfilter(const cv::ocl::oclMat &inputFrame, cv::ocl::oclMat &outputFrame, const unsigned int filterIndex = 0);
|
||||
void _spatiotemporalLPfilter_Irregular(const cv::ocl::oclMat &inputFrame, cv::ocl::oclMat &outputFrame, const unsigned int filterIndex = 0);
|
||||
void _localSquaringSpatioTemporalLPfilter(const cv::ocl::oclMat &inputFrame, cv::ocl::oclMat &LPfilterOutput, const unsigned int *integrationAreas, const unsigned int filterIndex = 0);
|
||||
void _localLuminanceAdaptation(const cv::ocl::oclMat &inputFrame, const cv::ocl::oclMat &localLuminance, cv::ocl::oclMat &outputFrame, const bool updateLuminanceMean = true);
|
||||
void _localLuminanceAdaptation(cv::ocl::oclMat &inputOutputFrame, const cv::ocl::oclMat &localLuminance);
|
||||
void _localLuminanceAdaptationPosNegValues(const cv::ocl::oclMat &inputFrame, const cv::ocl::oclMat &localLuminance, float *outputFrame);
|
||||
void _horizontalCausalFilter_addInput(const cv::ocl::oclMat &inputFrame, cv::ocl::oclMat &outputFrame);
|
||||
void _horizontalAnticausalFilter(cv::ocl::oclMat &outputFrame);
|
||||
void _verticalCausalFilter(cv::ocl::oclMat &outputFrame);
|
||||
void _horizontalAnticausalFilter_Irregular(cv::ocl::oclMat &outputFrame, const cv::ocl::oclMat &spatialConstantBuffer);
|
||||
void _verticalCausalFilter_Irregular(cv::ocl::oclMat &outputFrame, const cv::ocl::oclMat &spatialConstantBuffer);
|
||||
void _verticalAnticausalFilter_multGain(cv::ocl::oclMat &outputFrame);
|
||||
};
|
||||
|
||||
class MagnoRetinaFilter: public BasicRetinaFilter
|
||||
{
|
||||
public:
|
||||
MagnoRetinaFilter(const unsigned int NBrows, const unsigned int NBcolumns);
|
||||
virtual ~MagnoRetinaFilter();
|
||||
void clearAllBuffers();
|
||||
void resize(const unsigned int NBrows, const unsigned int NBcolumns);
|
||||
void setCoefficientsTable(const float parasolCells_beta, const float parasolCells_tau, const float parasolCells_k, const float amacrinCellsTemporalCutFrequency, const float localAdaptIntegration_tau, const float localAdaptIntegration_k);
|
||||
|
||||
const cv::ocl::oclMat &runFilter(const cv::ocl::oclMat &OPL_ON, const cv::ocl::oclMat &OPL_OFF);
|
||||
|
||||
inline const cv::ocl::oclMat &getMagnoON() const
|
||||
{
|
||||
return _magnoXOutputON;
|
||||
};
|
||||
inline const cv::ocl::oclMat &getMagnoOFF() const
|
||||
{
|
||||
return _magnoXOutputOFF;
|
||||
};
|
||||
inline const cv::ocl::oclMat &getMagnoYsaturated() const
|
||||
{
|
||||
return _magnoYsaturated;
|
||||
};
|
||||
inline void normalizeGrayOutputNearZeroCentreredSigmoide()
|
||||
{
|
||||
ocl::normalizeGrayOutputNearZeroCentreredSigmoide(_magnoYOutput, _magnoYsaturated);
|
||||
};
|
||||
inline float getTemporalConstant()
|
||||
{
|
||||
return this->_filteringCoeficientsTable[2];
|
||||
};
|
||||
private:
|
||||
cv::ocl::oclMat _previousInput_ON;
|
||||
cv::ocl::oclMat _previousInput_OFF;
|
||||
cv::ocl::oclMat _amacrinCellsTempOutput_ON;
|
||||
cv::ocl::oclMat _amacrinCellsTempOutput_OFF;
|
||||
cv::ocl::oclMat _magnoXOutputON;
|
||||
cv::ocl::oclMat _magnoXOutputOFF;
|
||||
cv::ocl::oclMat _localProcessBufferON;
|
||||
cv::ocl::oclMat _localProcessBufferOFF;
|
||||
cv::ocl::oclMat _magnoYOutput;
|
||||
cv::ocl::oclMat _magnoYsaturated;
|
||||
|
||||
float _temporalCoefficient;
|
||||
void _amacrineCellsComputing(const cv::ocl::oclMat &OPL_ON, const cv::ocl::oclMat &OPL_OFF);
|
||||
};
|
||||
|
||||
class ParvoRetinaFilter: public BasicRetinaFilter
|
||||
{
|
||||
public:
|
||||
ParvoRetinaFilter(const unsigned int NBrows = 480, const unsigned int NBcolumns = 640);
|
||||
virtual ~ParvoRetinaFilter();
|
||||
void resize(const unsigned int NBrows, const unsigned int NBcolumns);
|
||||
void clearAllBuffers();
|
||||
void setOPLandParvoFiltersParameters(const float beta1, const float tau1, const float k1, const float beta2, const float tau2, const float k2);
|
||||
|
||||
inline void setGanglionCellsLocalAdaptationLPfilterParameters(const float tau, const float k)
|
||||
{
|
||||
BasicRetinaFilter::setLPfilterParameters(0, tau, k, 2);
|
||||
};
|
||||
const cv::ocl::oclMat &runFilter(const cv::ocl::oclMat &inputFrame, const bool useParvoOutput = true);
|
||||
|
||||
inline const cv::ocl::oclMat &getPhotoreceptorsLPfilteringOutput() const
|
||||
{
|
||||
return _photoreceptorsOutput;
|
||||
};
|
||||
|
||||
inline const cv::ocl::oclMat &getHorizontalCellsOutput() const
|
||||
{
|
||||
return _horizontalCellsOutput;
|
||||
};
|
||||
|
||||
inline const cv::ocl::oclMat &getParvoON() const
|
||||
{
|
||||
return _parvocellularOutputON;
|
||||
};
|
||||
|
||||
inline const cv::ocl::oclMat &getParvoOFF() const
|
||||
{
|
||||
return _parvocellularOutputOFF;
|
||||
};
|
||||
|
||||
inline const cv::ocl::oclMat &getBipolarCellsON() const
|
||||
{
|
||||
return _bipolarCellsOutputON;
|
||||
};
|
||||
|
||||
inline const cv::ocl::oclMat &getBipolarCellsOFF() const
|
||||
{
|
||||
return _bipolarCellsOutputOFF;
|
||||
};
|
||||
|
||||
inline float getPhotoreceptorsTemporalConstant()
|
||||
{
|
||||
return this->_filteringCoeficientsTable[2];
|
||||
};
|
||||
|
||||
inline float getHcellsTemporalConstant()
|
||||
{
|
||||
return this->_filteringCoeficientsTable[5];
|
||||
};
|
||||
private:
|
||||
cv::ocl::oclMat _photoreceptorsOutput;
|
||||
cv::ocl::oclMat _horizontalCellsOutput;
|
||||
cv::ocl::oclMat _parvocellularOutputON;
|
||||
cv::ocl::oclMat _parvocellularOutputOFF;
|
||||
cv::ocl::oclMat _bipolarCellsOutputON;
|
||||
cv::ocl::oclMat _bipolarCellsOutputOFF;
|
||||
cv::ocl::oclMat _localAdaptationOFF;
|
||||
cv::ocl::oclMat _localAdaptationON;
|
||||
cv::ocl::oclMat _parvocellularOutputONminusOFF;
|
||||
void _OPL_OnOffWaysComputing();
|
||||
};
|
||||
class RetinaColor: public BasicRetinaFilter
|
||||
{
|
||||
public:
|
||||
RetinaColor(const unsigned int NBrows, const unsigned int NBcolumns, const int samplingMethod = RETINA_COLOR_DIAGONAL);
|
||||
virtual ~RetinaColor();
|
||||
|
||||
void clearAllBuffers();
|
||||
void resize(const unsigned int NBrows, const unsigned int NBcolumns);
|
||||
inline void runColorMultiplexing(const cv::ocl::oclMat &inputRGBFrame)
|
||||
{
|
||||
runColorMultiplexing(inputRGBFrame, _multiplexedFrame);
|
||||
};
|
||||
void runColorMultiplexing(const cv::ocl::oclMat &demultiplexedInputFrame, cv::ocl::oclMat &multiplexedFrame);
|
||||
void runColorDemultiplexing(const cv::ocl::oclMat &multiplexedColorFrame, const bool adaptiveFiltering = false, const float maxInputValue = 255.0);
|
||||
|
||||
void setColorSaturation(const bool saturateColors = true, const float colorSaturationValue = 4.0)
|
||||
{
|
||||
_saturateColors = saturateColors;
|
||||
_colorSaturationValue = colorSaturationValue;
|
||||
};
|
||||
|
||||
void setChrominanceLPfilterParameters(const float beta, const float tau, const float k)
|
||||
{
|
||||
setLPfilterParameters(beta, tau, k);
|
||||
};
|
||||
|
||||
bool applyKrauskopfLMS2Acr1cr2Transform(cv::ocl::oclMat &result);
|
||||
bool applyLMS2LabTransform(cv::ocl::oclMat &result);
|
||||
inline const cv::ocl::oclMat &getMultiplexedFrame() const
|
||||
{
|
||||
return _multiplexedFrame;
|
||||
};
|
||||
|
||||
inline const cv::ocl::oclMat &getDemultiplexedColorFrame() const
|
||||
{
|
||||
return _demultiplexedColorFrame;
|
||||
};
|
||||
|
||||
inline const cv::ocl::oclMat &getLuminance() const
|
||||
{
|
||||
return _luminance;
|
||||
};
|
||||
inline const cv::ocl::oclMat &getChrominance() const
|
||||
{
|
||||
return _chrominance;
|
||||
};
|
||||
void clipRGBOutput_0_maxInputValue(cv::ocl::oclMat &inputOutputBuffer, const float maxOutputValue = 255.0);
|
||||
void normalizeRGBOutput_0_maxOutputValue(const float maxOutputValue = 255.0);
|
||||
inline void setDemultiplexedColorFrame(const cv::ocl::oclMat &demultiplexedImage)
|
||||
{
|
||||
_demultiplexedColorFrame = demultiplexedImage;
|
||||
};
|
||||
protected:
|
||||
inline unsigned int bayerSampleOffset(unsigned int index)
|
||||
{
|
||||
return index + ((index / getNBcolumns()) % 2) * getNBpixels() + ((index % getNBcolumns()) % 2) * getNBpixels();
|
||||
}
|
||||
inline Rect getROI(int idx)
|
||||
{
|
||||
return Rect(0, idx * _NBrows, _NBcols, _NBrows);
|
||||
}
|
||||
int _samplingMethod;
|
||||
bool _saturateColors;
|
||||
float _colorSaturationValue;
|
||||
cv::ocl::oclMat _luminance;
|
||||
cv::ocl::oclMat _multiplexedFrame;
|
||||
cv::ocl::oclMat _RGBmosaic;
|
||||
cv::ocl::oclMat _tempMultiplexedFrame;
|
||||
cv::ocl::oclMat _demultiplexedTempBuffer;
|
||||
cv::ocl::oclMat _demultiplexedColorFrame;
|
||||
cv::ocl::oclMat _chrominance;
|
||||
cv::ocl::oclMat _colorLocalDensity;
|
||||
cv::ocl::oclMat _imageGradient;
|
||||
|
||||
float _pR, _pG, _pB;
|
||||
bool _objectInit;
|
||||
|
||||
void _initColorSampling();
|
||||
void _adaptiveSpatialLPfilter(const cv::ocl::oclMat &inputFrame, const cv::ocl::oclMat &gradient, cv::ocl::oclMat &outputFrame);
|
||||
void _adaptiveHorizontalCausalFilter_addInput(const cv::ocl::oclMat &inputFrame, const cv::ocl::oclMat &gradient, cv::ocl::oclMat &outputFrame);
|
||||
void _adaptiveVerticalAnticausalFilter_multGain(const cv::ocl::oclMat &gradient, cv::ocl::oclMat &outputFrame);
|
||||
void _computeGradient(const cv::ocl::oclMat &luminance, cv::ocl::oclMat &gradient);
|
||||
void _normalizeOutputs_0_maxOutputValue(void);
|
||||
void _applyImageColorSpaceConversion(const cv::ocl::oclMat &inputFrame, cv::ocl::oclMat &outputFrame, const float *transformTable);
|
||||
};
|
||||
class RetinaFilter
|
||||
{
|
||||
public:
|
||||
RetinaFilter(const unsigned int sizeRows, const unsigned int sizeColumns, const bool colorMode = false, const int samplingMethod = RETINA_COLOR_BAYER, const bool useRetinaLogSampling = false, const double reductionFactor = 1.0, const double samplingStrenght = 10.0);
|
||||
~RetinaFilter();
|
||||
|
||||
void clearAllBuffers();
|
||||
void resize(const unsigned int NBrows, const unsigned int NBcolumns);
|
||||
bool checkInput(const cv::ocl::oclMat &input, const bool colorMode);
|
||||
bool runFilter(const cv::ocl::oclMat &imageInput, const bool useAdaptiveFiltering = true, const bool processRetinaParvoMagnoMapping = false, const bool useColorMode = false, const bool inputIsColorMultiplexed = false);
|
||||
|
||||
void setGlobalParameters(const float OPLspatialResponse1 = 0.7, const float OPLtemporalresponse1 = 1, const float OPLassymetryGain = 0, const float OPLspatialResponse2 = 5, const float OPLtemporalresponse2 = 1, const float LPfilterSpatialResponse = 5, const float LPfilterGain = 0, const float LPfilterTemporalresponse = 0, const float MovingContoursExtractorCoefficient = 5, const bool normalizeParvoOutput_0_maxOutputValue = false, const bool normalizeMagnoOutput_0_maxOutputValue = false, const float maxOutputValue = 255.0, const float maxInputValue = 255.0, const float meanValue = 128.0);
|
||||
|
||||
inline void setPhotoreceptorsLocalAdaptationSensitivity(const float V0CompressionParameter)
|
||||
{
|
||||
_photoreceptorsPrefilter.setV0CompressionParameter(1 - V0CompressionParameter);
|
||||
_setInitPeriodCount();
|
||||
};
|
||||
|
||||
inline void setParvoGanglionCellsLocalAdaptationSensitivity(const float V0CompressionParameter)
|
||||
{
|
||||
_ParvoRetinaFilter.setV0CompressionParameter(V0CompressionParameter);
|
||||
_setInitPeriodCount();
|
||||
};
|
||||
|
||||
inline void setGanglionCellsLocalAdaptationLPfilterParameters(const float spatialResponse, const float temporalResponse)
|
||||
{
|
||||
_ParvoRetinaFilter.setGanglionCellsLocalAdaptationLPfilterParameters(temporalResponse, spatialResponse);
|
||||
_setInitPeriodCount();
|
||||
};
|
||||
|
||||
inline void setMagnoGanglionCellsLocalAdaptationSensitivity(const float V0CompressionParameter)
|
||||
{
|
||||
_MagnoRetinaFilter.setV0CompressionParameter(V0CompressionParameter);
|
||||
_setInitPeriodCount();
|
||||
};
|
||||
|
||||
void setOPLandParvoParameters(const float beta1, const float tau1, const float k1, const float beta2, const float tau2, const float k2, const float V0CompressionParameter)
|
||||
{
|
||||
_ParvoRetinaFilter.setOPLandParvoFiltersParameters(beta1, tau1, k1, beta2, tau2, k2);
|
||||
_ParvoRetinaFilter.setV0CompressionParameter(V0CompressionParameter);
|
||||
_setInitPeriodCount();
|
||||
};
|
||||
|
||||
void setMagnoCoefficientsTable(const float parasolCells_beta, const float parasolCells_tau, const float parasolCells_k, const float amacrinCellsTemporalCutFrequency, const float V0CompressionParameter, const float localAdaptintegration_tau, const float localAdaptintegration_k)
|
||||
{
|
||||
_MagnoRetinaFilter.setCoefficientsTable(parasolCells_beta, parasolCells_tau, parasolCells_k, amacrinCellsTemporalCutFrequency, localAdaptintegration_tau, localAdaptintegration_k);
|
||||
_MagnoRetinaFilter.setV0CompressionParameter(V0CompressionParameter);
|
||||
_setInitPeriodCount();
|
||||
};
|
||||
|
||||
inline void activateNormalizeParvoOutput_0_maxOutputValue(const bool normalizeParvoOutput_0_maxOutputValue)
|
||||
{
|
||||
_normalizeParvoOutput_0_maxOutputValue = normalizeParvoOutput_0_maxOutputValue;
|
||||
};
|
||||
|
||||
inline void activateNormalizeMagnoOutput_0_maxOutputValue(const bool normalizeMagnoOutput_0_maxOutputValue)
|
||||
{
|
||||
_normalizeMagnoOutput_0_maxOutputValue = normalizeMagnoOutput_0_maxOutputValue;
|
||||
};
|
||||
|
||||
inline void setMaxOutputValue(const float maxOutputValue)
|
||||
{
|
||||
_maxOutputValue = maxOutputValue;
|
||||
};
|
||||
|
||||
void setColorMode(const bool desiredColorMode)
|
||||
{
|
||||
_useColorMode = desiredColorMode;
|
||||
};
|
||||
inline void setColorSaturation(const bool saturateColors = true, const float colorSaturationValue = 4.0)
|
||||
{
|
||||
_colorEngine.setColorSaturation(saturateColors, colorSaturationValue);
|
||||
};
|
||||
inline const cv::ocl::oclMat &getLocalAdaptation() const
|
||||
{
|
||||
return _photoreceptorsPrefilter.getOutput();
|
||||
};
|
||||
inline const cv::ocl::oclMat &getPhotoreceptors() const
|
||||
{
|
||||
return _ParvoRetinaFilter.getPhotoreceptorsLPfilteringOutput();
|
||||
};
|
||||
|
||||
inline const cv::ocl::oclMat &getHorizontalCells() const
|
||||
{
|
||||
return _ParvoRetinaFilter.getHorizontalCellsOutput();
|
||||
};
|
||||
inline bool areContoursProcessed()
|
||||
{
|
||||
return _useParvoOutput;
|
||||
};
|
||||
bool getParvoFoveaResponse(cv::ocl::oclMat &parvoFovealResponse);
|
||||
inline void activateContoursProcessing(const bool useParvoOutput)
|
||||
{
|
||||
_useParvoOutput = useParvoOutput;
|
||||
};
|
||||
|
||||
const cv::ocl::oclMat &getContours();
|
||||
|
||||
inline const cv::ocl::oclMat &getContoursON() const
|
||||
{
|
||||
return _ParvoRetinaFilter.getParvoON();
|
||||
};
|
||||
|
||||
inline const cv::ocl::oclMat &getContoursOFF() const
|
||||
{
|
||||
return _ParvoRetinaFilter.getParvoOFF();
|
||||
};
|
||||
|
||||
inline bool areMovingContoursProcessed()
|
||||
{
|
||||
return _useMagnoOutput;
|
||||
};
|
||||
|
||||
inline void activateMovingContoursProcessing(const bool useMagnoOutput)
|
||||
{
|
||||
_useMagnoOutput = useMagnoOutput;
|
||||
};
|
||||
|
||||
inline const cv::ocl::oclMat &getMovingContours() const
|
||||
{
|
||||
return _MagnoRetinaFilter.getOutput();
|
||||
};
|
||||
|
||||
inline const cv::ocl::oclMat &getMovingContoursSaturated() const
|
||||
{
|
||||
return _MagnoRetinaFilter.getMagnoYsaturated();
|
||||
};
|
||||
|
||||
inline const cv::ocl::oclMat &getMovingContoursON() const
|
||||
{
|
||||
return _MagnoRetinaFilter.getMagnoON();
|
||||
};
|
||||
|
||||
inline const cv::ocl::oclMat &getMovingContoursOFF() const
|
||||
{
|
||||
return _MagnoRetinaFilter.getMagnoOFF();
|
||||
};
|
||||
|
||||
inline const cv::ocl::oclMat &getRetinaParvoMagnoMappedOutput() const
|
||||
{
|
||||
return _retinaParvoMagnoMappedFrame;
|
||||
};
|
||||
|
||||
inline const cv::ocl::oclMat &getParvoContoursChannel() const
|
||||
{
|
||||
return _colorEngine.getLuminance();
|
||||
};
|
||||
|
||||
inline const cv::ocl::oclMat &getParvoChrominance() const
|
||||
{
|
||||
return _colorEngine.getChrominance();
|
||||
};
|
||||
inline const cv::ocl::oclMat &getColorOutput() const
|
||||
{
|
||||
return _colorEngine.getDemultiplexedColorFrame();
|
||||
};
|
||||
|
||||
inline bool isColorMode()
|
||||
{
|
||||
return _useColorMode;
|
||||
};
|
||||
bool getColorMode()
|
||||
{
|
||||
return _useColorMode;
|
||||
};
|
||||
|
||||
inline bool isInitTransitionDone()
|
||||
{
|
||||
if (_ellapsedFramesSinceLastReset < _globalTemporalConstant)
|
||||
{
|
||||
return false;
|
||||
}
|
||||
return true;
|
||||
};
|
||||
inline float getRetinaSamplingBackProjection(const float projectedRadiusLength)
|
||||
{
|
||||
return projectedRadiusLength;
|
||||
};
|
||||
|
||||
inline unsigned int getInputNBrows()
|
||||
{
|
||||
return _photoreceptorsPrefilter.getNBrows();
|
||||
};
|
||||
|
||||
inline unsigned int getInputNBcolumns()
|
||||
{
|
||||
return _photoreceptorsPrefilter.getNBcolumns();
|
||||
};
|
||||
|
||||
inline unsigned int getInputNBpixels()
|
||||
{
|
||||
return _photoreceptorsPrefilter.getNBpixels();
|
||||
};
|
||||
|
||||
inline unsigned int getOutputNBrows()
|
||||
{
|
||||
return _photoreceptorsPrefilter.getNBrows();
|
||||
};
|
||||
|
||||
inline unsigned int getOutputNBcolumns()
|
||||
{
|
||||
return _photoreceptorsPrefilter.getNBcolumns();
|
||||
};
|
||||
|
||||
inline unsigned int getOutputNBpixels()
|
||||
{
|
||||
return _photoreceptorsPrefilter.getNBpixels();
|
||||
};
|
||||
private:
|
||||
bool _useParvoOutput;
|
||||
bool _useMagnoOutput;
|
||||
|
||||
unsigned int _ellapsedFramesSinceLastReset;
|
||||
unsigned int _globalTemporalConstant;
|
||||
|
||||
cv::ocl::oclMat _retinaParvoMagnoMappedFrame;
|
||||
BasicRetinaFilter _photoreceptorsPrefilter;
|
||||
ParvoRetinaFilter _ParvoRetinaFilter;
|
||||
MagnoRetinaFilter _MagnoRetinaFilter;
|
||||
RetinaColor _colorEngine;
|
||||
|
||||
bool _useMinimalMemoryForToneMappingONLY;
|
||||
bool _normalizeParvoOutput_0_maxOutputValue;
|
||||
bool _normalizeMagnoOutput_0_maxOutputValue;
|
||||
float _maxOutputValue;
|
||||
bool _useColorMode;
|
||||
|
||||
void _setInitPeriodCount();
|
||||
void _processRetinaParvoMagnoMapping();
|
||||
void _runGrayToneMapping(const cv::ocl::oclMat &grayImageInput, cv::ocl::oclMat &grayImageOutput , const float PhotoreceptorsCompression = 0.6, const float ganglionCellsCompression = 0.6);
|
||||
};
|
||||
|
||||
} /* namespace ocl */
|
||||
} /* namespace bioinspired */
|
||||
} /* namespace cv */
|
||||
|
||||
#endif /* HAVE_OPENCV_OCL */
|
||||
#endif /* __OCL_RETINA_HPP__ */
|
139
modules/bioinspired/test/test_retina_ocl.cpp
Normal file
139
modules/bioinspired/test/test_retina_ocl.cpp
Normal file
@ -0,0 +1,139 @@
|
||||
/*M///////////////////////////////////////////////////////////////////////////////////////
|
||||
//
|
||||
// IMPORTANT: READ BEFORE DOWNLOADING, COPYING, INSTALLING OR USING.
|
||||
//
|
||||
// By downloading, copying, installing or using the software you agree to this license.
|
||||
// If you do not agree to this license, do not download, install,
|
||||
// copy or use the software.
|
||||
//
|
||||
//
|
||||
// License Agreement
|
||||
// For Open Source Computer Vision Library
|
||||
//
|
||||
// Copyright (C) 2010-2013, Multicoreware, Inc., all rights reserved.
|
||||
// Copyright (C) 2010-2013, Advanced Micro Devices, Inc., all rights reserved.
|
||||
// Third party copyrights are property of their respective owners.
|
||||
//
|
||||
// @Authors
|
||||
// Peng Xiao, pengxiao@multicorewareinc.com
|
||||
//
|
||||
// Redistribution and use in source and binary forms, with or without modification,
|
||||
// are permitted provided that the following conditions are met:
|
||||
//
|
||||
// * Redistribution's of source code must retain the above copyright notice,
|
||||
// this list of conditions and the following disclaimer.
|
||||
//
|
||||
// * Redistribution's in binary form must reproduce the above copyright notice,
|
||||
// this list of conditions and the following disclaimer in the documentation
|
||||
// and/or other oclMaterials provided with the distribution.
|
||||
//
|
||||
// * The name of the copyright holders may not be used to endorse or promote products
|
||||
// derived from this software without specific prior written permission.
|
||||
//
|
||||
// This software is provided by the copyright holders and contributors as is and
|
||||
// any express or implied warranties, including, but not limited to, the implied
|
||||
// warranties of merchantability and fitness for a particular purpose are disclaimed.
|
||||
// In no event shall the Intel Corporation or contributors be liable for any direct,
|
||||
// indirect, incidental, special, exemplary, or consequential damages
|
||||
// (including, but not limited to, procurement of substitute goods or services;
|
||||
// loss of use, data, or profits; or business interruption) however caused
|
||||
// and on any theory of liability, whether in contract, strict liability,
|
||||
// or tort (including negligence or otherwise) arising in any way out of
|
||||
// the use of this software, even if advised of the possibility of such damage.
|
||||
//
|
||||
//M*/
|
||||
|
||||
#include "test_precomp.hpp"
|
||||
#include "opencv2/opencv_modules.hpp"
|
||||
#include "opencv2/bioinspired.hpp"
|
||||
#include "opencv2/imgproc.hpp"
|
||||
#include "opencv2/highgui.hpp"
|
||||
#include "opencv2/ocl.hpp"
|
||||
|
||||
#if defined(HAVE_OPENCV_OCL) && defined(HAVE_OPENCL)
|
||||
|
||||
#define RETINA_ITERATIONS 5
|
||||
|
||||
static double checkNear(const cv::Mat &m1, const cv::Mat &m2)
|
||||
{
|
||||
return cv::norm(m1, m2, cv::NORM_INF);
|
||||
}
|
||||
|
||||
#define PARAM_TEST_CASE(name, ...) struct name : testing::TestWithParam< std::tr1::tuple< __VA_ARGS__ > >
|
||||
#define GET_PARAM(k) std::tr1::get< k >(GetParam())
|
||||
|
||||
|
||||
PARAM_TEST_CASE(Retina_OCL, bool, int, bool, double, double)
|
||||
{
|
||||
bool colorMode;
|
||||
int colorSamplingMethod;
|
||||
bool useLogSampling;
|
||||
double reductionFactor;
|
||||
double samplingStrength;
|
||||
|
||||
std::vector<cv::ocl::Info> infos;
|
||||
|
||||
virtual void SetUp()
|
||||
{
|
||||
colorMode = GET_PARAM(0);
|
||||
colorSamplingMethod = GET_PARAM(1);
|
||||
useLogSampling = GET_PARAM(2);
|
||||
reductionFactor = GET_PARAM(3);
|
||||
samplingStrength = GET_PARAM(4);
|
||||
|
||||
cv::ocl::getDevice(infos);
|
||||
std::cout << "Device name:" << infos[0].DeviceName[0] << std::endl;
|
||||
}
|
||||
};
|
||||
|
||||
TEST_P(Retina_OCL, Accuracy)
|
||||
{
|
||||
using namespace cv;
|
||||
Mat input = imread(cvtest::TS::ptr()->get_data_path() + "shared/lena.png", colorMode);
|
||||
CV_Assert(!input.empty());
|
||||
ocl::oclMat ocl_input(input);
|
||||
|
||||
Ptr<bioinspired::Retina> ocl_retina = bioinspired::createRetina_OCL(
|
||||
input.size(),
|
||||
colorMode,
|
||||
colorSamplingMethod,
|
||||
useLogSampling,
|
||||
reductionFactor,
|
||||
samplingStrength);
|
||||
|
||||
Ptr<bioinspired::Retina> gold_retina = bioinspired::createRetina(
|
||||
input.size(),
|
||||
colorMode,
|
||||
colorSamplingMethod,
|
||||
useLogSampling,
|
||||
reductionFactor,
|
||||
samplingStrength);
|
||||
|
||||
Mat gold_parvo;
|
||||
Mat gold_magno;
|
||||
ocl::oclMat ocl_parvo;
|
||||
ocl::oclMat ocl_magno;
|
||||
|
||||
for(int i = 0; i < RETINA_ITERATIONS; i ++)
|
||||
{
|
||||
ocl_retina->run(ocl_input);
|
||||
gold_retina->run(input);
|
||||
|
||||
gold_retina->getParvo(gold_parvo);
|
||||
gold_retina->getMagno(gold_magno);
|
||||
|
||||
ocl_retina->getParvo(ocl_parvo);
|
||||
ocl_retina->getMagno(ocl_magno);
|
||||
|
||||
EXPECT_LE(checkNear(gold_parvo, (Mat)ocl_parvo), 1.0);
|
||||
EXPECT_LE(checkNear(gold_magno, (Mat)ocl_magno), 1.0);
|
||||
}
|
||||
}
|
||||
|
||||
INSTANTIATE_TEST_CASE_P(Contrib, Retina_OCL, testing::Combine(
|
||||
testing::Values(false, true),
|
||||
testing::Values((int)cv::bioinspired::RETINA_COLOR_BAYER),
|
||||
testing::Values(false/*,true*/),
|
||||
testing::Values(1.0, 0.5),
|
||||
testing::Values(10.0, 5.0)));
|
||||
#endif
|
Loading…
Reference in New Issue
Block a user