QTrk
Classes | Public Member Functions | Private Member Functions | Static Private Member Functions | Private Attributes | List of all members
QueuedCPUTracker Class Reference

CPU implementation of the QueuedTracker interface. More...

#include <QueuedCPUTracker.h>

Inheritance diagram for QueuedCPUTracker:
QueuedTracker

Classes

struct  Job
 Structure around LocalizationJob to link data and job. More...
 
struct  Thread
 Structure to facilitate the use of threads. More...
 

Public Member Functions

 QueuedCPUTracker (const QTrkComputedConfig &cc)
 Create a QueuedCPUTracker instance with specified configuration. More...
 
 ~QueuedCPUTracker ()
 Delete the instance and free related memory. More...
 
void Start ()
 Start up the threads and initiate their worker function. More...
 
void Break (bool pause)
 Enable or disable handling of jobs. More...
 
void GenerateTestImage (float *dst, float xp, float yp, float size, float SNratio)
 Generate a test image based on a the point spread function. More...
 
int NumThreads ()
 Get the used number of threads. More...
 
void SetLocalizationMode (LocMode_t lt) override
 Select which algorithm is to be used. More...
 
void SetRadialZLUT (float *data, int num_zluts, int planes) override
 Set the radial lookup tables to be used for z tracking. More...
 
void GetRadialZLUT (float *zlut) override
 Get the radial lookup tables used for z tracking. More...
 
void GetRadialZLUTSize (int &count, int &planes, int &rsteps) override
 Get the dimensions of the radial lookup table data. More...
 
void SetRadialWeights (float *rweights) override
 Set radial weights used for comparing LUT profiles. More...
 
void SetRadialWeights (std::vector< float > weights)
 Overload of SetRadialWeights to use a vector. More...
 
void ScheduleLocalization (void *data, int pitch, QTRK_PixelDataType pdt, const LocalizationJob *jobInfo) override
 Add a job to the queue to be processed. A job entails running the required algorithms on a single region of interest. More...
 
void EnableRadialZLUTCompareProfile (bool enabled)
 Set a flag to enable saving of error curves. More...
 
void GetRadialZLUTCompareProfile (float *dst)
 Get saved error curve. More...
 
void ZLUTSelfTest ()
 Test the planes saved in the lookup table against itself. More...
 
void BeginLUT (uint flags)
 Setup to begin building a lookup table. More...
 
void BuildLUT (void *data, int pitch, QTRK_PixelDataType pdt, int plane, vector2f *known_pos=0) override
 Add a new lookup table plane. More...
 
void FinalizeLUT () override
 Finalize the lookup tables in memory. More...
 
void GetImageZLUTSize (int *dims)
 Get the dimensions of the image lookup table data. More...
 
void GetImageZLUT (float *dst)
 Get the image lookup tables used. More...
 
bool SetImageZLUT (float *dst, float *radial_zlut, int *dims)
 Set the image lookup tables to be used for z tracking. More...
 
void SetPixelCalibrationImages (float *offset, float *gain) override
 Set the pixel calibration images. More...
 
void SetPixelCalibrationFactors (float offsetFactor, float gainFactor) override
 Set the pixel calibration factors. More...
 
int GetQueueLength (int *maxQueueLength=0) override
 Get the lengths of the queue of jobs to be handled. More...
 
int FetchResults (LocalizationResult *dstResult, int maxResults) override
 Fetch available results. More...
 
void ClearResults () override
 Clear results. More...
 
void Flush () override
 Stop waiting for more jobs to do, and just process the current batch. More...
 
bool IsIdle () override
 Test to see if the tracker is idle. More...
 
int GetResultCount () override
 Get the number of finished localization jobs (=results) available in memory. More...
 
bool GetDebugImage (int id, int *w, int *h, float **data)
 Get the debug image of a tracker. More...
 
ConfigValueMap GetConfigValues () override
 Get the used additional configurations. More...
 
void SetConfigValue (std::string name, std::string value) override
 Set an additional setting. More...
 
std::string GetProfileReport ()
 Get the output of performance profiling. More...
 
float * GetZLUTByIndex (int index)
 Get a pointer to the starting pixel of the specified LUT. More...
 
- Public Member Functions inherited from QueuedTracker
 QueuedTracker ()
 
virtual ~QueuedTracker ()
 
void ScheduleImageData (ImageData *data, const LocalizationJob *jobInfo)
 Quick function to schedule a single ROI from an ImageData object. More...
 
virtual int ScheduleFrame (void *imgptr, int pitch, int width, int height, ROIPosition *positions, int numROI, QTRK_PixelDataType pdt, const LocalizationJob *jobInfo)
 Schedule an entire frame at once, allowing for further optimizations. More...
 
virtual std::string GetWarnings ()
 Get a report of encountered errors. More...
 
ImageData DebugImage (int ID)
 Get the debug image as an ImageData object. More...
 
void ScheduleLocalization (uchar *data, int pitch, QTRK_PixelDataType pdt, uint frame, uint timestamp, vector3f *initial, uint zlutIndex)
 Add an image to the queue to be processed. Creates a job. More...
 
void ComputeZBiasCorrection (int bias_planes, CImageData *result, int smpPerPixel, bool useSplineInterp)
 
float ZLUTBiasCorrection (float z, int zlut_planes, int bead)
 
void SetZLUTBiasCorrection (const CImageData &data)
 
CImageDataGetZLUTBiasCorrection ()
 

Private Member Functions

void UpdateZLUTs ()
 Update the ZLUTs for all threads. More...
 
int ImageLUTNumBeads ()
 Return the number of beads in the image LUT. More...
 
int ImageLUTWidth ()
 Return the height of the image LUT. More...
 
int ImageLUTHeight ()
 Return the width of the image LUT. More...
 
float * GetImageLUTByIndex (int index, int plane=0)
 Get a specific image LUT or plane from an image LUT. More...
 
void JobFinished (Job *j)
 Flag job memory for reuse. More...
 
JobGetNextJob ()
 Get the next job in the queue. More...
 
JobAllocateJob ()
 Allocate job memory. More...
 
void AddJob (Job *j)
 Add a job to the queue. More...
 
void ProcessJob (Thread *th, Job *j)
 Process a job. More...
 
void SetTrackerImage (CPUTracker *trk, Job *job)
 Copy the ROI to the tracker's memory. More...
 
void ApplyOffsetGain (CPUTracker *trk, int beadIndex)
 Calibrate an image. More...
 

Static Private Member Functions

static void WorkerThreadMain (void *arg)
 The loop executed by the threads. More...
 

Private Attributes

LocMode_t localizeMode
 Localization settings. More...
 
Threads::Mutex jobs_mutex
 Mutex for the job queue. More...
 
Threads::Mutex jobs_buffer_mutex
 Mutex for the job buffer. More...
 
Threads::Mutex results_mutex
 Mutex for the results. More...
 
std::deque< Job * > jobs
 Queue to hold the jobs. More...
 
int jobCount
 Number of jobs in the queue. More...
 
std::vector< Job * > jobs_buffer
 Stores memory. Enables reuse of allocated Job memory after a job has been completed. More...
 
std::deque< LocalizationResultresults
 Queue to store the localization results. More...
 
int resultCount
 Number of results available. More...
 
int maxQueueSize
 Maximum number of jobs in the queue. More...
 
int jobsInProgress
 Number of jobs currently being processed. More...
 
Threads::Mutex gc_mutex
 Image calibration mutex. More...
 
float * calib_gain
 Vector with gains used for pixel correction. More...
 
float * calib_offset
 Vector with offsets used for pixel correction. More...
 
float gc_gainFactor
 Factor by which to scale the gain. More...
 
float gc_offsetFactor
 Factor by which to scale the pixel calibration offset. More...
 
int downsampleWidth
 Width of ROIs after downsampling. More...
 
int downsampleHeight
 Height of ROIs after downsampling. More...
 
std::vector< Threadthreads
 Vector with active threads. More...
 
float * zluts
 Pointer to the first pixel of the Z lookup tables. More...
 
float * zlut_cmpprofiles
 Array in which to save errorcurves if enabled with EnableRadialZLUTCompareProfile. More...
 
bool zlut_enablecmpprof
 Flag to save errorcurves. See EnableRadialZLUTCompareProfile. More...
 
int zlut_count
 Number of ZLUTs (= number of beads). More...
 
int zlut_planes
 Number of planes per ZLUT. More...
 
uint zlut_buildflags
 ZLUT build flags, set through BeginLUT. More...
 
std::vector< float > zcmp
 Scaling factors for the ZLUT algorithm. More...
 
std::vector< float > qi_radialbinweights
 Scaling factors for the QI algorithm. More...
 
int image_lut_dims [4]
 Image LUT dimensions, 4D. [numBeads][numPlanes][height][width]. More...
 
int image_lut_nElem_per_bead
 Image LUT number of pixels in LUT per bead. Is Width * Height * Planes. More...
 
float * image_lut
 Image LUT data pointer. More...
 
float * image_lut_dz
 Image LUT first derivative. (???) More...
 
float * image_lut_dz2
 Image LUT second derivative. (???) More...
 
bool quitWork
 Signal threads to stop their work. More...
 
bool processJobs
 Flag for threads to continue processing jobs. More...
 
bool dbgPrintResults
 Flag to enable easy printing of intermediate data. More...
 

Additional Inherited Members

- Public Types inherited from QueuedTracker
typedef std::map< std::string, std::string > ConfigValueMap
 Datastructure used to return additional settings in a string-string key-value pair mapping. More...
 
- Public Attributes inherited from QueuedTracker
QTrkComputedConfig cfg
 The settings used by this instance of QueuedTracker. More...
 
- Protected Attributes inherited from QueuedTracker
CImageDatazlut_bias_correction
 

Detailed Description

CPU implementation of the QueuedTracker interface.

Creates and maintains multiple threads (one per available core by default) and schedules localizations over them. Each thread has its own CPUTracker instance which provides the actual calculations.

Definition at line 16 of file QueuedCPUTracker.h.

Constructor & Destructor Documentation

§ QueuedCPUTracker()

QueuedCPUTracker::QueuedCPUTracker ( const QTrkComputedConfig cc)

Create a QueuedCPUTracker instance with specified configuration.

Parameters
[in]ccComputed config to use.
Returns
A QueuedCPUTracker instance.

Copy the configuration.

If QTrkComputedConfig::numThreads is not specified (<0), select one thread per CPU.

Queue length is 250 per thread, min 500.

Calculate the radial profile weights from the amount of steps.

Initialize to only COM localizations.

Invoke the Start function.

Definition at line 96 of file QueuedCPUTracker.cpp.

