UsTK : Ultrasound ToolKit  version 2.0.1 under development (2023-12-07)
usRFToPreScan2DConverter.cpp
1 /****************************************************************************
2  *
3  * This file is part of the ustk software.
4  * Copyright (C) 2016 - 2017 by Inria. All rights reserved.
5  *
6  * This software is free software; you can redistribute it and/or
7  * modify it under the terms of the GNU General Public License
8  * ("GPL") version 2 as published by the Free Software Foundation.
9  * See the file LICENSE.txt at the root directory of this source
10  * distribution for additional information about the GNU GPL.
11  *
12  * For using ustk with software that can not be combined with the GNU
13  * GPL, please contact Inria about acquiring a ViSP Professional
14  * Edition License.
15  *
16  * This software was developed at:
17  * Inria Rennes - Bretagne Atlantique
18  * Campus Universitaire de Beaulieu
19  * 35042 Rennes Cedex
20  * France
21  *
22  * If you have questions regarding the use of this file, please contact
23  * Inria at ustk@inria.fr
24  *
25  * This file is provided AS IS with NO WARRANTY OF ANY KIND, INCLUDING THE
26  * WARRANTY OF DESIGN, MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE.
27  *
28  * Authors:
29  * Pedro Patlan
30  * Marc Pouliquen
31  *
32  *****************************************************************************/
33 
39 #include <visp3/ustk_core/usRFToPreScan2DConverter.h>
40 
41 #if defined(USTK_HAVE_FFTW)
42 
48  : m_logCompressor(), m_decimationFactor(decimationFactor), m_isInit(false)
49 {
50 }
51 
56 {
57  if (m_isInit) {
58  fftw_free(m_fft_in);
59  fftw_free(m_fft_out);
60  fftw_free(m_fft_conv);
61  fftw_free(m_fft_out_inv);
62  fftw_destroy_plan(m_p);
63  fftw_destroy_plan(m_pinv);
64  delete m_env;
65  delete m_comp;
66  }
67 }
68 
74 void usRFToPreScan2DConverter::init(int widthRF, int heigthRF)
75 {
76  if (m_isInit && (m_signalSize != heigthRF || m_scanLineNumber != widthRF)) {
77  fftw_free(m_fft_in);
78  fftw_free(m_fft_out);
79  fftw_free(m_fft_conv);
80  fftw_free(m_fft_out_inv);
81  fftw_destroy_plan(m_p);
82  fftw_destroy_plan(m_pinv);
83  delete m_env;
84  delete m_comp;
85  } else if (m_signalSize == heigthRF && m_scanLineNumber == heigthRF)
86  return;
87 
88  // for FFT
89  m_fft_in = (fftw_complex *)fftw_malloc(sizeof(fftw_complex) * heigthRF);
90  m_fft_out = (fftw_complex *)fftw_malloc(sizeof(fftw_complex) * heigthRF);
91  m_fft_conv = (fftw_complex *)fftw_malloc(sizeof(fftw_complex) * heigthRF);
92  m_fft_out_inv = (fftw_complex *)fftw_malloc(sizeof(fftw_complex) * heigthRF);
93 
94  // log compression
95  m_env = new double[heigthRF * widthRF];
96  m_comp = new unsigned char[heigthRF * widthRF];
97 
98  m_signalSize = heigthRF;
99  m_scanLineNumber = widthRF;
100 
101  m_p = fftw_plan_dft_1d(m_signalSize, m_fft_in, m_fft_out, FFTW_FORWARD, FFTW_ESTIMATE);
102  m_pinv = fftw_plan_dft_1d(m_signalSize, m_fft_out, m_fft_out_inv, FFTW_BACKWARD, FFTW_ESTIMATE);
103 
104  m_isInit = true;
105 }
106 
113 void usRFToPreScan2DConverter::enveloppeDetection(const short int *s, double *out)
114 {
116  // const clock_t begin_time = clock();
117  int N = m_signalSize;
118 
119  // Put signal data s into in
120  for (int i = 0; i < N; i++) {
121  m_fft_in[i][0] = (double)s[i];
122  m_fft_in[i][1] = 0.0;
123  }
124 
125  // Obtain the FFT
126  fftw_execute(m_p);
127 
128  for (int i = 0; i < N; i++) {
129 
130  if (i < N / 2)
131  m_fft_out[i][1] = -m_fft_out[i][1];
132  else if (i == N / 2) {
133  m_fft_out[i][0] = 0;
134  m_fft_out[i][1] = 0;
135  } else if (i > N / 2)
136  m_fft_out[i][0] = -m_fft_out[i][0];
137  if (i == 0) {
138  m_fft_out[i][0] = 0;
139  m_fft_out[i][1] = 0;
140  }
141 
142  double z = m_fft_out[i][1];
143  m_fft_out[i][1] = m_fft_out[i][0];
144  m_fft_out[i][0] = z;
145  }
146 
147  // FFT Inverse
148 
149  // Obtain the IFFT
150  fftw_execute(m_pinv);
151  // Put the iFFT output in sa
152  double Ndouble = (double)N;
153  for (int i = 0; i < N; i++) {
154  out[i] = (unsigned char)sqrt(abs(std::complex<double>(m_fft_in[i][0], -m_fft_out_inv[i][0] / Ndouble)));
155  }
156 }
157 
168  usImagePreScan2D<unsigned char> &preScanImage)
169 {
170 
171  if (!m_isInit || ((int)rfImage.getWidth()) != m_scanLineNumber || ((int)rfImage.getHeight()) != m_signalSize) {
172  init(rfImage.getWidth(), rfImage.getHeight());
173  }
174  preScanImage.resize(rfImage.getHeight() / m_decimationFactor, rfImage.getWidth());
175 
176  // First we copy the transducer settings
177  preScanImage.setImagePreScanSettings(rfImage);
178 
179  unsigned int w = rfImage.getWidth();
180  unsigned int h = rfImage.getHeight();
181 
182  unsigned int frameSize = w * h;
183 
184  // prepare signals corresponding to each scanline or column (the bitmap is storing row after row)
185  std::vector<short *> rfSignals;
186  for (unsigned int j = 0; j < w; j++) {
187  short *s = new short[h];
188  for (unsigned int i = 0; i < h; i++) {
189  s[i] = rfImage(i, j);
190  }
191  rfSignals.push_back(s);
192  }
193 
194  // Run envelope detector
195  for (unsigned int i = 0; i < w; ++i) {
196  enveloppeDetection(rfImage.getSignal(i), m_env + i * h);
197  }
198 
199  // Log-compress
200  m_logCompressor.run(m_comp, m_env, frameSize);
201 
202  // find min & max values
203  double min = 1e8;
204  double max = -1e8;
205  for (unsigned int i = 0; i < frameSize; ++i) {
206  if (m_comp[i] < min)
207  min = m_comp[i];
208  if (m_comp[i] > max)
209  max = m_comp[i];
210  }
211  // max-min computation
212  double maxMinDiff = max - min;
213 
214  // Decimate and normalize
215  unsigned int k = 0;
216  for (unsigned int i = 0; i < h; i += m_decimationFactor) {
217  for (unsigned int j = 0; j < w; ++j) {
218  unsigned int vcol = (unsigned int)(((m_comp[i + h * j] - min) / maxMinDiff) * 255);
219  preScanImage[k][j] = (vcol > 255) ? 255 : vcol;
220  }
221  k++;
222  if (k == preScanImage.getHeight()) // prevent overflow for k at last iteration (index given to operator[] on vpImage
223  // goes from 0 to N-1)
224  return;
225  }
226 }
227 
232 int usRFToPreScan2DConverter::getDecimationFactor() { return m_decimationFactor; }
233 
238 void usRFToPreScan2DConverter::setDecimationFactor(int decimationFactor) { m_decimationFactor = decimationFactor; }
239 
240 #endif
void resize(const unsigned int h, const unsigned int w)
void setImagePreScanSettings(const usImagePreScanSettings &preScanSettings)
unsigned int getHeight() const
Definition: usImageRF2D.h:412
const Type * getSignal(unsigned int scanlineIndex) const
Definition: usImageRF2D.h:437
unsigned int getWidth() const
Definition: usImageRF2D.h:424
void run(unsigned char *dst, const double *src, unsigned int size)
void convert(const usImageRF2D< short int > &rfImage, usImagePreScan2D< unsigned char > &preScanImage)
usRFToPreScan2DConverter(int decimationFactor=10)
void setDecimationFactor(int decimationFactor)
Class containing a set of static methods to compute various processes on RF images (gradients,...