UsTK : Ultrasound ToolKit  version 2.0.1 under development (2025-02-01)
tutorial-elastography-BMA-2D.cpp
1 
3 #include <iostream>
4 #include <visp3/ustk_core/usConfig.h>
5 
6 #if (defined(USTK_HAVE_QT5) || defined(USTK_HAVE_VTK_QT)) && (defined(VISP_HAVE_X11) || defined(VISP_HAVE_GDI)) && \
7  defined(USTK_HAVE_FFTW)
8 
9 #include <QApplication>
10 #include <QtCore/QThread>
11 
12 #include <visp3/ustk_core/usRFToPreScan2DConverter.h>
13 #include <visp3/ustk_elastography/usElastography.h>
14 #include <visp3/ustk_grabber/usNetworkGrabberRF2D.h>
15 
16 #include <visp3/gui/vpDisplayGDI.h>
17 #include <visp3/gui/vpDisplayX.h>
18 #include <visp3/io/vpImageIo.h>
19 
20 int main(int argc, char **argv)
21 {
22  // QT application
23  QApplication app(argc, argv);
24 
25  usElastography *elastography = new usElastography;
27  elastography->setROI(40, 3200, 50, 500);
28 
29  usNetworkGrabberRF2D *qtGrabber = new usNetworkGrabberRF2D();
30  qtGrabber->setIPAddress("127.0.0.1");
31  qtGrabber->connectToServer();
32 
33  // setting acquisition parameters
35  header.probeId = 15; // 4DC7 id = 15
36  header.slotId = 0; // top slot id = 0
37  header.imagingMode = 12; // B-mode = 0, RF = 12
38 
39  // prepare image;
42 
43  usRFToPreScan2DConverter converter;
44 
45  // prepare converter
46  vpImage<unsigned char> strainImage;
47 
48 // Prepare display
49 #if defined(VISP_HAVE_X11)
50  vpDisplayX *displayEcho = NULL;
51  vpDisplayX *displayElas = NULL;
52 #elif defined(VISP_HAVE_GDI)
53  vpDisplayGDI *displayEcho = NULL;
54  vpDisplayGDI *displayElas = NULL;
55 #endif
56  bool displayInit = false;
57  bool captureRunning = true;
58 
59  // sending acquisition parameters
60  qtGrabber->initAcquisition(header);
61 
62  qtGrabber->runAcquisition();
63 
64  std::cout << "waiting ultrasound initialisation..." << std::endl;
65 
66  // our local grabbing loop
67  do {
68  if (qtGrabber->isFirstFrameAvailable()) {
69  grabbedFrame = qtGrabber->acquire();
70 
71  std::cout << "MAIN THREAD received frame No : " << grabbedFrame->getFrameCount() << std::endl;
72 
73  elastography->updateRF(*grabbedFrame);
74  strainImage = elastography->run();
75 
76  converter.convert(*grabbedFrame, preScanImage);
77 
78  // init display
79  if (!displayInit && strainImage.getHeight() != 0 && strainImage.getWidth() != 0) {
80 #if defined(VISP_HAVE_X11)
81  displayEcho = new vpDisplayX(preScanImage);
82  displayElas = new vpDisplayX(strainImage);
83 #elif defined(VISP_HAVE_GDI)
84  displayEcho = new vpDisplayGDI(preScanImage);
85  displayElas = new vpDisplayGDI(strainImage);
86 #endif
87  displayInit = true;
88  }
89 
90  // processing display
91  if (displayInit) {
92  vpDisplay::display(preScanImage);
93  vpDisplay::displayRectangle(preScanImage, 320, 40, 50, 50, vpColor::red);
94  vpDisplay::display(strainImage);
95  vpDisplay::flush(preScanImage);
96  vpDisplay::flush(strainImage);
97  }
98  }
99 
100  else {
101  vpTime::wait(10);
102  }
103  } while (captureRunning);
104 
105  if (displayInit) {
106  delete displayElas;
107  delete displayEcho;
108  }
109  return app.exec();
110 }
111 
112 #else
113 int main()
114 {
115  std::cout << "You should intall Qt5 (with wigdets and network modules), FFTW and a display graphic system (GDI or "
116  "X11) to run this tutorial"
117  << std::endl;
118  return 0;
119 }
120 
121 #endif
Computation of a strain map using two sucessive RF images acquired at different compressions of the p...
vpImage< unsigned char > run()
void setROI(int tx, int ty, int tw, int th)
void updateRF(const usImageRF2D< short int > &image)
void setMotionEstimator(MotionEstimator t_mest)
Class to store additionnal informations arriving on the network with ultrasound images grabbed,...
quint32 getFrameCount() const
Specific class to grab RF frames from the ultrasound station on the network.
usFrameGrabbedInfo< usImageRF2D< short int > > * acquire()
bool initAcquisition(const usNetworkGrabber::usInitHeaderSent &header)
void setIPAddress(const std::string &s_ip)
2D conversion from RF signal to pre-scan image
void convert(const usImageRF2D< short int > &rfImage, usImagePreScan2D< unsigned char > &preScanImage)