97  : jobs_mutex("jobs"), jobs_buffer_mutex("jobs_buffer"), results_mutex("results")
98 {
99  if(0){//cc.testRun){
100  std::string folder = GetCurrentOutputPath();
101  if(GetFileAttributesA(folder.c_str()) & FILE_ATTRIBUTE_DIRECTORY)
102  CreateDirectory((LPCTSTR)folder.c_str(),NULL);
103  }
104 
106  cfg = cc;
107  quitWork = false;
108 
110  if (cfg.numThreads < 0) {
111  cfg.numThreads = Threads::GetCPUCount();
112  dbgprintf("Using %d threads\n", cfg.numThreads);
113  }
114 
116  maxQueueSize = std::max(2,cfg.numThreads) * 250;
117  jobCount = 0;
118  resultCount = 0;
119 
120  zlut_cmpprofiles = 0;
121  zlut_enablecmpprof = false;
122  zluts = 0;
123  zlut_count = zlut_planes = 0;
124  processJobs = false;
125  jobsInProgress = 0;
126  dbgPrintResults = false;
127 
129  qi_radialbinweights = ComputeRadialBinWindow(cfg.qi_radialsteps);
130 
131  calib_gain = calib_offset = 0;
135 
136  for (int i=0;i<4;i++) image_lut_dims[i]=0;
138  image_lut = 0;
140 
141  downsampleWidth = cfg.width >> cfg.downsample;
142  downsampleHeight = cfg.height >> cfg.downsample;
143 
145  Start();
146 }
int jobsInProgress
Number of jobs currently being processed.
Threads::Mutex jobs_mutex
Mutex for the job queue.
float gc_gainFactor
Factor by which to scale the gain.
float * image_lut
Image LUT data pointer.
std::vector< float > qi_radialbinweights
Scaling factors for the QI algorithm.
Threads::Mutex jobs_buffer_mutex
Mutex for the job buffer.
float * zlut_cmpprofiles
Array in which to save errorcurves if enabled with EnableRadialZLUTCompareProfile.
int zlut_count
Number of ZLUTs (= number of beads).
void Start()
Start up the threads and initiate their worker function.
int maxQueueSize
Maximum number of jobs in the queue.
bool processJobs
Flag for threads to continue processing jobs.
float * calib_gain
Vector with gains used for pixel correction.
int jobCount
Number of jobs in the queue.
LocMode_t localizeMode
Localization settings.
Threads::Mutex results_mutex
Mutex for the results.
float * image_lut_dz
Image LUT first derivative. (???)
float * zluts
Pointer to the first pixel of the Z lookup tables.
std::vector< float > ComputeRadialBinWindow(int rsteps)
Calculate the radial weights for ZLUT profile comparisons.
Definition: utils.cpp:706
int resultCount
Number of results available.
std::string GetCurrentOutputPath(bool ext)
Definition: utils.cpp:19
static int GetCPUCount()
Definition: threads.h:138
QTrkComputedConfig cfg
The settings used by this instance of QueuedTracker.
bool quitWork
Signal threads to stop their work.
int zlut_planes
Number of planes per ZLUT.
int image_lut_dims[4]
Image LUT dimensions, 4D. [numBeads][numPlanes][height][width].
float * image_lut_dz2
Image LUT second derivative. (???)
void dbgprintf(const char *fmt,...)
Definition: utils.cpp:149
Use only COM.
Definition: qtrk_c_api.h:7
float * calib_offset
Vector with offsets used for pixel correction.
bool zlut_enablecmpprof
Flag to save errorcurves. See EnableRadialZLUTCompareProfile.
int downsampleHeight
Height of ROIs after downsampling.
int downsampleWidth
Width of ROIs after downsampling.
int image_lut_nElem_per_bead
Image LUT number of pixels in LUT per bead. Is Width * Height * Planes.
float gc_offsetFactor
Factor by which to scale the pixel calibration offset.
bool dbgPrintResults
Flag to enable easy printing of intermediate data.

§ ~QueuedCPUTracker()

QueuedCPUTracker::~QueuedCPUTracker ( )

Delete the instance and free related memory.

Definition at line 148 of file QueuedCPUTracker.cpp.

149 {
150  // wait for threads to finish
151  quitWork = true;
152 
153  for (int k=0;k<threads.size();k++) {
154  delete threads[k].mutex;
155  Threads::WaitAndClose(threads[k].thread);
156  delete threads[k].tracker;
157  }
158 
159  // free job memory
162 
163  delete[] calib_gain;
164  delete[] calib_offset;
165  if (zluts) delete[] zluts;
166  if (zlut_cmpprofiles) delete[] zlut_cmpprofiles;
167  if (image_lut) delete[] image_lut;
168  if (image_lut_dz) delete[] image_lut_dz;
169  if (image_lut_dz2) delete[] image_lut_dz2;
170 }
float * image_lut
Image LUT data pointer.
float * zlut_cmpprofiles
Array in which to save errorcurves if enabled with EnableRadialZLUTCompareProfile.
std::vector< Job * > jobs_buffer
Stores memory. Enables reuse of allocated Job memory after a job has been completed.
float * calib_gain
Vector with gains used for pixel correction.
float * image_lut_dz
Image LUT first derivative. (???)
float * zluts
Pointer to the first pixel of the Z lookup tables.
std::deque< Job * > jobs
Queue to hold the jobs.
void DeleteAllElems(T &c)
Definition: utils.h:19
static void WaitAndClose(Handle *h)
Definition: threads.h:128
bool quitWork
Signal threads to stop their work.
std::vector< Thread > threads
Vector with active threads.
float * image_lut_dz2
Image LUT second derivative. (???)
float * calib_offset
Vector with offsets used for pixel correction.

Member Function Documentation

§ AddJob()

void QueuedCPUTracker::AddJob ( Job j)
private

Add a job to the queue.

Parameters
[in]jThe job to add.

Definition at line 75 of file QueuedCPUTracker.cpp.

76 {
77  jobs_mutex.lock();
78  jobs.push_back(j);
79  jobCount++;
81 }
Threads::Mutex jobs_mutex
Mutex for the job queue.
void lock()
Definition: threads.h:74
int jobCount
Number of jobs in the queue.
std::deque< Job * > jobs
Queue to hold the jobs.
void unlock()
Definition: threads.h:79

§ AllocateJob()

QueuedCPUTracker::Job * QueuedCPUTracker::AllocateJob ( )
private

Allocate job memory.

Tries to re-use memory of finished jobs before allocating new memory.

Returns
Pointer to a usable job instance.

Definition at line 62 of file QueuedCPUTracker.cpp.

63 {
66  if (!jobs_buffer.empty()) {
67  j = jobs_buffer.back();
68  jobs_buffer.pop_back();
69  } else
70  j = new Job;
72  return j;
73 }
Threads::Mutex jobs_buffer_mutex
Mutex for the job buffer.
std::vector< Job * > jobs_buffer
Stores memory. Enables reuse of allocated Job memory after a job has been completed.
void lock()
Definition: threads.h:74
Structure around LocalizationJob to link data and job.
void unlock()
Definition: threads.h:79

§ ApplyOffsetGain()

void QueuedCPUTracker::ApplyOffsetGain ( CPUTracker trk,
int  beadIndex 
)
private

Calibrate an image.

Calibrates the image set on a tracker through SetTrackerImage.

Parameters
[in]trkThe tracker whose image to calibrate.
[in]beadIndexThe bead this image belongs to.

Definition at line 265 of file QueuedCPUTracker.cpp.

266 {
267  if (calib_offset || calib_gain) {
268  int index = cfg.width*cfg.height*beadIndex;
269 
270  //gc_mutex.lock();
271  //float gf = gc_gainFactor, of = gc_offsetFactor;
272  //gc_mutex.unlock();
273  float gf=1,of=1;
274 
275  trk->ApplyOffsetGain(calib_offset ? &calib_offset[index] : 0, calib_gain ? &calib_gain[index] : 0, of, gf);
276 // if (j->job.frame%100==0)
277  }
278 }
int width
Width of regions of interest to be handled. Typically equals height (square ROI). ...
Definition: qtrk_c_api.h:106
int height
Height of regions of interest to be handled. Typically equals width (square ROI). ...
Definition: qtrk_c_api.h:107
float * calib_gain
Vector with gains used for pixel correction.
void ApplyOffsetGain(float *offset, float *gain, float offsetFactor, float gainFactor)
Scale the input image with the background calibration values.
Definition: cpu_tracker.cpp:95
QTrkComputedConfig cfg
The settings used by this instance of QueuedTracker.
float * calib_offset
Vector with offsets used for pixel correction.

§ BeginLUT()

void QueuedCPUTracker::BeginLUT ( uint  flags)
virtual

Setup to begin building a lookup table.

Sets the flags by which the lookup table is built. Flags are defined in a uint bitmask format as:

Name Value Description
0 Default, no special flags.
BUILDLUT_IMAGELUT 1 (2^0) Build an image LUT. An image LUT seems to be a work in progress to save ROIs rather than profiles in the LUT.
BUILDLUT_FOURIER 2 (2^1) Build a fourier LUT.
BUILDLUT_NORMALIZE 4 (2^2) Normalize radial profiles. Irrelevant, since FinalizeLUT always normalizes.
BUILDLUT_BIASCORRECT8 (2^3) Enable bias correction. Currently only partly implemented, don't use.
Parameters
[in]flagsUINT interpreted as a series of bits to set settings.

Implements QueuedTracker.

Definition at line 584 of file QueuedCPUTracker.cpp.

585 {
586  zlut_buildflags = flags;
587 }
uint zlut_buildflags
ZLUT build flags, set through BeginLUT.

§ Break()

void QueuedCPUTracker::Break ( bool  pause)

Enable or disable handling of jobs.

Parameters
[in]pauseTrue stops new jobs from being started.

Definition at line 214 of file QueuedCPUTracker.cpp.

215 {
216  processJobs = !brk;
217 }
bool processJobs
Flag for threads to continue processing jobs.

§ BuildLUT()

void QueuedCPUTracker::BuildLUT ( void *  data,
int  pitch,
QTRK_PixelDataType  pdt,
int  plane,
vector2f known_pos = 0 
)
overridevirtual

Add a new lookup table plane.

Takes a stack of ROI images through data. Determines and adds the profile for each ROI to its respective LUT.

Parameters
[in]dataPointer to the start of an image stack.
[in]pitchWidth of the data in memory so offsets can be calculated.
[in]pdtPixel data type for the data. See QTRK_PixelDataType.
[in]planeThe plane number the ROIs are taken for.
[in]known_posCenter position from which to start making the radial profile. A standard QI localization with applied settings is performed if 0 (default).

Implements QueuedTracker.

Definition at line 589 of file QueuedCPUTracker.cpp.

