iRoCS Toolbox  1.1.0
ProfileSampler-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_PROFILESAMPLER_INL_HH
24 #define LIBSEGMENTATION_SRC_PROFILESAMPLER_INL_HH
25 
26 #ifdef HAVE_CONFIG_H
27 #include <config.hh>
28 #endif
29 
30 #include "ProfileSampler.hh"
31 
32 #include <omp.h>
33 
36 
37 namespace segmentation
38 {
39 
40  template<class T, int Dim>
41  ProfileSampler<T,Dim>::ProfileSampler(int numSamples, T samplingDistUm)
42  : _numSamples(numSamples), _samplingDistUm(samplingDistUm)
43  {}
44 
45  template<class T, int Dim>
46  blitz::Array<T,1> ProfileSampler<T, Dim>::sample(
47  blitz::Array<T,Dim> const &data,
48  blitz::TinyVector<double,Dim> const &elSize,
49  blitz::TinyVector<double,Dim> const &centerUm,
50  blitz::TinyVector<double,Dim> const &positionUm) const
51  {
52  // sampling direction, from center to edge
53  blitz::TinyVector<double,Dim> samplingStepUm(positionUm - centerUm);
54 
55  // scale with sampling dist to optain stampling step vector
56  // sampling step points away from the center
57  samplingStepUm *= _samplingDistUm /
58  std::sqrt(blitz::dot(samplingStepUm, samplingStepUm));
59 
60  // half total sampling vector, used as offset to sample in both
61  // directions around position_um staring from the point closest to the
62  // center
63  blitz::TinyVector<double,Dim> samplingOffsetUm(
64  0.5 * samplingStepUm * _numSamples);
65  blitz::TinyVector<double,Dim> samplingStartUm(
66  positionUm - samplingOffsetUm);
67 
69  blitz::Array<T, 1> edgeSamples(_numSamples);
70  for (int i = 0; i < _numSamples; ++i)
71  {
72  blitz::TinyVector<double,Dim> samplePosUm(
73  samplingStartUm + samplingStepUm * static_cast<double>(i));
74  blitz::TinyVector<double,Dim> samplePosPx(samplePosUm / elSize);
75  edgeSamples(i) = ip.get(data, samplePosPx);
76  }
77 
78  // forward finite difference approximation of edge_samples
79  blitz::Array<T,1> dEdge(_numSamples - 1);
80  for (int i = 0; i < _numSamples - 1; ++i)
81  dEdge(i) = edgeSamples(i + 1) - edgeSamples(i);
82 
83  // normalize vector
84  T normalizer = blitz::max(blitz::abs(dEdge));
85 
86  // we don't like our results to be screwed up by inf values
87  if(normalizer > 0) dEdge /= normalizer;
88 
89  return dEdge;
90  }
91 
92  template<class T, int Dim>
94  blitz::Array<T,Dim> const &data,
95  blitz::TinyVector<double,Dim> const &elSize,
96  blitz::TinyVector<double,Dim> const &centerUm,
97  std::vector< blitz::TinyVector<double,Dim> > const &candidates,
98  blitz::Array<T,2> &profiles,
99  iRoCS::ProgressReporter *progress) const
100  {
101  profiles.resize(candidates.size(), _numSamples - 1);
102 
103  ptrdiff_t progressCnt = 0;
104  ptrdiff_t progressUpdate = 0;
105  if (progress != NULL)
106  progressUpdate = candidates.size() /
107  (progress->taskProgressMax() - progress->taskProgressMin());
108 
109 #ifdef _OPENMP
110 #pragma omp parallel for
111 #endif
112  for (ptrdiff_t i = 0; i < static_cast<ptrdiff_t>(candidates.size()); ++i)
113  {
114  if (progress != NULL)
115  {
116  if (progress->isAborted()) continue;
117  if (progressCnt % progressUpdate == 0)
118  progress->updateProgress(
119  static_cast<int>(progressCnt / progressUpdate +
120  progress->taskProgressMin()));
121 #ifdef _OPENMP
122 #pragma omp atomic
123 #endif
124  ++progressCnt;
125  }
126  profiles(i, blitz::Range::all()) =
127  sample(data, elSize, centerUm, candidates[i]);
128  }
129  }
130 
131  template<class T, int Dim>
133  blitz::Array<T,Dim> const &data,
134  blitz::TinyVector<double,Dim> const &elSize,
135  blitz::TinyVector<double,Dim> const &centerUm,
136  blitz::Array<T,Dim> const &candidates, T threshold,
137  blitz::Array<T,2> &profiles,
138  std::vector< blitz::TinyVector<double,Dim> > &usedCandidates,
139  blitz::TinyVector<double,Dim> const &lbUm,
140  blitz::TinyVector<double,Dim> const &ubUm,
141  iRoCS::ProgressReporter *progress) const
142  {
143  usedCandidates.clear();
144 
145  blitz::TinyVector<ptrdiff_t,Dim> lb(lbUm / elSize);
146  blitz::TinyVector<ptrdiff_t,Dim> ub(ubUm / elSize);
147  for (int d = 0; d < Dim; ++d)
148  {
149  if (lb(d) < 0) lb(d) = 0;
150  if (ub(d) < 0) ub(d) = 0;
151  if (lb(d) >= data.extent(d)) lb(d) = data.extent(d) - 1;
152  if (ub(d) >= data.extent(d)) ub(d) = data.extent(d) - 1;
153  }
154  blitz::RectDomain<Dim> roi(lb, ub);
155 
156  typename blitz::Array<T, Dim>::const_iterator it;
157  for (it = candidates(roi).begin(); it != candidates(roi).end(); ++it)
158  {
159  if (*it < threshold) continue;
160  blitz::TinyVector<double,Dim> positionUm(it.position() + lb);
161  positionUm *= elSize;
162  usedCandidates.push_back(positionUm);
163  }
164  sample(data, elSize, centerUm, usedCandidates, profiles, progress);
165  }
166 
167 } //namespace
168 #endif //LIBSEGMENTATION_SRC_PROFILESAMPLER_INL_HH
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.
ProfileSampler(int numSamples, T samplingDistUm)
The LinearInterpolator class provides sub-pixel access to blitz++ Arrays using linear interpolation...
virtual bool updateProgress(int progress)=0
int taskProgressMax() const
int taskProgressMin() const
Classes and functions for sub-pixel Array access with different interpolation strategies.
bool all(blitz::TinyMatrix< bool, NRows, NColumns > const &matrix)
all() reduction for boolean blitz::TinyMatrix.
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.
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