UsTK : Ultrasound ToolKit  version 2.0.1 under development (2023-12-07)
usRFToPreScan3DConverter.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/usRFToPreScan3DConverter.h>
40 
41 #if defined(USTK_HAVE_FFTW)
42 
43 #ifdef VISP_HAVE_OPENMP
44 #include <omp.h>
45 #endif
46 
50 usRFToPreScan3DConverter::usRFToPreScan3DConverter() : m_frameNumber(), m_isInit(false), m_decimationFactor(10) {}
51 
56 
67  usImagePreScan3D<unsigned char> &preScanImage)
68 {
69  if (!m_isInit || (((int)rfImage.getNumberOfFrames()) != m_frameNumber) ||
70  (((int)rfImage.getHeight()) != m_heightRF) || (((int)rfImage.getWidth()) != m_widthRF)) {
71  init(rfImage.getHeight(), rfImage.getWidth(), rfImage.getNumberOfFrames());
72  }
73  preScanImage.resize(rfImage.getHeight() / getDecimationFactor(), rfImage.getWidth(), rfImage.getNumberOfFrames());
74  // First we copy the transducer/motor settings
75  preScanImage.setImagePreScanSettings(rfImage);
76  preScanImage.setAxialResolution(rfImage.getDepth() / preScanImage.getHeight());
77  preScanImage.setMotorSettings(rfImage);
78  std::vector<usImagePreScan2D<unsigned char> > preScanFrame;
79  preScanFrame.resize(m_frameNumber);
80  std::vector<usImageRF2D<short int> > frameRF;
81  frameRF.resize(m_frameNumber);
82  // loop to convert each frame of the volume
83 
84  for (int i = 0; i < m_frameNumber; i++) {
85  rfImage.getFrame(frameRF.at(i), i);
86  }
87 
88  if (!m_isInit) {
89  init(rfImage.getHeight(), rfImage.getWidth(), rfImage.getNumberOfFrames());
90  }
91 
92 #ifdef VISP_HAVE_OPENMP
93 #pragma omp parallel for
94 #endif
95  for (int i = 0; i < m_frameNumber; i++) {
96  m_converter[i].convert(frameRF[i], preScanFrame.at(i));
97  }
98 
99  for (int i = 0; i < m_frameNumber; i++) {
100  preScanImage.insertFrame(preScanFrame[i], i);
101  }
102 }
103 
108 int usRFToPreScan3DConverter::getDecimationFactor() { return m_decimationFactor; }
109 
115 {
116  for (int i = 0; i < m_frameNumber; i++)
117  m_converter[i].setDecimationFactor(decimationFactor);
118  m_decimationFactor = decimationFactor;
119 }
120 
124 void usRFToPreScan3DConverter::init(int heightRF, int widthRF, int frameNumber)
125 {
126  // first init
127  if (!m_isInit) {
128  m_converter = new usRFToPreScan2DConverter[frameNumber];
129  for (int i = 0; i < frameNumber; i++) {
130  m_converter[i].init(widthRF, heightRF);
131  m_converter[i].setDecimationFactor(m_decimationFactor);
132  }
133 
134  m_isInit = true;
135  m_frameNumber = frameNumber;
136  m_heightRF = heightRF;
137  m_widthRF = widthRF;
138  }
139  // update image size
140  else if (m_frameNumber != frameNumber || m_heightRF != heightRF || m_widthRF != widthRF) {
141  for (int i = 0; i < m_frameNumber; i++)
142  delete &(m_converter[i]);
143 
144  m_frameNumber = frameNumber;
145  m_heightRF = heightRF;
146  m_widthRF = widthRF;
147  m_converter = new usRFToPreScan2DConverter[m_frameNumber];
148  for (int i = 0; i < m_frameNumber; i++) {
149  m_converter[i].init(m_widthRF, m_heightRF);
150  m_converter[i].setDecimationFactor(m_decimationFactor);
151  }
152  }
153 }
154 #endif
unsigned int getHeight() const
Definition: usImage3D.h:131
void resize(unsigned int height, unsigned int width, unsigned int numberOfFrames)
void insertFrame(const usImagePreScan2D< Type > &frame, unsigned int index)
void setImagePreScanSettings(const usImagePreScanSettings &preScanSettings)
void setAxialResolution(const double axialResolution)
void getFrame(usImageRF2D< Type > &image, unsigned int index) const
Definition: usImageRF3D.h:395
unsigned int getWidth() const
Definition: usImageRF3D.h:471
unsigned int getHeight() const
Definition: usImageRF3D.h:477
unsigned int getNumberOfFrames() const
Definition: usImageRF3D.h:483
void setMotorSettings(const usMotorSettings &other)
2D conversion from RF signal to pre-scan image
void convert(const usImageRF2D< short int > &rfImage, usImagePreScan2D< unsigned char > &preScanImage)
void setDecimationFactor(int decimationFactor)
void init(int heightRF, int widthRF, int frameNumber)
void setDecimationFactor(int decimationFactor)
void convert(const usImageRF3D< short int > &rfImage, usImagePreScan3D< unsigned char > &preScanImage)