590 {
591  parallel_for(zlut_count,[&] (int i) {
592 // for(int i=0;i<zlut_count;i++) {
593  CPUTracker trk (cfg.width,cfg.height);
594  void *img_data = (uchar*)data + pitch * cfg.height * i;
595 
596  if (pdt == QTrkFloat) {
597  trk.SetImage((float*)img_data, pitch);
598  } else if (pdt == QTrkU8) {
599  trk.SetImage8Bit((uchar*)img_data, pitch);
600  } else {
601  trk.SetImage16Bit((ushort*)img_data,pitch);
602  }
603  ApplyOffsetGain(&trk, i);
604 
605  vector2f pos;
606 
607  if (known_pos) {
608  pos = known_pos[i];
609  } else {
610  vector2f com = trk.ComputeMeanAndCOM();
611  bool bhit;
613  //dbgprintf("BuildLUT() COMPos: %f,%f, QIPos: x=%f, y=%f\n", com.x,com.y, pos.x, pos.y);
614  }
616  int h=ImageLUTHeight(), w=ImageLUTWidth();
617  float* lut_dst = &image_lut[ i * image_lut_nElem_per_bead + w*h* plane ];
618 
619  vector2f ilut_scale(1,1);
620  float startx = pos.x - w/2*ilut_scale.x;
621  float starty = pos.y - h/2*ilut_scale.y;
622 
623  for (int y=0;y<h;y++) {
624  for (int x=0;x<w;x++) {
625 
626  float px = startx + x*ilut_scale.x;
627  float py = starty + y*ilut_scale.y;
628 
629  bool outside=false;
630  float v = Interpolate(trk.srcImage, trk.width, trk.height, px, py, &outside);
631  lut_dst[y*w+x] += v - trk.mean;
632  }
633  }
634  }
635 
636  float *bead_zlut=GetZLUTByIndex(i);
637  float *tmp = new float[cfg.zlut_radialsteps];
638 
640  trk.FourierRadialProfile(tmp, cfg.zlut_radialsteps, cfg.zlut_angularsteps, cfg.zlut_minradius, cfg.zlut_maxradius);
641  if (plane==0) {
642  for (int i=0;i<trk.width*trk.height;i++)
643  trk.srcImage[i]=sqrtf(trk.srcImage[i]);
644  trk.SaveImage("freqimg.jpg");
645  }
646  } else {
647  trk.ComputeRadialProfile(tmp, cfg.zlut_radialsteps, cfg.zlut_angularsteps, cfg.zlut_minradius, cfg.zlut_maxradius, pos, false);
648  }
649  // WriteArrayAsCSVRow("rlut-test.csv", tmp, cfg.zlut_radialsteps, plane>0);
650  for(int i=0;i<cfg.zlut_radialsteps;i++)
651  bead_zlut[plane*cfg.zlut_radialsteps+i] += tmp[i];
652  delete[] tmp;
653  } );
654 }
float zlut_minradius
Distance in pixels from the bead center from which to start sampling profiles. Default 1...
Definition: qtrk_c_api.h:135
void ApplyOffsetGain(CPUTracker *trk, int beadIndex)
Calibrate an image.
64 bit float
Definition: qtrk_c_api.h:37
float * image_lut
Image LUT data pointer.
#define BUILDLUT_FOURIER
float qi_minradius
Distance in pixels from the bead center from which to start sampling profiles. Default 1...
Definition: qtrk_c_api.h:141
int width
Width of regions of interest to be handled. Typically equals height (square ROI). ...
Definition: qtrk_c_api.h:106
int zlut_count
Number of ZLUTs (= number of beads).
Class with all CPU algorithm implementations.
Definition: cpu_tracker.h:38
int height
Height of regions of interest to be handled. Typically equals width (square ROI). ...
Definition: qtrk_c_api.h:107
unsigned short ushort
Definition: std_incl.h:128
8 bit unsigned int
Definition: qtrk_c_api.h:35
int zlut_radialsteps
Number of radial steps to sample on.
Definition: qtrk_c_api.h:198
void parallel_for(int count, TF f)
Definition: threads.h:251
float qi_maxradius
Max radius in pixels of the sampling circle.
Definition: qtrk_c_api.h:204
float * GetZLUTByIndex(int index)
Get a pointer to the starting pixel of the specified LUT.
int ImageLUTWidth()
Return the height of the image LUT.
QTrkComputedConfig cfg
The settings used by this instance of QueuedTracker.
uint zlut_buildflags
ZLUT build flags, set through BeginLUT.
int qi_radialsteps
Number of radial steps to sample on.
Definition: qtrk_c_api.h:202
int ImageLUTHeight()
Return the width of the image LUT.
T Interpolate(T *image, int width, int height, float x, float y, bool *outside=0)
Definition: utils.h:43
int qi_angstepspq
Number of angular steps to sample on per quadrant.
Definition: qtrk_c_api.h:203
int qi_iterations
Number of times to run the QI algorithm, sampling around the last found position. ...
Definition: qtrk_c_api.h:140
float qi_angstep_factor
Factor to reduce angular steps on lower iterations. Default 1.0 (no effect).
Definition: qtrk_c_api.h:157
int image_lut_nElem_per_bead
Image LUT number of pixels in LUT per bead. Is Width * Height * Planes.
#define BUILDLUT_IMAGELUT
int zlut_angularsteps
Number of angular steps to sample on.
Definition: qtrk_c_api.h:199
float zlut_maxradius
Max radius in pixels of the sampling circle.
Definition: qtrk_c_api.h:200
unsigned char uchar
Definition: std_incl.h:130

§ ClearResults()

void QueuedCPUTracker::ClearResults ( )
overridevirtual

Clear results.

Implements QueuedTracker.

Definition at line 29 of file QueuedCPUTracker.cpp.

30 {
32  resultCount = 0;
33  results.clear();
35 }
void lock()
Definition: threads.h:74
Threads::Mutex results_mutex
Mutex for the results.
std::deque< LocalizationResult > results
Queue to store the localization results.
int resultCount
Number of results available.
void unlock()
Definition: threads.h:79

§ EnableRadialZLUTCompareProfile()

void QueuedCPUTracker::EnableRadialZLUTCompareProfile ( bool  enabled)
virtual

Set a flag to enable saving of error curves.

Errors obtained by comparing a radial profile to a ZLUT will be kept in memory rather than destroyed. Only saves for one localization. Error curve can be retreived by GetRadialZLUTCompareProfile.

Note
Not implemented for CUDA.
Parameters
[in]enabledFlag (boolean) to save error curves. Default is false.

Implements QueuedTracker.

Definition at line 544 of file QueuedCPUTracker.cpp.

545 {
546  zlut_enablecmpprof=enabled;
547 }
bool zlut_enablecmpprof
Flag to save errorcurves. See EnableRadialZLUTCompareProfile.

§ FetchResults()

int QueuedCPUTracker::FetchResults ( LocalizationResult results,
int  maxResults 
)
overridevirtual

Fetch available results.

Note
Removes results from internal QueuedTracker memory.
Parameters
[in]resultsArray of pre-allocated LocalizationResult to which to add the results.
[in]maxResultsMaximum number of results to fetch. Corresponds to maximum size of dstResult.
Returns
Number of fetched results.

Implements QueuedTracker.

Definition at line 473 of file QueuedCPUTracker.cpp.

474 {
475  int numResults = 0;
477  while (numResults < maxResults && !results.empty()) {
478  dstResults[numResults++] = results.back();
479  results.pop_back();
480  resultCount--;
481  }
483  return numResults;
484 }
void lock()
Definition: threads.h:74
Threads::Mutex results_mutex
Mutex for the results.
int resultCount
Number of results available.
void unlock()
Definition: threads.h:79

§ FinalizeLUT()

void QueuedCPUTracker::FinalizeLUT ( )
overridevirtual

Finalize the lookup tables in memory.

Normalizes the profiles for radial lookup tables and calculates derivates and adds boundary conditions for image LUTs.

Implements QueuedTracker.

Definition at line 656 of file QueuedCPUTracker.cpp.

657 {
658  // normalize radial LUT?
659 
660  if (zluts) {
661  for (int i=0;i<zlut_count*zlut_planes;i++) {
662 
663  // WriteArrayAsCSVRow("finalize-lut.csv", &zluts[cfg.zlut_radialsteps*i], cfg.zlut_radialsteps, i>0);
665 
666  }
667  }
668 
669  int w = ImageLUTWidth();
670  int h = ImageLUTHeight();
671 
672  if (w * h > 0) {
673 
676 
677  // Compute 1st and 2nd order derivatives
678  for (int i=0;i<image_lut_dims[0];i++) {
679  for (int z=1;z<image_lut_dims[1]-1;z++) {
680  float *img = &image_lut[ image_lut_nElem_per_bead * i + w*h*z ]; // current plane
681  float *imgL = &image_lut[ image_lut_nElem_per_bead * i + w*h*(z-1) ]; // one plane below
682  float *imgU = &image_lut[ image_lut_nElem_per_bead * i + w*h*(z+1) ]; // one plane above
683 
684  float *img_dz = &image_lut_dz[ image_lut_nElem_per_bead * i + w*h*z ];
685  float *img_dz2 = &image_lut_dz2[ image_lut_nElem_per_bead * i + w*h*z ];
686 
687  // Numerical approx of derivatives..
688  for (int y=0;y<h;y++) {
689  for (int x=0;x<w;x++) {
690  const float h = 1.0f;
691  img_dz[y*w+x] = 0.5f * ( imgU[y*w+x] - imgL[y*w+x] ) / h;
692  img_dz2[y*w+x] = (imgU[y*w+x] - 2*img[y*w+x] + imgL[y*w+x]) / (h*h);
693  }
694  }
695  }
696  for (int k=0;k<w*h;k++) {
697  // Top and bottom planes are simply copied from the neighbouring planes
698  image_lut_dz[ image_lut_nElem_per_bead * i + w*h*0 + k] = image_lut_dz[ image_lut_nElem_per_bead * i + w*h*1 + k];
699  image_lut_dz[ image_lut_nElem_per_bead * i + w*h*(image_lut_dims[1]-1) + k] = image_lut_dz[ image_lut_nElem_per_bead * i + w*h*(image_lut_dims[1]-2) + k];
701  image_lut_dz2[ image_lut_nElem_per_bead * i + w*h*(image_lut_dims[1]-1) + k] = image_lut_dz2[ image_lut_nElem_per_bead * i + w*h*(image_lut_dims[1]-2) + k];
702  }
703  }
704  }
705 }
float * image_lut
Image LUT data pointer.
int zlut_count
Number of ZLUTs (= number of beads).
int ImageLUTNumBeads()
Return the number of beads in the image LUT.
void NormalizeRadialProfile(scalar_t *prof, int rsteps)
Definition: utils.cpp:257
int zlut_radialsteps
Number of radial steps to sample on.
Definition: qtrk_c_api.h:198
float * image_lut_dz
Image LUT first derivative. (???)
float * zluts
Pointer to the first pixel of the Z lookup tables.
int ImageLUTWidth()
Return the height of the image LUT.
QTrkComputedConfig cfg
The settings used by this instance of QueuedTracker.
int zlut_planes
Number of planes per ZLUT.
int image_lut_dims[4]
Image LUT dimensions, 4D. [numBeads][numPlanes][height][width].
float * image_lut_dz2
Image LUT second derivative. (???)
int ImageLUTHeight()
Return the width of the image LUT.
int image_lut_nElem_per_bead
Image LUT number of pixels in LUT per bead. Is Width * Height * Planes.

