iRoCS Toolbox  1.1.0
EdgeFilter-inl.hh
Go to the documentation of this file.
1 /**************************************************************************
2  *
3  * Copyright (C) 2015 Margret Keuper, Thorsten Falk
4  *
5  * Image Analysis Lab, University of Freiburg, Germany
6  *
7  * This program is free software; you can redistribute it and/or modify
8  * it under the terms of the GNU General Public License as published by
9  * the Free Software Foundation; either version 3 of the License, or
10  * (at your option) any later version.
11  *
12  * This program is distributed in the hope that it will be useful,
13  * but WITHOUT ANY WARRANTY; without even the implied warranty of
14  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15  * GNU General Public License for more details.
16  *
17  * You should have received a copy of the GNU General Public License
18  * along with this program; if not, write to the Free Software Foundation,
19  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA
20  *
21  **************************************************************************/
22 
23 #ifndef LIBMARGRET_SRC_EDGEFILTER_INL_HH
24 #define LIBMARGRET_SRC_EDGEFILTER_INL_HH
25 
26 #ifdef HAVE_CONFIG_H
27 #include <config.hh>
28 #endif
29 
30 #include "EdgeFilter.hh"
31 
34 #include <omp.h>
35 #include <list>
36 
37 namespace segmentation
38 {
39 
40  template<class DataT>
41  void edgeFilter(
42  blitz::Array<DataT,3> const &data,
43  blitz::TinyVector<double,3> const &elSize,
44  blitz::Array<DataT,3> &result,
45  blitz::TinyVector<double,3> const &lbUm,
46  blitz::TinyVector<double,3> const &ubUm,
47  iRoCS::ProgressReporter *progress)
48  {
49  // Cut out ROI
50  blitz::TinyVector<ptrdiff_t,3> lb(lbUm / elSize);
51  blitz::TinyVector<ptrdiff_t,3> ub(ubUm / elSize);
52  for (int d = 0; d < 3; ++d)
53  {
54  if (lb(d) < 0) lb(d) = 0;
55  if (ub(d) < 0) ub(d) = 0;
56  if (lb(d) >= data.extent(d)) lb(d) = data.extent(d) - 1;
57  if (ub(d) >= data.extent(d)) ub(d) = data.extent(d) - 1;
58  }
59  blitz::RectDomain<3> roi(lb, ub);
60  blitz::Array<DataT,3> cropped(ub - lb + 1);
61  cropped = data(roi);
62 
63  // Smooth ROI
64  blitz::TinyVector<double,3> sigma(elSize(1));
66  cropped, elSize, cropped, sigma, atb::BlitzIndexT(0), atb::RepeatBT);
67 
68  // Compute gradient
69  blitz::Array<blitz::TinyVector<DataT,3>,3> dCropped;
71  cropped, elSize, dCropped,
73  blitz::Array<DataT,3> dCroppedMag(cropped.shape());
74 
75  // Split up into magnitude and direction
76 #ifdef _OPENMP
77 #pragma omp parallel for
78 #endif
79  for (ptrdiff_t i = 0; i < dCropped.size(); ++i)
80  {
81  dCroppedMag.data()[i] = static_cast<DataT>(
82  std::sqrt(blitz::dot(dCropped.data()[i], dCropped.data()[i])));
83  if (dCroppedMag.data()[i] != 0.0)
84  dCropped.data()[i] /= dCroppedMag.data()[i];
85  }
86 
87  // Non maximum suppression
88  result.resize(data.shape());
89  result = 0;
90 
91  ptrdiff_t progressCnt = 0;
92  ptrdiff_t progressUpdate = 0;
93  if (progress != NULL)
94  progressUpdate = cropped.size() /
95  (progress->taskProgressMax() - progress->taskProgressMin());
96 
98 
99 #ifdef _OPENMP
100 #pragma omp parallel for
101 #endif
102  for (ptrdiff_t i = 0; i < cropped.size(); ++i)
103  {
104 
105  if (progress != NULL)
106  {
107  if (progress->isAborted()) continue;
108  if (progressCnt % progressUpdate == 0)
109  progress->updateProgress(
110  static_cast<int>(progressCnt / progressUpdate +
111  progress->taskProgressMin()));
112 #ifdef _OPENMP
113 #pragma omp atomic
114 #endif
115  ++progressCnt;
116  }
117 
118  blitz::TinyVector<ptrdiff_t,3> p;
119  ptrdiff_t tmp = i;
120  for (int d = 2; d >= 0; --d)
121  {
122  p(d) = tmp % cropped.extent(d);
123  tmp /= cropped.extent(d);
124  }
125  blitz::TinyVector<double,3> pl(
126  blitz::TinyVector<double,3>(p) - dCropped(p));
127  blitz::TinyVector<double,3> pu(
128  blitz::TinyVector<double,3>(p) + dCropped(p));
129 
130  if (dCroppedMag.data()[i] > ip.get(dCroppedMag, pl) &&
131  dCroppedMag.data()[i] > ip.get(dCroppedMag, pu))
132  result(blitz::TinyVector<ptrdiff_t,3>(p + lb)) =
133  dCroppedMag.data()[i];
134  }
135  }
136 
137  template<class DataT>
138  blitz::Array<DataT, 3> _calculateGradients(
139  blitz::Array<DataT, 3> const &data,
140  blitz::TinyVector<double,3> const &elSize,
141  int direction, double scaling,
142  iRoCS::ProgressReporter *progress)
143  {
144  //mit Randbehandlung
145  blitz::TinyVector<ptrdiff_t,3> dataShape(data.shape());
146 
147  DataT maxVal = blitz::max(data);
148  DataT minVal = blitz::min(data);
149 
150  blitz::Array<DataT,3> edgeNormedIm(dataShape);
151 
152  if (maxVal != minVal)
153  edgeNormedIm = (data - minVal) / (maxVal - minVal);
154 
155  blitz::Array<DataT,3> dx(dataShape);
156  dx = 0;
157 
158  DataT* data_in = edgeNormedIm.data();
159  DataT* data_out_x = dx.data();
160 
161  ptrdiff_t r2r_offs = dataShape(2);
162  ptrdiff_t l2l_offs = dataShape(1) * dataShape(2);
163 
164  switch (direction)
165  {
166  case 0:
167  {
168  // compute gradX
169  for (ptrdiff_t lev = 0; lev < dataShape(0); ++lev)
170  {
171  for (ptrdiff_t row = 0; row < dataShape(1); ++row)
172  {
173  ptrdiff_t offset = lev * l2l_offs + row * r2r_offs + 1;
174  for (ptrdiff_t col = 1; col < dataShape(2) - 1; ++col)
175  {
176  data_out_x[offset] =
177  static_cast<DataT>(
178  0.5 * (data_in[offset+1] - data_in[offset-1]) * scaling);
179  ++offset;
180  }
181  }
182  }
183  for (ptrdiff_t lev = 0; lev < dataShape(0); ++lev)
184  {
185  for (ptrdiff_t row = 0; row < dataShape(1); ++row)
186  {
187  dx(lev, row, static_cast<ptrdiff_t>(0)) =
188  dx(lev, row, static_cast<ptrdiff_t>(1));
189  dx(lev, row, dataShape(2) - 1) = dx(lev, row, dataShape(2) - 2);
190  }
191  }
192  return dx;
193  break;
194  }
195  case 1: {
196  // compute gradY
197  for (ptrdiff_t lev = 0; lev < dataShape(0); ++lev)
198  {
199  for (ptrdiff_t row = 1; row < dataShape(1) - 1; ++row)
200  {
201  ptrdiff_t offset = lev * l2l_offs + row * r2r_offs;
202  for (ptrdiff_t col = 0; col < dataShape(2); ++col)
203  {
204  data_out_x[offset] = static_cast<DataT>(
205  0.5 * (data_in[offset+r2r_offs]
206  - data_in[offset-r2r_offs]) * scaling);
207  ++offset;
208  }
209  }
210  }
211  for (ptrdiff_t lev = 0; lev < dataShape(0); ++lev)
212  {
213  for (ptrdiff_t col = 0; col < dataShape(2); ++col)
214  {
215  dx(lev, static_cast<ptrdiff_t>(0), col) =
216  dx(lev, static_cast<ptrdiff_t>(1), col);
217  dx(lev, dataShape(1) - 1, col) = dx(lev, dataShape(1) - 2, col);
218  }
219  }
220  return dx;
221  break;
222  }
223  case 2: {
224  // compute gradZ
225  for (ptrdiff_t lev = 1; lev < dataShape(0) - 1; ++lev)
226  {
227  for (ptrdiff_t row = 0; row < dataShape(1); ++row)
228  {
229  ptrdiff_t offset = lev * l2l_offs + row * r2r_offs;
230  for (ptrdiff_t col = 0; col < dataShape(2); ++col)
231  {
232  data_out_x[offset] = static_cast<DataT>(
233  0.5 * (data_in[offset+l2l_offs]
234  - data_in[offset-l2l_offs]) * scaling);
235  ++offset;
236  }
237  }
238  }
239  for (ptrdiff_t row = 0; row < dataShape(1); ++row)
240  {
241  for (ptrdiff_t col = 0; col < dataShape(2); ++col)
242  {
243  dx(static_cast<ptrdiff_t>(0), row, col) =
244  dx(static_cast<ptrdiff_t>(1), row, col);
245  dx(dataShape(0) - 1, row, col) = dx(dataShape(0) - 2, row, col);
246  }
247  }
248  return dx;
249  break;
250  }
251  default:
252  return data;
253  }
254  }
255 } //namespace
256 
257 #endif //guard
virtual bool isAborted() const =0
DataT get(blitz::Array< DataT, Dim > const &data, blitz::TinyVector< double, Dim > const &pos) const
Get the Array value at the given subpixel position.
The Implementation of the filter computing the gradient of the input data using central differences...
blitz::Array< DataT, 3 > _calculateGradients(blitz::Array< DataT, 3 > const &data, blitz::TinyVector< double, 3 > const &elSize, int direction, double scaling, iRoCS::ProgressReporter *progress)
The Implementation of the convolution of the input data with a Gaussian.
The LinearInterpolator class provides sub-pixel access to blitz++ Arrays using linear interpolation...
virtual bool updateProgress(int progress)=0
int taskProgressMax() const
int BlitzIndexT
The native integer type for indexing blitz++ Arrays.
Definition: TypeTraits.hh:56
int taskProgressMin() const
Mirror out-of-Array coordinates at the Array boundaries to get in-Array indices and values...
virtual void apply(blitz::Array< DataT, Dim > const &data, blitz::TinyVector< double, Dim > const &elementSizeUm, blitz::Array< ResultT, Dim > &filtered, iRoCS::ProgressReporter *pr=NULL) const
Apply the filter to the given Array.
Always pick the closest in-Array position to the position passed.
void edgeFilter(blitz::Array< DataT, 3 > const &data, blitz::TinyVector< double, 3 > const &elSize, blitz::Array< DataT, 3 > &result, blitz::TinyVector< double, 3 > const &lbUm, blitz::TinyVector< double, 3 > const &ubUm, iRoCS::ProgressReporter *progress)
calculate gradient magnitude with non-maximum-supression
double dot(blitz::Array< DataT, 1 > const &vec1, blitz::Array< DataT, 1 > const &vec2)
Calculation of the inner product (dot product) of the given vectors.
The CentralGradientFilter class implements the SeparableFilter interface and provides gradient comput...