iRoCS Toolbox  1.1.0
ProfileFilter-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_SRC_PROFILEFILTER_INL_HH
24 #define LIBSEGMENTATION_SRC_PROFILEFILTER_INL_HH
25 
26 #ifdef HAVE_CONFIG_H
27 #include <config.hh>
28 #endif
29 
30 #include "NormalPDF.hh"
31 
32 #include "ProfileFilter.hh"
33 #include "ProfileSampler.hh"
34 
35 #ifdef HAVE_BLITZ_V9
36 #include <blitz/tinyvec-et.h>
37 #endif
38 
39 namespace segmentation
40 {
41 
42  template<class T, int Dim>
43  blitz::RectDomain<Dim>
45  ProfileSampler<T,Dim> const &sampler,
46  std::vector<NormalPDF<T>*> const &pdfs,
47  blitz::TinyVector<double,Dim> centerUm,
48  blitz::Array<T,Dim> const &data,
49  blitz::TinyVector<double,Dim> const &elSize,
50  blitz::Array<T,Dim> &result,
51  blitz::TinyVector<double,Dim> lbUm,
52  blitz::TinyVector<double,Dim> const &ubUm,
53  iRoCS::ProgressReporter *progress) const
54  {
55  blitz::TinyVector<ptrdiff_t,Dim> lbPx(lbUm / elSize);
56  blitz::TinyVector<ptrdiff_t,Dim> ubPx(ubUm / elSize);
57 
58  for (int d = 0; d < Dim; ++d)
59  {
60  if (lbPx(d) < 0) lbPx(d) = 0;
61  if (ubPx(d) < 0) ubPx(d) = 0;
62  if (lbPx(d) >= data.extent(d)) lbPx(d) = data.extent(d) - 1;
63  if (ubPx(d) >= data.extent(d)) ubPx(d) = data.extent(d) - 1;
64  }
65  blitz::RectDomain<Dim> roi(lbPx, ubPx);
66  lbUm = blitz::TinyVector<double,Dim>(lbPx) * elSize;
67 
68  typename blitz::Array<T,Dim>::iterator itDst;
69  blitz::TinyVector<ptrdiff_t,Dim> outShape(ubPx - lbPx + 1);
70  result.resize(outShape);
71 
73  T maxVal = -std::numeric_limits<T>::infinity();
74 
75  ptrdiff_t progressCnt = 0;
76  ptrdiff_t progressUpdate = 0;
77  if (progress != NULL)
78  progressUpdate = result.size() /
79  (progress->taskProgressMax() - progress->taskProgressMin());
80 
81 #ifdef _OPENMP
82 #pragma omp parallel for
83 #endif
84 #if defined _OPENMP && defined _WIN32
85  for (ptrdiff_t i = 0; i < static_cast<ptrdiff_t>(result.size()); ++i) {
86 #else
87  for (size_t i = 0; i < result.size(); ++i) {
88 #endif
89  if (progress != NULL)
90  {
91  if (progress->isAborted()) continue;
92  if (progressCnt % progressUpdate == 0)
93  progress->updateProgress(
94  static_cast<int>(progressCnt / progressUpdate +
95  progress->taskProgressMin()));
96 #ifdef _OPENMP
97 #pragma omp atomic
98 #endif
99  ++progressCnt;
100  }
101 
102  blitz::TinyVector<ptrdiff_t,Dim> srcPosPx;
103  ptrdiff_t tmp = i;
104  for (int d = Dim - 1; d >= 0; --d)
105  {
106  srcPosPx(d) = tmp % result.extent(d);
107  tmp /= result.extent(d);
108  }
109  srcPosPx += lbPx;
110 
111  blitz::TinyVector<double,Dim> srcPosUm(srcPosPx);
112  srcPosUm *= elSize;
113 
114  // distance to nearest pdf
115  blitz::Array<T,1> profile(
116  sampler.sample(data, elSize, centerUm, srcPosUm));
117 
118  T minDist = std::numeric_limits<T>::infinity();
119  for (size_t j = 0; j < pdfs.size(); ++j)
120  {
121  T d = dist(*(pdfs[j]), profile);
122  if (d < minDist) minDist = d;
123  }
124 
125  if(minDist == std::numeric_limits<T>::infinity())
126  {
127  if (blitz::any(srcPosUm != centerUm))
128  {
129 #ifdef _OPENMP
130 #pragma omp critical
131 #endif
132  std::cerr << "ProfileFilter: WARNING: nan value detected and set to "
133  << "infinity at position" << srcPosUm << std::endl;
134  }
135  minDist = -1;
136  }
137 
138  // filter response distance to nearest pdf
139  result.data()[i] = minDist;
140 #ifdef _OPENMP
141 #pragma omp critical
142 #endif
143  {
144  if (minDist > maxVal) maxVal = minDist;
145  }
146  }
147 #ifdef _OPENMP
148 #pragma omp parallel for
149 #endif
150 #if defined _OPENMP && defined _WIN32
151  for (ptrdiff_t i = 0; i < static_cast<ptrdiff_t>(result.size()); ++i)
152 #else
153  for (size_t i = 0; i < result.size(); ++i)
154 #endif
155  result.data()[i] = (result.data()[i] == -1) ? maxVal : result.data()[i];
156 
157  std::cout << "profile filter max val " << maxVal << std::endl;
158 
159  // invert grayvalues and normalize
160 #ifdef _OPENMP
161 #pragma omp parallel for
162 #endif
163 #if defined _OPENMP && defined _WIN32
164  for (ptrdiff_t i = 0; i < static_cast<ptrdiff_t>(result.size()); ++i)
165 #else
166  for (size_t i = 0; i < result.size(); ++i)
167 #endif
168  result.data()[i] = static_cast<T>(1) - result.data()[i] / maxVal;
169 
170  return roi;
171  }
172 
173 }
174 #endif //LIBSEGMENTATION_SRC_PROFILEFILTER_INL_HH
virtual bool isAborted() const =0
virtual bool updateProgress(int progress)=0
int taskProgressMax() const
Encapsulates a normal pdf.
Definition: NormalPDF.hh:40
bool any(blitz::TinyMatrix< bool, NRows, NColumns > const &matrix)
any() reduction for boolean blitz::TinyMatrix.
int taskProgressMin() const
blitz::RectDomain< Dim > operator()(ProfileSampler< T, Dim > const &sampler, std::vector< NormalPDF< T > *> const &pdfs, blitz::TinyVector< double, Dim > centerUm, blitz::Array< T, Dim > const &data, blitz::TinyVector< double, Dim > const &elSize, blitz::Array< T, Dim > &result, blitz::TinyVector< double, Dim > lbUm=0.0, blitz::TinyVector< double, Dim > const &ubUm=(std::numeric_limits< double >::infinity()), iRoCS::ProgressReporter *progress=NULL) const
Filters input data.
blitz::Array< T, 1 > sample(blitz::Array< T, Dim > const &data, blitz::TinyVector< double, Dim > const &elSize, blitz::TinyVector< double, Dim > const &centerUm, blitz::TinyVector< double, Dim > const &positionUm) const