§ Flush()

void QueuedCPUTracker::Flush ( )
inlineoverridevirtual

Stop waiting for more jobs to do, and just process the current batch.

Implements QueuedTracker.

Definition at line 94 of file QueuedCPUTracker.h.

94 { };

§ GenerateTestImage()

void QueuedCPUTracker::GenerateTestImage ( float *  dst,
float  xp,
float  yp,
float  size,
float  SNratio 
)

Generate a test image based on a the point spread function.

Creates an ImageData object and calls GenerateTestImage. Automatically uses ROI size from QtrkComputedConfig::width and QtrkComputedConfig::height.

Parameters
[out]dstPre-allocated float array in which to save the image.
[in]xpX coordinate of the PSF center.
[in]ypY coordinate of the PSF center.
[in]sizeSpread of the PSF. Can't be zero. Higher size increases radius of created diffraction pattern.
[in]SNratioSignal to noise ratio.

Definition at line 486 of file QueuedCPUTracker.cpp.

487 {
488  ImageData img(dst,cfg.width,cfg.height);
489  ::GenerateTestImage(img,xp,yp,size,SNratio);
490 }
int width
Width of regions of interest to be handled. Typically equals height (square ROI). ...
Definition: qtrk_c_api.h:106
int height
Height of regions of interest to be handled. Typically equals width (square ROI). ...
Definition: qtrk_c_api.h:107
QTrkComputedConfig cfg
The settings used by this instance of QueuedTracker.
void GenerateTestImage(float *dst, float xp, float yp, float size, float SNratio)
Generate a test image based on a the point spread function.

§ GetConfigValues()

QueuedTracker::ConfigValueMap QueuedCPUTracker::GetConfigValues ( )
overridevirtual

Get the used additional configurations.

Implements QueuedTracker.

Definition at line 510 of file QueuedCPUTracker.cpp.

511 {
512  ConfigValueMap cvm;
513  cvm["trace"] = dbgPrintResults ? "1" : "0";
514  return cvm;
515 }
std::map< std::string, std::string > ConfigValueMap
Datastructure used to return additional settings in a string-string key-value pair mapping...
bool dbgPrintResults
Flag to enable easy printing of intermediate data.

§ GetDebugImage()

bool QueuedCPUTracker::GetDebugImage ( int  id,
int *  w,
int *  h,
float **  data 
)
virtual

Get the debug image of a tracker.

Parameters
[in]idID of the tracker (= threadID).
[out]wReference to an int in which the image width will be stored.
[out]hReference to an int in which the image height will be stored.
[out]dataPointer to a float array that will be filled with the image.
Returns
True if the call was valid.

Reimplemented from QueuedTracker.

Definition at line 492 of file QueuedCPUTracker.cpp.

493 {
494  if (id >= 0 && id < threads.size()) {
495  threads[id].lock();
496 
497  *w = cfg.width;
498  *h = cfg.height;
499 
500  *data = new float [cfg.width*cfg.height];
501  memcpy(*data, threads[id].tracker->GetDebugImage(), sizeof(float)* cfg.width*cfg.height);
502 
503  threads[id].unlock();
504  return true;
505  }
506 
507  return false;
508 }
int width
Width of regions of interest to be handled. Typically equals height (square ROI). ...
Definition: qtrk_c_api.h:106
int height
Height of regions of interest to be handled. Typically equals width (square ROI). ...
Definition: qtrk_c_api.h:107
QTrkComputedConfig cfg
The settings used by this instance of QueuedTracker.
std::vector< Thread > threads
Vector with active threads.

§ GetImageLUTByIndex()

float* QueuedCPUTracker::GetImageLUTByIndex ( int  index,
int  plane = 0 
)
inlineprivate

Get a specific image LUT or plane from an image LUT.

Parameters
[in]indexThe beadindex of the requested LUT.
[in]planeIndex of the requested plane. 0 by default.
Returns
Pointer to the first element of the requested LUT or plane.

Definition at line 210 of file QueuedCPUTracker.h.

210  {
211  return &image_lut [ index * image_lut_nElem_per_bead + plane * (image_lut_dims[2]*image_lut_dims[3]) ];
212  }
float * image_lut
Image LUT data pointer.
int image_lut_dims[4]
Image LUT dimensions, 4D. [numBeads][numPlanes][height][width].
int image_lut_nElem_per_bead
Image LUT number of pixels in LUT per bead. Is Width * Height * Planes.

§ GetImageZLUT()

void QueuedCPUTracker::GetImageZLUT ( float *  dst)
virtual

Get the image lookup tables used.

Note
Use of image LUT is currently not clear. Radial ZLUT is always used.
Parameters
[out]dstPointer to the pre-allocated memory in which to save the data.

Reimplemented from QueuedTracker.

Definition at line 537 of file QueuedCPUTracker.cpp.

538 {
539  if (image_lut) {
540  memcpy(dst, image_lut, sizeof(float)*image_lut_nElem_per_bead*image_lut_dims[0]);
541  }
542 }
float * image_lut
Image LUT data pointer.
int image_lut_dims[4]
Image LUT dimensions, 4D. [numBeads][numPlanes][height][width].
int image_lut_nElem_per_bead
Image LUT number of pixels in LUT per bead. Is Width * Height * Planes.

§ GetImageZLUTSize()

void QueuedCPUTracker::GetImageZLUTSize ( int *  dims)
virtual

Get the dimensions of the image lookup table data.

Note
Use of image LUT is currently not clear. Radial ZLUT is always used.
Parameters
[out]dimsReference to pre-allocated int array. Returns [ count, planes, height, width ].

Reimplemented from QueuedTracker.

Definition at line 530 of file QueuedCPUTracker.cpp.

531 {
532  for (int i=0;i<4;i++)
533  dims[i]=image_lut_dims[i];
534 }
int image_lut_dims[4]
Image LUT dimensions, 4D. [numBeads][numPlanes][height][width].

§ GetNextJob()

QueuedCPUTracker::Job * QueuedCPUTracker::GetNextJob ( )
private

Get the next job in the queue.

Returns
The next job in the queue.

Definition at line 48 of file QueuedCPUTracker.cpp.

49 {
50  QueuedCPUTracker::Job *j = 0;
51  jobs_mutex.lock();
52  if (!jobs.empty()) {
53  j = jobs.front();
54  jobs.pop_front();
55  jobCount --;
56  jobsInProgress ++;
57  }
59  return j;
60 }
int jobsInProgress
Number of jobs currently being processed.
Threads::Mutex jobs_mutex
Mutex for the job queue.
void lock()
Definition: threads.h:74
int jobCount
Number of jobs in the queue.
std::deque< Job * > jobs
Queue to hold the jobs.
Structure around LocalizationJob to link data and job.
void unlock()
Definition: threads.h:79

§ GetProfileReport()

std::string QueuedCPUTracker::GetProfileReport ( )
inlinevirtual

Get the output of performance profiling.

Note
Profiling is only implemented in CUDA at the moment.
Returns
String with the parsed profiling output.

Reimplemented from QueuedTracker.

Definition at line 113 of file QueuedCPUTracker.h.

113 { return "CPU tracker currently has no profile reporting"; }

§ GetQueueLength()

int QueuedCPUTracker::GetQueueLength ( int *  maxQueueLen = 0)
overridevirtual

Get the lengths of the queue of jobs to be handled.

Parameters
[out]maxQueueLenPre-allocated integer that returns the maximum size of the queue if nonzero.
Returns
Number of jobs currently being handled and in the queue.

Implements QueuedTracker.

Definition at line 83 of file QueuedCPUTracker.cpp.

84 {
85  int jc;
86  jobs_mutex.lock();
87  jc = jobCount + jobsInProgress;
89 
90  if (maxQueueLength)
91  *maxQueueLength = this->maxQueueSize;
92 
93  return jc;
94 }
int jobsInProgress
Number of jobs currently being processed.
Threads::Mutex jobs_mutex
Mutex for the job queue.
void lock()
Definition: threads.h:74
int maxQueueSize
Maximum number of jobs in the queue.
int jobCount
Number of jobs in the queue.
void unlock()
Definition: threads.h:79

§ GetRadialZLUT()

void QueuedCPUTracker::GetRadialZLUT ( float *  dst)
overridevirtual

Get the radial lookup tables used for z tracking.

Parameters
[out]dstPointer to the pre-allocated memory in which to save the data.

Implements QueuedTracker.

Definition at line 440 of file QueuedCPUTracker.cpp.

441 {
443  if (nElem>0) {
445  memcpy(zlut, zluts, sizeof(float)* nElem);
447  }
448 }
int zlut_count
Number of ZLUTs (= number of beads).
void lock()
Definition: threads.h:74
Threads::Mutex results_mutex
Mutex for the results.
int zlut_radialsteps
Number of radial steps to sample on.
Definition: qtrk_c_api.h:198
float * zluts
Pointer to the first pixel of the Z lookup tables.
QTrkComputedConfig cfg
The settings used by this instance of QueuedTracker.
int zlut_planes
Number of planes per ZLUT.
void unlock()
Definition: threads.h:79

§ GetRadialZLUTCompareProfile()

void QueuedCPUTracker::GetRadialZLUTCompareProfile ( float *  dst)
virtual

Get saved error curve.

See EnableRadialZLUTCompareProfile.

Note
Not implemented for CUDA.
Parameters
[in]dstPointer to the pre-allocated memory in which to save the error curve. Size is count * planes.

Implements QueuedTracker.

Definition at line 550 of file QueuedCPUTracker.cpp.

551 {
552  if (zlut_cmpprofiles){
553  for (int i=0;i<zlut_count*zlut_planes;i++)
554  dst[i]=zlut_cmpprofiles[i];
555  }
556 }
float * zlut_cmpprofiles
Array in which to save errorcurves if enabled with EnableRadialZLUTCompareProfile.
int zlut_count
Number of ZLUTs (= number of beads).
int zlut_planes
Number of planes per ZLUT.

§ GetRadialZLUTSize()

void QueuedCPUTracker::GetRadialZLUTSize ( int &  count,
int &  planes,
int &  radialsteps 
)
overridevirtual

Get the dimensions of the radial lookup table data.

Parameters
[out]countReference to pre-allocated int. Returns number of lookup tables.
[out]planesReference to pre-allocated int. Returns number of planes per lookup table.
[out]radialstepsReference to pre-allocated int. Returns number of steps per plane.

Implements QueuedTracker.

Definition at line 433 of file QueuedCPUTracker.cpp.

