iRoCS Toolbox  1.1.0
gvf-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 LIBSEGMENTATION_GVF_INL_HH
24 #define LIBSEGMENTATION_GVF_INL_HH
25 
26 #ifdef HAVE_CONFIG_H
27 #include <config.hh>
28 #endif
29 
30 #include "gvf.hh"
31 
33 
34 #ifdef HAVE_BLITZ_V9
35 #include <blitz/tinyvec-et.h>
36 #endif
37 
38 namespace segmentation
39 {
40 
50  template<typename T>
52  blitz::Array<blitz::TinyVector<T, 3>, 3> &gradient,
53  blitz::TinyVector<T,3> const &el_size_um, T mu, T nu, int max_iter,
54  iRoCS::ProgressReporter *progress)
55  {
56  blitz::TinyVector<atb::BlitzIndexT,3> shape(gradient.shape());
57  // normalize gradient and gradient magnitude
58  blitz::Array<T,3> gradient_norm_sq(shape);
59 
60  T max_norm_sq = T(0);
61 #ifdef _OPENMP
62 #pragma omp parallel for
63 #endif
64 #if defined _OPENMP && defined _WIN32
65  for (ptrdiff_t i = 0; i < static_cast<ptrdiff_t>(gradient.size()); ++i) {
66 #else
67  for (size_t i = 0; i < gradient.size(); ++i) {
68 #endif
69  gradient_norm_sq.data()[i] =
70  gradient.data()[i](0) * gradient.data()[i](0)
71  + gradient.data()[i](1) * gradient.data()[i](1)
72  + gradient.data()[i](2) * gradient.data()[i](2);
73 #ifdef _OPENMP
74 #pragma omp critical
75 #endif
76  if (gradient_norm_sq.data()[i] > max_norm_sq)
77  max_norm_sq = gradient_norm_sq.data()[i];
78  }
79 
80 #ifdef _OPENMP
81 #pragma omp parallel for
82 #endif
83 #if defined _OPENMP && defined _WIN32
84  for (ptrdiff_t i = 0; i < static_cast<ptrdiff_t>(gradient.size()); ++i) {
85 #else
86  for (size_t i = 0; i < gradient.size(); ++i) {
87 #endif
88  gradient.data()[i] /= std::sqrt(max_norm_sq);
89  gradient_norm_sq.data()[i] /= max_norm_sq;
90  }
91 
92  // normalizer
93 
94  T max_el_size = blitz::max(el_size_um);
95  blitz::TinyVector<T,3> normalizer(
96  (el_size_um / max_el_size) * (el_size_um / max_el_size));
97  normalizer = 1;
98  std::cout << "GVF - normalizer " << normalizer << std::endl;
99 
100  // const part in linear equation
101  blitz::Array<blitz::TinyVector<T, 3>,3> b_arr(shape);
102  b_arr = gradient_norm_sq * gradient;
103 
104  // short notation
105  blitz::Array<blitz::TinyVector<T, 3>, 3>& x = gradient;
106 
107  for (int iter = 1; iter <= max_iter; ++iter)
108  {
109  if (progress != NULL && !progress->updateProgress(
110  static_cast<double>(iter) / static_cast<double>(max_iter) *
111  progress->taskProgressMax() - progress->taskProgressMin()) +
112  progress->taskProgressMin()) break;
113 #ifdef _OPENMP
114 #pragma omp parallel for
115 #endif
116  for (int lev = 1; lev < shape(0) - 1; ++lev)
117  {
118  for (int row = 1; row < shape(1) - 1; ++row)
119  {
120  for (int col = 1; col < shape(2) - 1; ++col)
121  {
122  // own coefficients
123  // laplacian - squared gradient magnitude
124  T a_ii = -6 - gradient_norm_sq(lev, row, col)/mu;
125  // const part
126  const blitz::TinyVector<T, 3>& b = -b_arr(lev, row, col)/mu;
127 
128  // other coefficients * x
129  // this is just the finite-difference laplacian
130  blitz::TinyVector<T, 3> a_ij_x =
131  (x(lev - 1, row, col) +
132  x(lev + 1, row, col)) * normalizer(0) +
133  (x(lev, row - 1, col) +
134  x(lev, row + 1, col)) * normalizer(1) +
135  (x(lev, row, col - 1) +
136  x(lev, row, col + 1)) * normalizer(2);
137  // SOR iteration
138  x(lev, row, col) = (1. - nu) * x(lev, row, col) +
139  nu/a_ii * (b - a_ij_x);
140  }
141  }
142  }
143  } // for maxiter
144  }
145 
146  template<typename T>
148  blitz::Array<blitz::TinyVector<T, 3>, 3> &gradient,
149  blitz::TinyVector<T,3> const &el_size_um, T mu, T hs, T hr, int max_iter,
150  iRoCS::ProgressReporter *progress)
151  {
152  std::cout << "MSGVF - Initializing" << std::endl;
153  T scaling=el_size_um(2)/el_size_um(0);
154  T hz=hs*scaling;
155  blitz::TinyVector<atb::BlitzIndexT,3> Shape(gradient.shape());
156  hr*=hr;
157  hz=std::floor(hz);
158  atb::BlitzIndexT hz_int = static_cast<atb::BlitzIndexT>(hz);
159  atb::BlitzIndexT hs_int = static_cast<atb::BlitzIndexT>(hs);
160  blitz::Array<T,3> dx2(
161  Shape(0) + 2 * hz_int, Shape(1) + 2 * hs_int, Shape(2) + 2 * hs_int);
162  blitz::Array<T,3> dy2(dx2.shape());
163  blitz::Array<T,3> dz2(dx2.shape());
164 
165  blitz::Array<T,3> gradient_norm(Shape);
166  gradient_norm =
167  blitz::sqrt(blitz::pow2(gradient[0]) +
168  blitz::pow2(gradient[1]) +
169  blitz::pow2(gradient[2]));
170  T maxgrad=blitz::max(gradient_norm);
171 #ifdef _OPENMP
172 #pragma omp parallel for
173 #endif
174  for (ptrdiff_t i = 0; i < static_cast<ptrdiff_t>(gradient.size()); ++i)
175  {
176  gradient.data()[i] /= maxgrad;
177  gradient_norm.data()[i] /= maxgrad;
178  }
179 
180  dx2=0;
181  dy2=0;
182  dz2=0;
183 
184  dx2(blitz::Range(hz_int, hz_int + Shape(0) - 1),
185  blitz::Range(hs_int, hs_int + Shape(1) - 1),
186  blitz::Range(hs_int, hs_int + Shape(2) - 1)) = gradient[0];
187  dy2(blitz::Range(hz_int, hz_int + Shape(0) - 1),
188  blitz::Range(hs_int, hs_int + Shape(1) - 1),
189  blitz::Range(hs_int, hs_int + Shape(2) - 1)) = gradient[1];
190  dz2(blitz::Range(hz_int, hz_int + Shape(0) - 1),
191  blitz::Range(hs_int, hs_int + Shape(1) - 1),
192  blitz::Range(hs_int, hs_int + Shape(2) - 1)) = gradient[2];
193 
194  blitz::TinyVector<atb::BlitzIndexT,3> ShapeL(dx2.shape());
195  blitz::Array<T,3> b(ShapeL);
196 
197  b = dx2 * dx2 + dy2 * dy2 + dz2 * dz2;
198 
199  blitz::Array<T,3> c1(ShapeL);
200  c1 = b * dx2;
201  blitz::Array<T,3> c2(ShapeL);
202  c2 = b * dy2;
203  blitz::Array<T,3> c3(ShapeL);
204  c3 = b * dz2;
205  b = static_cast<T>(1) - b;
206  blitz::Array<blitz::TinyVector<T,3>,3> ms_f(ShapeL);
207  ms_f[0]=0;
208  ms_f[1]=0;
209  ms_f[2]=0;
210 
211  T hssq = hs * hs;
212  scaling *= scaling;
213  T invscaling = static_cast<T>(1.0) / scaling;
214 
215  std::cout << "MSGVF - Running" << std::endl;
216 
217  int count;
218  /************* Solve GVF = (u,v) iteratively ***************/
219  for (count = 1; count <= max_iter; count++)
220  {
221  if (progress != NULL && !progress->updateProgress(
222  static_cast<double>(count) / static_cast<double>(max_iter) *
223  (progress->taskProgressMax() - progress->taskProgressMin()) +
224  progress->taskProgressMin())) return;
225 
226  // /* IV: Compute Laplace operator using Neuman condition */
227  // blitz::Array<float,2> Lu = laplacian(u);
228  // blitz::Array<float,2> Lv = laplacian(v);
229  //
230  // /******** V: Update GVF ************/
231  // u = (1-b) * u + mu * Lu + c1;
232  // v = (1-b) * v + mu * Lv + c2;
233 
234 #ifdef _OPENMP
235 #pragma omp parallel for
236 #endif
237  for(atb::BlitzIndexT lev = hz_int; lev < ShapeL(0) - hz_int; ++lev)
238  {
239  if (progress != NULL && progress->isAborted()) continue;
240  for(atb::BlitzIndexT row = hs_int; row < ShapeL(1) - hs_int; ++row)
241  {
242  for(atb::BlitzIndexT col = hs_int; col < ShapeL(2) - hs_int; ++col)
243  {
244  ms_f(lev, row, col)[0]=0;
245  ms_f(lev, row, col)[1]=0;
246  ms_f(lev, row, col)[2]=0;
247 
248  T count=0;
249 
250  T dx_l=dx2(lev, row, col);
251  T dy_l=dy2(lev, row, col);
252  T dz_l=dz2(lev, row, col);
253 
254  T distv=0;
255 
256  //float dweight=1;
257  //for spherical kernel
258  for(atb::BlitzIndexT lev2 = lev - hz_int;
259  lev2 < lev + hz_int + 1; ++lev2)
260  {
261  for(atb::BlitzIndexT row2 = row - hs_int;
262  row2 < row + hs_int + 1; ++row2)
263  {
264  for(atb::BlitzIndexT col2 = col - hs_int;
265  col2 < col + hs_int + 1; ++col2)
266  {
267  T dx_l2 = dx2(lev2, row2, col2);
268  T dy_l2 = dy2(lev2, row2, col2);
269  T dz_l2 = dz2(lev2, row2, col2);
270  distv = ((dx_l2 - dx_l) * (dx_l2 - dx_l) +
271  (dy_l2 - dy_l) * (dy_l2 - dy_l) +
272  (dz_l2 - dz_l) * (dz_l2 - dz_l)) / hr +
273  ((invscaling * static_cast<T>(lev - lev2) *
274  static_cast<T>(lev - lev2)) +
275  (static_cast<T>(row - row2) *
276  static_cast<T>(row - row2) +
277  static_cast<T>(col - col2) *
278  static_cast<T>(col - col2))) / hssq;
279 
280  if(distv < 1)
281  {
282  //float dweight=1-distv;
283  //count+=dweight;
284  count++;
285  ms_f(lev, row, col)[0] += dx_l2;
286  ms_f(lev, row, col)[1] += dy_l2;
287  ms_f(lev, row, col)[2] += dz_l2;
288  }
289  }
290  }
291  }
292  if(count > 0.5)
293  {
294  ms_f(lev, row, col)[0] /= count;
295  ms_f(lev, row, col)[1] /= count;
296  ms_f(lev, row, col)[2] /= count;
297  }
298  }
299  }
300  }
301  if (progress != NULL && progress->isAborted()) return;
302 
303  dx2 = ms_f[0];
304 #ifdef _OPENMP
305 #pragma omp parallel for
306 #endif
307  for(atb::BlitzIndexT lev = 1; lev < ShapeL(0) - 1; ++lev)
308  {
309  if (progress != NULL && progress->isAborted()) continue;
310  for(atb::BlitzIndexT row = 1; row < ShapeL(1) - 1; ++row)
311  {
312  for(atb::BlitzIndexT col = 1; col < ShapeL(2) - 1; ++col)
313  {
314  dx2(lev, row, col) =
315  ((b(lev, row, col)) * dx2(lev, row, col) +
316  mu * ((dx2(lev - 1, row, col) +
317  dx2(lev + 1, row, col) -
318  2 * dx2(lev, row, col)) * scaling +
319  dx2(lev, row - 1, col) +
320  dx2(lev, row + 1, col) +
321  dx2(lev, row, col - 1) +
322  dx2(lev, row, col + 1) -
323  4 * dx2(lev, row, col)) + c1(lev, row, col));
324  }
325  }
326  }
327  if (progress != NULL && progress->isAborted()) return;
328 
329  dy2 = ms_f[1];
330 #ifdef _OPENMP
331 #pragma omp parallel for
332 #endif
333  for(atb::BlitzIndexT lev = 1; lev < ShapeL(0) - 1; ++lev)
334  {
335  if (progress != NULL && progress->isAborted()) continue;
336  for(atb::BlitzIndexT row = 1; row < ShapeL(1) - 1; ++row)
337  {
338  for(atb::BlitzIndexT col = 1; col < ShapeL(2) - 1; ++col)
339  {
340  dy2(lev, row, col) =
341  ((b(lev, row, col)) * dy2(lev, row, col) +
342  mu * ((dy2(lev - 1, row, col) +
343  dy2(lev + 1, row, col) -
344  2 * dy2(lev, row, col)) * scaling +
345  dy2(lev, row - 1, col) +
346  dy2(lev, row + 1, col) +
347  dy2(lev, row, col - 1) +
348  dy2(lev, row, col + 1) -
349  4 * dy2(lev, row, col)) + c2(lev, row, col));
350  }
351  }
352  }
353  if (progress != NULL && progress->isAborted()) return;
354 
355  dz2 = ms_f[2];
356 #ifdef _OPENMP
357 #pragma omp parallel for
358 #endif
359  for(atb::BlitzIndexT lev = 1; lev < ShapeL(0) - 1; ++lev)
360  {
361  if (progress != NULL && progress->isAborted()) continue;
362  for(atb::BlitzIndexT row = 1; row < ShapeL(1) - 1; ++row)
363  {
364  for(atb::BlitzIndexT col = 1; col < ShapeL(2) - 1; ++col)
365  {
366  dz2(lev, row, col) =
367  ((b(lev, row, col)) * dz2(lev, row, col) +
368  mu * ((dz2(lev - 1, row, col) +
369  dz2(lev + 1, row, col) -
370  2 * dz2(lev, row, col)) * scaling +
371  dz2(lev, row - 1, col) +
372  dz2(lev, row + 1, col) +
373  dz2(lev, row, col - 1) +
374  dz2(lev, row, col + 1) -
375  4 * dz2(lev, row, col)) + c3(lev, row, col));
376  }
377  }
378  }
379  }
380  if (progress != NULL && progress->isAborted()) return;
381 
382  gradient[0] = dx2(
383  blitz::Range(hz_int, hz_int + Shape(0) - 1),
384  blitz::Range(hs_int, hs_int + Shape(1) - 1),
385  blitz::Range(hs_int, hs_int + Shape(2) - 1));
386  gradient[1] = dy2(
387  blitz::Range(hz_int, hz_int + Shape(0) - 1),
388  blitz::Range(hs_int, hs_int + Shape(1) - 1),
389  blitz::Range(hs_int, hs_int + Shape(2) - 1));
390  gradient[2] = dz2(
391  blitz::Range(hz_int, hz_int + Shape(0) - 1),
392  blitz::Range(hs_int, hs_int + Shape(1) - 1),
393  blitz::Range(hs_int, hs_int + Shape(2) - 1));
394 
395  std::cout << "MSGVF done after " << count - 1 << " of " << max_iter
396  << " iterations" << std::endl;
397  }
398 
399 } // namespace
400 
401 #endif //guard
virtual bool isAborted() const =0
virtual bool updateProgress(int progress)=0
void msGradientVectorFlow(blitz::Array< blitz::TinyVector< T, 3 >, 3 > &gradient, blitz::TinyVector< T, 3 > const &el_size_um, T mu, T hs, T hr, int max_iter, iRoCS::ProgressReporter *progress)
Definition: gvf-inl.hh:147
int taskProgressMax() const
Query specific information about different data types.
int BlitzIndexT
The native integer type for indexing blitz++ Arrays.
Definition: TypeTraits.hh:56
int taskProgressMin() const
void gradientVectorFlowSOR(blitz::Array< blitz::TinyVector< T, 3 >, 3 > &gradient, blitz::TinyVector< T, 3 > const &el_size_um, T mu, T nu, int max_iter, iRoCS::ProgressReporter *progress)
Solve Euler-Lagrange equation for gradient vector flow using successive over-relaxation u : output g...
Definition: gvf-inl.hh:51