23 #ifndef LIBSEGMENTATION_GVF_INL_HH 24 #define LIBSEGMENTATION_GVF_INL_HH 35 #include <blitz/tinyvec-et.h> 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,
56 blitz::TinyVector<atb::BlitzIndexT,3> shape(gradient.shape());
58 blitz::Array<T,3> gradient_norm_sq(shape);
62 #pragma omp parallel for 64 #if defined _OPENMP && defined _WIN32 65 for (ptrdiff_t i = 0; i < static_cast<ptrdiff_t>(gradient.size()); ++i) {
67 for (
size_t i = 0; i < gradient.size(); ++i) {
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);
76 if (gradient_norm_sq.data()[i] > max_norm_sq)
77 max_norm_sq = gradient_norm_sq.data()[i];
81 #pragma omp parallel for 83 #if defined _OPENMP && defined _WIN32 84 for (ptrdiff_t i = 0; i < static_cast<ptrdiff_t>(gradient.size()); ++i) {
86 for (
size_t i = 0; i < gradient.size(); ++i) {
88 gradient.data()[i] /= std::sqrt(max_norm_sq);
89 gradient_norm_sq.data()[i] /= max_norm_sq;
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));
98 std::cout <<
"GVF - normalizer " << normalizer << std::endl;
101 blitz::Array<blitz::TinyVector<T, 3>,3> b_arr(shape);
102 b_arr = gradient_norm_sq * gradient;
105 blitz::Array<blitz::TinyVector<T, 3>, 3>& x = gradient;
107 for (
int iter = 1; iter <= max_iter; ++iter)
110 static_cast<double>(iter) / static_cast<double>(max_iter) *
114 #pragma omp parallel for 116 for (
int lev = 1; lev < shape(0) - 1; ++lev)
118 for (
int row = 1; row < shape(1) - 1; ++row)
120 for (
int col = 1; col < shape(2) - 1; ++col)
124 T a_ii = -6 - gradient_norm_sq(lev, row, col)/mu;
126 const blitz::TinyVector<T, 3>& b = -b_arr(lev, row, col)/mu;
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);
138 x(lev, row, col) = (1. - nu) * x(lev, row, col) +
139 nu/a_ii * (b - a_ij_x);
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,
152 std::cout <<
"MSGVF - Initializing" << std::endl;
153 T scaling=el_size_um(2)/el_size_um(0);
155 blitz::TinyVector<atb::BlitzIndexT,3> Shape(gradient.shape());
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());
165 blitz::Array<T,3> gradient_norm(Shape);
167 blitz::sqrt(blitz::pow2(gradient[0]) +
168 blitz::pow2(gradient[1]) +
169 blitz::pow2(gradient[2]));
170 T maxgrad=blitz::max(gradient_norm);
172 #pragma omp parallel for 174 for (ptrdiff_t i = 0; i < static_cast<ptrdiff_t>(gradient.size()); ++i)
176 gradient.data()[i] /= maxgrad;
177 gradient_norm.data()[i] /= maxgrad;
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];
194 blitz::TinyVector<atb::BlitzIndexT,3> ShapeL(dx2.shape());
195 blitz::Array<T,3> b(ShapeL);
197 b = dx2 * dx2 + dy2 * dy2 + dz2 * dz2;
199 blitz::Array<T,3> c1(ShapeL);
201 blitz::Array<T,3> c2(ShapeL);
203 blitz::Array<T,3> c3(ShapeL);
205 b =
static_cast<T
>(1) - b;
206 blitz::Array<blitz::TinyVector<T,3>,3> ms_f(ShapeL);
213 T invscaling =
static_cast<T
>(1.0) / scaling;
215 std::cout <<
"MSGVF - Running" << std::endl;
219 for (count = 1; count <= max_iter; count++)
222 static_cast<double>(count) / static_cast<double>(max_iter) *
235 #pragma omp parallel for 239 if (progress != NULL && progress->
isAborted())
continue;
244 ms_f(lev, row, col)[0]=0;
245 ms_f(lev, row, col)[1]=0;
246 ms_f(lev, row, col)[2]=0;
250 T dx_l=dx2(lev, row, col);
251 T dy_l=dy2(lev, row, col);
252 T dz_l=dz2(lev, row, col);
259 lev2 < lev + hz_int + 1; ++lev2)
262 row2 < row + hs_int + 1; ++row2)
265 col2 < col + hs_int + 1; ++col2)
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;
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;
294 ms_f(lev, row, col)[0] /= count;
295 ms_f(lev, row, col)[1] /= count;
296 ms_f(lev, row, col)[2] /= count;
301 if (progress != NULL && progress->
isAborted())
return;
305 #pragma omp parallel for 309 if (progress != NULL && progress->
isAborted())
continue;
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));
327 if (progress != NULL && progress->
isAborted())
return;
331 #pragma omp parallel for 335 if (progress != NULL && progress->
isAborted())
continue;
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));
353 if (progress != NULL && progress->
isAborted())
return;
357 #pragma omp parallel for 361 if (progress != NULL && progress->
isAborted())
continue;
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));
380 if (progress != NULL && progress->
isAborted())
return;
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));
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));
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));
395 std::cout <<
"MSGVF done after " << count - 1 <<
" of " << max_iter
396 <<
" iterations" << std::endl;
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)
int taskProgressMax() const
Query specific information about different data types.
int BlitzIndexT
The native integer type for indexing blitz++ Arrays.
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...