434 {
435  count = zlut_count;
436  planes = zlut_planes;
437  rsteps = cfg.zlut_radialsteps;
438 }
int zlut_count
Number of ZLUTs (= number of beads).
int zlut_radialsteps
Number of radial steps to sample on.
Definition: qtrk_c_api.h:198
QTrkComputedConfig cfg
The settings used by this instance of QueuedTracker.
int zlut_planes
Number of planes per ZLUT.

§ GetResultCount()

int QueuedCPUTracker::GetResultCount ( )
overridevirtual

Get the number of finished localization jobs (=results) available in memory.

Returns
The number of available results.

Implements QueuedTracker.

Definition at line 21 of file QueuedCPUTracker.cpp.

22 {
24  int rc = resultCount;
26  return rc;
27 }
void lock()
Definition: threads.h:74
Threads::Mutex results_mutex
Mutex for the results.
int resultCount
Number of results available.
void unlock()
Definition: threads.h:79

§ GetZLUTByIndex()

float* QueuedCPUTracker::GetZLUTByIndex ( int  index)
inline

Get a pointer to the starting pixel of the specified LUT.

LUTs are saved in one large contiguous memory region. This function returns a reference to a specific LUT.

Parameters
[in]indexIndex of the requested LUT.
Returns
Pointer to the LUT starting pixel

Definition at line 123 of file QueuedCPUTracker.h.

123 { return &zluts[ index * (zlut_planes*cfg.zlut_radialsteps) ]; }
int zlut_radialsteps
Number of radial steps to sample on.
Definition: qtrk_c_api.h:198
float * zluts
Pointer to the first pixel of the Z lookup tables.
QTrkComputedConfig cfg
The settings used by this instance of QueuedTracker.
int zlut_planes
Number of planes per ZLUT.

§ ImageLUTHeight()

int QueuedCPUTracker::ImageLUTHeight ( )
inlineprivate

Return the width of the image LUT.

Definition at line 198 of file QueuedCPUTracker.h.

§ ImageLUTNumBeads()

int QueuedCPUTracker::ImageLUTNumBeads ( )
inlineprivate

Return the number of beads in the image LUT.

Definition at line 196 of file QueuedCPUTracker.h.

§ ImageLUTWidth()

int QueuedCPUTracker::ImageLUTWidth ( )
inlineprivate

Return the height of the image LUT.

Definition at line 197 of file QueuedCPUTracker.h.

§ IsIdle()

bool QueuedCPUTracker::IsIdle ( )
overridevirtual

Test to see if the tracker is idle.

That is, GetQueueLength == 0.

Returns
Boolean indicating if the tracker is idle.

Implements QueuedTracker.

Definition at line 16 of file QueuedCPUTracker.cpp.

17 {
18  return GetQueueLength() == 0;
19 }
int GetQueueLength(int *maxQueueLength=0) override
Get the lengths of the queue of jobs to be handled.

§ JobFinished()

void QueuedCPUTracker::JobFinished ( QueuedCPUTracker::Job j)
private

Flag job memory for reuse.

See jobs_buffer.

Warning
Only call when a job is fully finished, including results copying.
Parameters
[in]jThe job to free.

Definition at line 37 of file QueuedCPUTracker.cpp.

38 {
40  jobs_buffer.push_back(j);
42 
43  jobs_mutex.lock();
46 }
int jobsInProgress
Number of jobs currently being processed.
Threads::Mutex jobs_mutex
Mutex for the job queue.
Threads::Mutex jobs_buffer_mutex
Mutex for the job buffer.
std::vector< Job * > jobs_buffer
Stores memory. Enables reuse of allocated Job memory after a job has been completed.
void lock()
Definition: threads.h:74
void unlock()
Definition: threads.h:79

§ NumThreads()

int QueuedCPUTracker::NumThreads ( )
inline

Get the used number of threads.

Returns
The number of threads.

Definition at line 54 of file QueuedCPUTracker.h.

54 { return cfg.numThreads; }
QTrkComputedConfig cfg
The settings used by this instance of QueuedTracker.
int numThreads
Number of threads/streams to use. Defaults differ between CPU and GPU implementations.
Definition: qtrk_c_api.h:114

§ ProcessJob()

void QueuedCPUTracker::ProcessJob ( QueuedCPUTracker::Thread th,
Job j 
)
private

Process a job.

Parameters
[in]thThe thread on which to execute.
[in]jThe job to execute.

Choose the CPUTracker instance based on the specified thread.

Set the tracker's image memory.

Calibrate the images if needed.

Always compute the COM.

Compute 1DXCor, QI or 2DGaussian based on settings.

Compute Z profile and Z position if Z localization is requested.

Add the results to the available results.

Definition at line 280 of file QueuedCPUTracker.cpp.

281 {
283  CPUTracker* trk = th->tracker;
284  th->lock();
285 
287  SetTrackerImage(trk, j);
288 
289  //FloatToJPEGFile("dbg.jpg",(float*)j->data,cfg.width,cfg.height);
290 
292  trk->srcImage[0]=trk->srcImage[1]=trk->srcImage[2]=trk->srcImage[3]=0;
293  }
294 
296  ApplyOffsetGain(trk, j->job.zlutIndex);
297 
298 // dbgprintf("Job: id %d, bead %d\n", j->id, j->zlut);
299 
300  LocalizationResult result={};
301  result.job = j->job;
302 
305  result.imageMean = trk->mean;
306 
307  if (_isnan(com.x) || _isnan(com.y))
308  com = vector2f(cfg.width/2,cfg.height/2);
309 
310  bool boundaryHit = false;
311 
313  if (localizeMode & LT_XCor1D) {
314  result.firstGuess = com;
315  vector2f resultPos = trk->ComputeXCorInterpolated(com, cfg.xc1_iterations, cfg.xc1_profileWidth, boundaryHit);
316  result.pos.x = resultPos.x;
317  result.pos.y = resultPos.y;
318  } else if (localizeMode & LT_QI ){
319  result.firstGuess = com;
321  result.pos.x = resultPos.x;
322  result.pos.y = resultPos.y;
323  } else if (localizeMode & LT_Gaussian2D) {
324  result.firstGuess = com;
326  vector2f xy = gr.pos;
327  result.pos = vector3f(xy.x,xy.y,0.0f);
328  } else {
329  result.firstGuess.x = result.pos.x = com.x;
330  result.firstGuess.y = result.pos.y = com.y;
331  }
332 
334  bool normalizeProfile = (localizeMode & LT_NormalizeProfile)!=0;
335  if(localizeMode & LT_LocalizeZ) {
336  float* prof=ALLOCA_ARRAY(float,cfg.zlut_radialsteps);
337 
338  for (int i=0;i< ((localizeMode & LT_ZLUTAlign) ? 5 : 1) ; i++) {
339  if (localizeMode & LT_FourierLUT) {
341  } else {
342  trk->ComputeRadialProfile(prof,cfg.zlut_radialsteps, cfg.zlut_angularsteps, cfg.zlut_minradius, cfg.zlut_maxradius, result.pos2D(), false, &boundaryHit, normalizeProfile );
343  }
344 
345  float *cmpprof = 0;
346 
347  if (zlut_enablecmpprof)
348  cmpprof = &zlut_cmpprofiles[j->job.zlutIndex*zlut_planes];
349 
350  if (i > 0) {
351  // update with Quadrant Align
352  result.pos = trk->QuadrantAlign(result.pos, j->job.zlutIndex, cfg.qi_angstepspq, boundaryHit);
353  }
354  result.pos.z = trk->LUTProfileCompare(prof, j->job.zlutIndex, cmpprof, CPUTracker::LUTProfMaxQuadraticFit,(float*)0,(int*)0,j->job.frame);
355  //dbgprintf("[%d] x=%f, y=%f, z=%f\n", i, result.pos.x,result.pos.y,result.pos.z);
356 
357  }
358 
360  result.pos.z = trk->LUTProfileCompareAdjustedWeights(prof, j->job.zlutIndex, result.pos.z);
361  }
362 
364  result.pos.z = ZLUTBiasCorrection(result.pos.z, zlut_planes, j->job.zlutIndex);
365  }
366 
367  if(dbgPrintResults)
368  dbgprintf("fr:%d, bead: %d: x=%f, y=%f, z=%f\n",result.job.frame, result.job.zlutIndex, result.pos.x, result.pos.y, result.pos.z);
369 
370  th->unlock();
371 
372  result.error = boundaryHit ? 1 : 0;
375  results.push_back(result);
376  resultCount++;
378 }
int gauss2D_iterations
Number of times to run the 2D gaussian algorithm.
Definition: qtrk_c_api.h:163
float zlut_minradius
Distance in pixels from the bead center from which to start sampling profiles. Default 1...
Definition: qtrk_c_api.h:135
Ignore the first four pixels of a frame. Some cameras save timestamps rather than image data in the f...
Definition: qtrk_c_api.h:23
Gauss2DResult Compute2DGaussianMLE(vector2f initial, int iterations, float sigma)
Calculate a 2D Gaussian fit on the image.
void ApplyOffsetGain(CPUTracker *trk, int beadIndex)
Calibrate an image.
void ComputeRadialProfile(float *dst, int radialSteps, int angularSteps, float minradius, float maxradius, vector2f center, bool crp, bool *boundaryHit=0, bool normalize=true)
Wrapper to compute a radial profile.
std::vector< float > qi_radialbinweights
Scaling factors for the QI algorithm.
float qi_minradius
Distance in pixels from the bead center from which to start sampling profiles. Default 1...
Definition: qtrk_c_api.h:141
float * zlut_cmpprofiles
Array in which to save errorcurves if enabled with EnableRadialZLUTCompareProfile.
int width
Width of regions of interest to be handled. Typically equals height (square ROI). ...
Definition: qtrk_c_api.h:106
int zlutIndex
Bead number of this ROI. Used to get the right ZLUT from memory.
Definition: qtrk_c_api.h:58
Class with all CPU algorithm implementations.
Definition: cpu_tracker.h:38
void lock()
Definition: threads.h:74
Structure to group results from the 2D Gaussian fit.
Definition: cpu_tracker.h:180
Structure for job results.
Definition: qtrk_c_api.h:67
int height
Height of regions of interest to be handled. Typically equals width (square ROI). ...
Definition: qtrk_c_api.h:107
CPUTracker * tracker
The tracker for this thread.
vector2f pos
Found position.
Definition: cpu_tracker.h:181
LocMode_t localizeMode
Localization settings.
vector2< float > vector2f
Definition: std_incl.h:39
Threads::Mutex results_mutex
Mutex for the results.
float LUTProfileCompare(float *profile, int zlutIndex, float *cmpProf, LUTProfileMaxComputeMode maxPosMethod, float *fitcurve=0, int *maxPos=0, int frameNum=0)
Compare a profile to a LUT and calculate the optimum Z value for it.
int zlut_radialsteps
Number of radial steps to sample on.
Definition: qtrk_c_api.h:198
float imageMean
Average pixel value of the ROI.
Definition: qtrk_c_api.h:73
std::deque< LocalizationResult > results
Queue to store the localization results.
uint frame
Frame number this ROI belongs to.
Definition: qtrk_c_api.h:56
COM+XCor1D.
Definition: qtrk_c_api.h:8
vector2f ComputeMeanAndCOM(float bgcorrection=0.0f)
Calculate the center of mass of the image.
vector2f pos2D()
Final 2D position found.
Definition: qtrk_c_api.h:70
float * srcImage
Memory region that holds the image on which tracking is to be performed.
Definition: cpu_tracker.h:46
Make a fourier based lookup table.
Definition: qtrk_c_api.h:24
vector3f pos
Final 3D position found. If no z localization was performed, the value of z will be 0...
Definition: qtrk_c_api.h:69
Enable z localization.
Definition: qtrk_c_api.h:21
float com_bgcorrection
Background correction factor for COM. Defines the number of standard deviations data needs to be away...
Definition: qtrk_c_api.h:133
int resultCount
Number of results available.
float mean
Mean intensity of the ROI. Calculated in ComputeMeanAndCOM.
Definition: cpu_tracker.h:48
float qi_maxradius
Max radius in pixels of the sampling circle.
Definition: qtrk_c_api.h:204
Normalize found radial profiles.
Definition: qtrk_c_api.h:22
LocalizationJob job
Job metadata. See LocalizationJob.
Definition: qtrk_c_api.h:68
QTrkComputedConfig cfg
The settings used by this instance of QueuedTracker.
void lock()
Lock this thread&#39;s mutex.
void unlock()
Unlock this thread&#39;s mutex.
vector2f ComputeXCorInterpolated(vector2f initial, int iterations, int profileWidth, bool &boundaryHit)
Compute the cross correlation offsets and resulting position.
CImageData * zlut_bias_correction
int zlut_planes
Number of planes per ZLUT.
#define ALLOCA_ARRAY(T, N)
Definition: std_incl.h:153
float ZLUTBiasCorrection(float z, int zlut_planes, int bead)
Default. Use a 7-point quadratic fit around the error curve&#39;s maximum.
Definition: cpu_tracker.h:379
int qi_radialsteps
Number of radial steps to sample on.
Definition: qtrk_c_api.h:202
vector2f firstGuess
(x,y) position found by the COM localization. Used as initial position for the subsequent algorithms...
Definition: qtrk_c_api.h:71
COM+QI.
Definition: qtrk_c_api.h:9
2D Gaussian localization
Definition: qtrk_c_api.h:10
vector3< float > vector3f
Definition: std_incl.h:114
void dbgprintf(const char *fmt,...)
Definition: utils.cpp:149
int qi_angstepspq
Number of angular steps to sample on per quadrant.
Definition: qtrk_c_api.h:203
int qi_iterations
Number of times to run the QI algorithm, sampling around the last found position. ...
Definition: qtrk_c_api.h:140
float LUTProfileCompareAdjustedWeights(float *rprof, int zlutIndex, float z_estim)
Compare a profile to a LUT and calculate the optimum Z value for it using an initial estimate...
int xc1_iterations
Number of times to run the cross correlation algorithm.
Definition: qtrk_c_api.h:161
Enable ZLUT align.
Definition: qtrk_c_api.h:19
float qi_angstep_factor
Factor to reduce angular steps on lower iterations. Default 1.0 (no effect).
Definition: qtrk_c_api.h:157
int xc1_profileWidth
Profile width for the cross correlation.
Definition: qtrk_c_api.h:160
float gauss2D_sigma
Standard deviation to use in the 2D gaussian algorithm.
Definition: qtrk_c_api.h:164
bool zlut_enablecmpprof
Flag to save errorcurves. See EnableRadialZLUTCompareProfile.
uint error
Flag (boolean) indicating whether the ROI boundary was hit during localization. A 1 indicates a hit...
Definition: qtrk_c_api.h:72
void unlock()
Definition: threads.h:79
void FourierRadialProfile(float *dst, int radialSteps, int angularSteps, float minradius, float maxradius)
Calculate the radial profile based on the fourier transform.
vector3f QuadrantAlign(vector3f initial, int beadIndex, int angularStepsPerQuadrant, bool &boundaryHit)
Perform the quadrant align algorithm.
vector2f ComputeQI(vector2f initial, int iterations, int radialSteps, int angularStepsPerQuadrant, float angStepIterationFactor, float minRadius, float maxRadius, bool &boundaryHit, float *radialweights=0)
Execute the quadrant interpolation algorithm.
int zlut_angularsteps
Number of angular steps to sample on.
Definition: qtrk_c_api.h:199
float zlut_maxradius
Max radius in pixels of the sampling circle.
Definition: qtrk_c_api.h:200
void SetTrackerImage(CPUTracker *trk, Job *job)
Copy the ROI to the tracker&#39;s memory.
Do a ZLUT lookup with adjusted weights, for testing purposes.
Definition: qtrk_c_api.h:25
bool dbgPrintResults
Flag to enable easy printing of intermediate data.

§ ScheduleLocalization()

void QueuedCPUTracker::ScheduleLocalization ( void *  data,
int  pitch,
QTRK_PixelDataType  pdt,
const LocalizationJob jobInfo 
)
overridevirtual

Add a job to the queue to be processed. A job entails running the required algorithms on a single region of interest.

If a localization can not be added to the queue because it is full, the thread will be put to sleep and periodically try again.

Parameters
[in]dataPointer to the data. Type specified by pdt.
[in]pitchDistance in bytes between two successive rows of pixels (e.g. address of (0,0) - address of (0,1)).
[in]pdtType of data, specified by QTRK_PixelDataType.
[in]jobInfoStructure with metadata for the ROI to be handled. See LocalizationJob.
Bug:
So if process is stopped, it'll just queue regardless of queue sizes?

Implements QueuedTracker.

Definition at line 450 of file QueuedCPUTracker.cpp.

451 {
452  if (processJobs) {
453  while(maxQueueSize != 0 && GetQueueLength () >= maxQueueSize)
454  Threads::Sleep(5);
455  }
456 
457  Job* j = AllocateJob();
458  int dstPitch = PDT_BytesPerPixel(pdt) * cfg.width;
459 
460  if(!j->data || j->dataType != pdt) {
461  if (j->data) delete[] j->data;
462  j->data = new uchar[dstPitch * cfg.height];
463  }
464  for (int y=0; y<cfg.height; y++)
465  memcpy(&j->data[dstPitch*y], & ((uchar*)data)[pitch*y], dstPitch);
466 
467  j->dataType = pdt;
468  j->job = *jobInfo;
469 
470  AddJob(j);
471 }
int width
Width of regions of interest to be handled. Typically equals height (square ROI). ...
Definition: qtrk_c_api.h:106
Job * AllocateJob()
Allocate job memory.
int maxQueueSize
Maximum number of jobs in the queue.
bool processJobs
Flag for threads to continue processing jobs.
void AddJob(Job *j)
Add a job to the queue.
int height
Height of regions of interest to be handled. Typically equals width (square ROI). ...
Definition: qtrk_c_api.h:107
QTrkComputedConfig cfg
The settings used by this instance of QueuedTracker.
static void Sleep(int ms)
Definition: threads.h:134
int GetQueueLength(int *maxQueueLength=0) override
Get the lengths of the queue of jobs to be handled.
unsigned char uchar
Definition: std_incl.h:130
int PDT_BytesPerPixel(QTRK_PixelDataType pdt)
sizeof() equivalent for the QTRK_PixelDataType.

§ SetConfigValue()

void QueuedCPUTracker::SetConfigValue ( std::string  name,
std::string  value 
)
overridevirtual

Set an additional setting.

Parameters
[in]nameName of the setting.
[in]valueValue of the setting.

Implements QueuedTracker.

Definition at line 518 of file QueuedCPUTracker.cpp.

519 {
520  if (name == "trace")
521  dbgPrintResults = !!atoi(value.c_str());
522 }
bool dbgPrintResults
Flag to enable easy printing of intermediate data.

§ SetImageZLUT()

bool QueuedCPUTracker::SetImageZLUT ( float *  src,
float *  radial_zlut,
int *  dims 
)
virtual

Set the image lookup tables to be used for z tracking.

Note
Use of image LUT is currently not clear. Radial ZLUT is always used.
Parameters
[in]srcPointer to the data for the image LUT.
[in]radial_zlutPointer to the data for the radial LUT.
[in]dimsArray of dimension sizes for the image LUT. See GetImageZLUTSize.

Reimplemented from QueuedTracker.

Definition at line 558 of file QueuedCPUTracker.cpp.

559 {
560  if (image_lut) {
561  delete[] image_lut;
562  image_lut=0;
563  }
564 
565  for (int i=0;i<4;i++)
566  image_lut_dims[i]=dims[i];
567 
568  image_lut_nElem_per_bead = dims[1]*dims[2]*dims[3];
569 
570  if (image_lut_nElem_per_bead > 0 && dims[0] > 0) {
572  memset(image_lut, 0, sizeof(float)*image_lut_nElem_per_bead*image_lut_dims[0]);
573  }
574 
575  if (src) {
576  memcpy(image_lut, src, sizeof(float)*image_lut_nElem_per_bead*image_lut_dims[0]);
577  }
578 
579  SetRadialZLUT(radial_lut, dims[0], dims[1]);
580 
581  return true; // returning true indicates this implementation supports ImageLUT
582 }
float * image_lut
Image LUT data pointer.
void SetRadialZLUT(float *data, int num_zluts, int planes) override
Set the radial lookup tables to be used for z tracking.
int image_lut_dims[4]
Image LUT dimensions, 4D. [numBeads][numPlanes][height][width].
int image_lut_nElem_per_bead
Image LUT number of pixels in LUT per bead. Is Width * Height * Planes.

§ SetLocalizationMode()

void QueuedCPUTracker::SetLocalizationMode ( LocMode_t  locType)
overridevirtual

Select which algorithm is to be used.

Parameters
[in]locTypeAn integer used as a bitmask for settings based on LocalizeModeEnum.

Implements QueuedTracker.

Definition at line 524 of file QueuedCPUTracker.cpp.

525 {
526  while (!IsIdle());
527  localizeMode = lt;
528 }
bool IsIdle() override
Test to see if the tracker is idle.
LocMode_t localizeMode
Localization settings.

§ SetPixelCalibrationFactors()

void QueuedCPUTracker::SetPixelCalibrationFactors ( float  offsetFactor,
float  gainFactor 
)
overridevirtual

Set the pixel calibration factors.

The factors can be used to increase or decrease the effects of the images supplied through SetPixelCalibrationImages for further finetuning. These only have an effect when an image is actually set through that function.

Parameters
[in]offsetFactorFactor by which to scale the offset values.
[in]gainFactorFactor by which to scale the gain values.

Implements QueuedTracker.

Definition at line 205 of file QueuedCPUTracker.cpp.

206 {
207  gc_mutex.lock();
208  gc_gainFactor = gainFactor;
209  gc_offsetFactor = offsetFactor;
210  gc_mutex.unlock();
211 }
float gc_gainFactor
Factor by which to scale the gain.
void lock()
Definition: threads.h:74
Threads::Mutex gc_mutex
Image calibration mutex.
void unlock()
Definition: threads.h:79
float gc_offsetFactor
Factor by which to scale the pixel calibration offset.

§ SetPixelCalibrationImages()

void QueuedCPUTracker::SetPixelCalibrationImages ( float *  offset,
float *  gain 
)
overridevirtual

Set the pixel calibration images.

These images are used to scale the input image to get rid of background influences in the image. The values are per-pixel-per-ROI. Result = gain*(pixel+offset).

Parameters
[in]offsetArray with the offset values to use per pixel. Size and order is [width*height*numbeads].
[in]gainArray with gain values to use per pixel. Size and order is [width*height*numbeads].

Implements QueuedTracker.

Definition at line 172 of file QueuedCPUTracker.cpp.

173 {
174  if (zlut_count > 0) {
175  int nelem = cfg.width*cfg.height*zlut_count;
176 
177  if (calib_gain == 0 && gain) {
178  calib_gain = new float[nelem];
179  memcpy(calib_gain, gain, sizeof(float)*nelem);
180  }
181  else if (calib_gain && gain == 0) {
182  delete[] calib_gain;
183  calib_gain = 0;
184  }
185 
186  if (calib_offset == 0 && offset) {
187  calib_offset = new float[nelem];
188  memcpy(calib_offset, offset, sizeof(float)*nelem);
189  }
190  else if (calib_offset && offset == 0) {
191  delete[] calib_offset;
192  calib_offset = 0;
193  }
194 
195 #ifdef _DEBUG
196  std::string path = GetLocalModulePath();
197  for (int i=0;i<zlut_count;i++) {
198  if(calib_gain) FloatToJPEGFile( SPrintf("%s/gain-bead%d.jpg", path.c_str(), i).c_str(), &calib_gain[cfg.width*cfg.height*i], cfg.width,cfg.height);
199  if(calib_offset) FloatToJPEGFile( SPrintf("%s/offset-bead%d.jpg", path.c_str(), i).c_str(), &calib_offset[cfg.width*cfg.height*i], cfg.width,cfg.height);
200  }
201 #endif
202  }
203 }
int width
Width of regions of interest to be handled. Typically equals height (square ROI). ...
Definition: qtrk_c_api.h:106
int zlut_count
Number of ZLUTs (= number of beads).
int height
Height of regions of interest to be handled. Typically equals width (square ROI). ...
Definition: qtrk_c_api.h:107
float * calib_gain
Vector with gains used for pixel correction.
std::string GetLocalModulePath()
Definition: utils.cpp:114
QTrkComputedConfig cfg
The settings used by this instance of QueuedTracker.
void FloatToJPEGFile(const char *name, const float *d, int w, int h)
Definition: fastjpg.cpp:189
float * calib_offset
Vector with offsets used for pixel correction.
std::string SPrintf(const char *fmt,...)
Definition: utils.cpp:132

§ SetRadialWeights() [1/2]

void QueuedCPUTracker::SetRadialWeights ( float *  zcmp)
overridevirtual

Set radial weights used for comparing LUT profiles.

Parameters
[in]zcmpArray of radial weights to use. zcmp has to have zlut_radialsteps elements.

Implements QueuedTracker.

Definition at line 409 of file QueuedCPUTracker.cpp.

410 {
411  if (rweights)
412  zcmp.assign(rweights, rweights + cfg.zlut_radialsteps);
413  else
414  zcmp.clear();
415 
416  for (int i=0;i<threads.size();i++){
417  threads[i].lock();
418  threads[i].tracker->SetRadialWeights(rweights);
419  threads[i].unlock();
420  }
421 
422 }
std::vector< float > zcmp
Scaling factors for the ZLUT algorithm.
int zlut_radialsteps
Number of radial steps to sample on.
Definition: qtrk_c_api.h:198
QTrkComputedConfig cfg
The settings used by this instance of QueuedTracker.
std::vector< Thread > threads
Vector with active threads.

§ SetRadialWeights() [2/2]

void QueuedCPUTracker::SetRadialWeights ( std::vector< float >  weights)
inline

Overload of SetRadialWeights to use a vector.

Parameters
[in]weightsFloat vector with the radial weights.

Definition at line 67 of file QueuedCPUTracker.h.

67 { SetRadialWeights(&weights[0]); }
void SetRadialWeights(float *rweights) override
Set radial weights used for comparing LUT profiles.

§ SetRadialZLUT()

void QueuedCPUTracker::SetRadialZLUT ( float *  data,
int  count,
int  planes 
)
overridevirtual

Set the radial lookup tables to be used for z tracking.

Data can be zero to allocate ZLUT data. LUTs should have been created before by BuildLUT, but not necessarily by the current instance as long as setting match.

Parameters
[in]dataPointer to the start of the ZLUT data.
[in]countNumber of ZLUTs in the dataset.
[in]planesNumber of planes per ZLUT.

Implements QueuedTracker.

Definition at line 380 of file QueuedCPUTracker.cpp.

381 {
382 // jobs_mutex.lock();
383 // results_mutex.lock();
385  delete zlut_bias_correction;
386 
387  if (zluts) delete[] zluts;
388  int res = cfg.zlut_radialsteps;
389  int total = num_zluts*res*planes;
390  if (total > 0) {
391  zluts = new float[planes*res*num_zluts];
392  std::fill(zluts,zluts+(planes*res*num_zluts), 0.0f);
393  zlut_planes = planes;
394  zlut_count = num_zluts;
395  if(data)
396  std::copy(data, data+(planes*res*num_zluts), zluts);
397  }
398  else
399  zluts = 0;
400 
401  if (zlut_cmpprofiles) delete[] zlut_cmpprofiles;
402  zlut_cmpprofiles = new float [num_zluts*planes];
403 
404  UpdateZLUTs();
405 // results_mutex.unlock();
406 // jobs_mutex.unlock();
407 }
float * zlut_cmpprofiles
Array in which to save errorcurves if enabled with EnableRadialZLUTCompareProfile.
int zlut_count
Number of ZLUTs (= number of beads).
int zlut_radialsteps
Number of radial steps to sample on.
Definition: qtrk_c_api.h:198
float * zluts
Pointer to the first pixel of the Z lookup tables.
void UpdateZLUTs()
Update the ZLUTs for all threads.
QTrkComputedConfig cfg
The settings used by this instance of QueuedTracker.
CImageData * zlut_bias_correction
int zlut_planes
Number of planes per ZLUT.

§ SetTrackerImage()

void QueuedCPUTracker::SetTrackerImage ( CPUTracker trk,
Job job 
)
private

Copy the ROI to the tracker's memory.

Parameters
[in]trkThe tracker for which to set the image.
[in]jobThe job this tracker will execute. Image data is taken from the job.

Definition at line 736 of file QueuedCPUTracker.cpp.

737 {
738  if (j->dataType == QTrkU8) {
739  trk->SetImage8Bit(j->data, cfg.width);
740  } else if (j->dataType == QTrkU16) {
741  trk->SetImage16Bit((ushort*)j->data, cfg.width*2);
742  } else {
743  trk->SetImageFloat((float*)j->data);
744  }
745 }
int width
Width of regions of interest to be handled. Typically equals height (square ROI). ...
Definition: qtrk_c_api.h:106
16 bit unsigned int
Definition: qtrk_c_api.h:36
unsigned short ushort
Definition: std_incl.h:128
8 bit unsigned int
Definition: qtrk_c_api.h:35
void SetImage8Bit(uchar *srcImage, uint srcpitch)
Set an image with 8 bit type.
Definition: cpu_tracker.h:254
QTrkComputedConfig cfg
The settings used by this instance of QueuedTracker.
void SetImageFloat(float *srcImage)
Set an image with float type.
Definition: cpu_tracker.cpp:88
void SetImage16Bit(ushort *srcImage, uint srcpitch)
Set an image with 16 bit type.
Definition: cpu_tracker.h:247

§ Start()

void QueuedCPUTracker::Start ( )

Start up the threads and initiate their worker function.

Definition at line 220 of file QueuedCPUTracker.cpp.

221 {
222  quitWork = false;
223 
224  threads.resize(cfg.numThreads);
225  for (int k=0;k<cfg.numThreads;k++) {
226 
227  threads[k].mutex = new Threads::Mutex();
228 #if 0
229  threads[k].mutex->name = SPrintf("thread%d", k);
230  threads[k].mutex->trace = true;
231 #endif
232  threads[k].tracker = new CPUTracker(downsampleWidth, downsampleHeight, cfg.xc1_profileLength, false);//cfg.testRun);
233  threads[k].manager = this;
234  threads[k].tracker->trackerID = k;
235  }
236 
237  for (int k=0;k<threads.size();k++) {
239  Threads::SetBackgroundPriority(threads[k].thread, true);
240  }
241 
242  processJobs = true;
243 }
Class with all CPU algorithm implementations.
Definition: cpu_tracker.h:38
bool processJobs
Flag for threads to continue processing jobs.
static Handle * Create(ThreadEntryPoint method, void *param)
Definition: threads.h:99
QTrkComputedConfig cfg
The settings used by this instance of QueuedTracker.
bool quitWork
Signal threads to stop their work.
std::vector< Thread > threads
Vector with active threads.
static void WorkerThreadMain(void *arg)
The loop executed by the threads.
static void SetBackgroundPriority(Handle *thread, bool bg)
Definition: threads.h:118
int downsampleHeight
Height of ROIs after downsampling.
int numThreads
Number of threads/streams to use. Defaults differ between CPU and GPU implementations.
Definition: qtrk_c_api.h:114
int downsampleWidth
Width of ROIs after downsampling.
int xc1_profileLength
Profile length for the cross correlation.
Definition: qtrk_c_api.h:159
std::string SPrintf(const char *fmt,...)
Definition: utils.cpp:132

§ UpdateZLUTs()

void QueuedCPUTracker::UpdateZLUTs ( )
private

Update the ZLUTs for all threads.

Definition at line 424 of file QueuedCPUTracker.cpp.

425 {
426  for (int i=0;i<threads.size();i++){
427  threads[i].lock();
428  threads[i].tracker->SetRadialZLUT(zluts, zlut_planes, cfg.zlut_radialsteps, zlut_count, cfg.zlut_minradius, cfg.zlut_maxradius, false, false);
429  threads[i].unlock();
430  }
431 }
float zlut_minradius
Distance in pixels from the bead center from which to start sampling profiles. Default 1...
Definition: qtrk_c_api.h:135
int zlut_count
Number of ZLUTs (= number of beads).
int zlut_radialsteps
Number of radial steps to sample on.
Definition: qtrk_c_api.h:198
float * zluts
Pointer to the first pixel of the Z lookup tables.
QTrkComputedConfig cfg
The settings used by this instance of QueuedTracker.
std::vector< Thread > threads
Vector with active threads.
int zlut_planes
Number of planes per ZLUT.
float zlut_maxradius
Max radius in pixels of the sampling circle.
Definition: qtrk_c_api.h:200

§ WorkerThreadMain()

void QueuedCPUTracker::WorkerThreadMain ( void *  arg)
staticprivate

The loop executed by the threads.

Will start execution of a new job whenever available and allowed.

Parameters
[in]argReference to the QueuedCPUTracker::Thread that will execute this loop.

Definition at line 245 of file QueuedCPUTracker.cpp.

246 {
247  Thread* th = (Thread*)arg;
248  QueuedCPUTracker* this_ = th->manager;
249 
250  while (!this_->quitWork) {
251  Job* j = 0;
252  if (this_->processJobs)
253  j = this_->GetNextJob();
254 
255  if (j) {
256  this_->ProcessJob(th, j);
257  this_->JobFinished(j);
258  } else {
259  Threads::Sleep(1);
260  }
261  }
262  //dbgprintf("Thread %p ending.\n", arg);
263 }
bool processJobs
Flag for threads to continue processing jobs.
void JobFinished(Job *j)
Flag job memory for reuse.
CPU implementation of the QueuedTracker interface.
bool quitWork
Signal threads to stop their work.
Job * GetNextJob()
Get the next job in the queue.
static void Sleep(int ms)
Definition: threads.h:134
void ProcessJob(Thread *th, Job *j)
Process a job.

§ ZLUTSelfTest()

void QueuedCPUTracker::ZLUTSelfTest ( )

Test the planes saved in the lookup table against itself.

Implemented for testing purposes. Writes the output to a CSV file. Can be used to inspect ZLUT quality.

Definition at line 707 of file QueuedCPUTracker.cpp.

708 {
709  if (zluts){
710  /*bool compEnabled = zlut_enablecmpprof;
711  if(!compEnabled)
712  EnableRadialZLUTCompareProfile(true);*/
713 
714 
715  threads[0].lock();
716  CPUTracker* trk = threads[0].tracker;
717 
718  for (int ii=0;ii<zlut_count;ii++) {
719  //trk.SetRadialZLUT(&zluts[cfg.zlut_radialsteps*zlut_planes*ii],zlut_planes,cfg.zlut_radialsteps,zlut_count,cfg.zlut_minradius,cfg.zlut_maxradius,false,false);
720 
721  float* curZLUT = GetZLUTByIndex(ii);
722  for(int plane=0;plane<zlut_planes;plane++){
723  float* cmpprof = new float[zlut_planes];
725  WriteArrayAsCSVRow("D:\\TestImages\\zlutSelfTest.csv",cmpprof,zlut_planes,true);
726  delete[] cmpprof;
727  }
728  }
729 
730  threads[0].unlock();
731 
732  //EnableRadialZLUTCompareProfile(compEnabled);
733  }
734 }
void WriteArrayAsCSVRow(const char *file, float *d, int len, bool append)
Definition: utils.cpp:537
int zlut_count
Number of ZLUTs (= number of beads).
Class with all CPU algorithm implementations.
Definition: cpu_tracker.h:38
float LUTProfileCompare(float *profile, int zlutIndex, float *cmpProf, LUTProfileMaxComputeMode maxPosMethod, float *fitcurve=0, int *maxPos=0, int frameNum=0)
Compare a profile to a LUT and calculate the optimum Z value for it.
int zlut_radialsteps
Number of radial steps to sample on.
Definition: qtrk_c_api.h:198
float * zluts
Pointer to the first pixel of the Z lookup tables.
float * GetZLUTByIndex(int index)
Get a pointer to the starting pixel of the specified LUT.
QTrkComputedConfig cfg
The settings used by this instance of QueuedTracker.
std::vector< Thread > threads
Vector with active threads.
int zlut_planes
Number of planes per ZLUT.
Default. Use a 7-point quadratic fit around the error curve&#39;s maximum.
Definition: cpu_tracker.h:379

Member Data Documentation

§ calib_gain

float* QueuedCPUTracker::calib_gain
private

Vector with gains used for pixel correction.

Definition at line 166 of file QueuedCPUTracker.h.

§ calib_offset

float* QueuedCPUTracker::calib_offset
private

Vector with offsets used for pixel correction.

Definition at line 167 of file QueuedCPUTracker.h.

§ dbgPrintResults

bool QueuedCPUTracker::dbgPrintResults
private

Flag to enable easy printing of intermediate data.

Definition at line 216 of file QueuedCPUTracker.h.

§ downsampleHeight

int QueuedCPUTracker::downsampleHeight
private

Height of ROIs after downsampling.

Definition at line 172 of file QueuedCPUTracker.h.

§ downsampleWidth

int QueuedCPUTracker::downsampleWidth
private

Width of ROIs after downsampling.

Definition at line 171 of file QueuedCPUTracker.h.

§ gc_gainFactor

float QueuedCPUTracker::gc_gainFactor
private

Factor by which to scale the gain.

Definition at line 168 of file QueuedCPUTracker.h.

§ gc_mutex

Threads::Mutex QueuedCPUTracker::gc_mutex
private

Image calibration mutex.

Definition at line 165 of file QueuedCPUTracker.h.

§ gc_offsetFactor

float QueuedCPUTracker::gc_offsetFactor
private

Factor by which to scale the pixel calibration offset.

Definition at line 169 of file QueuedCPUTracker.h.

§ image_lut

float* QueuedCPUTracker::image_lut
private

Image LUT data pointer.

Definition at line 199 of file QueuedCPUTracker.h.

§ image_lut_dims

int QueuedCPUTracker::image_lut_dims[4]
private

Image LUT dimensions, 4D. [numBeads][numPlanes][height][width].

Definition at line 194 of file QueuedCPUTracker.h.

§ image_lut_dz

float* QueuedCPUTracker::image_lut_dz
private

Image LUT first derivative. (???)

Definition at line 200 of file QueuedCPUTracker.h.

§ image_lut_dz2

float* QueuedCPUTracker::image_lut_dz2
private

Image LUT second derivative. (???)

Definition at line 201 of file QueuedCPUTracker.h.

§ image_lut_nElem_per_bead

int QueuedCPUTracker::image_lut_nElem_per_bead
private

Image LUT number of pixels in LUT per bead. Is Width * Height * Planes.

Definition at line 195 of file QueuedCPUTracker.h.

§ jobCount

int QueuedCPUTracker::jobCount
private

Number of jobs in the queue.

Note
Why maintain this manually and not just use jobs.size()?

Definition at line 158 of file QueuedCPUTracker.h.

§ jobs

std::deque<Job*> QueuedCPUTracker::jobs
private

Queue to hold the jobs.

Definition at line 152 of file QueuedCPUTracker.h.

§ jobs_buffer

std::vector<Job*> QueuedCPUTracker::jobs_buffer
private

Stores memory. Enables reuse of allocated Job memory after a job has been completed.

Definition at line 159 of file QueuedCPUTracker.h.

§ jobs_buffer_mutex

Threads::Mutex QueuedCPUTracker::jobs_buffer_mutex
private

Mutex for the job buffer.

Definition at line 150 of file QueuedCPUTracker.h.

§ jobs_mutex

Threads::Mutex QueuedCPUTracker::jobs_mutex
private

Mutex for the job queue.

Definition at line 149 of file QueuedCPUTracker.h.

§ jobsInProgress

int QueuedCPUTracker::jobsInProgress
private

Number of jobs currently being processed.

Definition at line 163 of file QueuedCPUTracker.h.

§ localizeMode

LocMode_t QueuedCPUTracker::localizeMode
private

Localization settings.

Definition at line 148 of file QueuedCPUTracker.h.

§ maxQueueSize

int QueuedCPUTracker::maxQueueSize
private

Maximum number of jobs in the queue.

Definition at line 162 of file QueuedCPUTracker.h.

§ processJobs

bool QueuedCPUTracker::processJobs
private

Flag for threads to continue processing jobs.

Definition at line 215 of file QueuedCPUTracker.h.

§ qi_radialbinweights

std::vector<float> QueuedCPUTracker::qi_radialbinweights
private

Scaling factors for the QI algorithm.

Definition at line 190 of file QueuedCPUTracker.h.

§ quitWork

bool QueuedCPUTracker::quitWork
private

Signal threads to stop their work.

Definition at line 214 of file QueuedCPUTracker.h.

§ resultCount

int QueuedCPUTracker::resultCount
private

Number of results available.

Definition at line 161 of file QueuedCPUTracker.h.

§ results

std::deque<LocalizationResult> QueuedCPUTracker::results
private

Queue to store the localization results.

Definition at line 160 of file QueuedCPUTracker.h.

§ results_mutex

Threads::Mutex QueuedCPUTracker::results_mutex
private

Mutex for the results.

Definition at line 151 of file QueuedCPUTracker.h.

§ threads

std::vector<Thread> QueuedCPUTracker::threads
private

Vector with active threads.

Definition at line 174 of file QueuedCPUTracker.h.

§ zcmp

std::vector<float> QueuedCPUTracker::zcmp
private

Scaling factors for the ZLUT algorithm.

Definition at line 189 of file QueuedCPUTracker.h.

§ zlut_buildflags

uint QueuedCPUTracker::zlut_buildflags
private

ZLUT build flags, set through BeginLUT.

Definition at line 187 of file QueuedCPUTracker.h.

§ zlut_cmpprofiles

float* QueuedCPUTracker::zlut_cmpprofiles
private

Array in which to save errorcurves if enabled with EnableRadialZLUTCompareProfile.

Definition at line 183 of file QueuedCPUTracker.h.

§ zlut_count

int QueuedCPUTracker::zlut_count
private

Number of ZLUTs (= number of beads).

Definition at line 185 of file QueuedCPUTracker.h.

§ zlut_enablecmpprof

bool QueuedCPUTracker::zlut_enablecmpprof
private

Flag to save errorcurves. See EnableRadialZLUTCompareProfile.

Definition at line 184 of file QueuedCPUTracker.h.

§ zlut_planes

int QueuedCPUTracker::zlut_planes
private

Number of planes per ZLUT.

Definition at line 186 of file QueuedCPUTracker.h.

§ zluts

float* QueuedCPUTracker::zluts
private

Pointer to the first pixel of the Z lookup tables.

All LUTs are saved in one big contiguous section of memory. Calculate specific LUTs or planes based on their indexes. Order is [beadIndex][plane][step]. See GetZLUTByIndex.

Definition at line 182 of file QueuedCPUTracker.h.


The documentation for this class was generated from the following files: