Updated Changelog
[libgpiv.git] / lib / valid_proc.c
bloba57b3cb06d39614970f507c58e7ae5f6d533c95e
1 /* -*- Mode: C; indent-tabs-mode: nil; c-basic-offset: 4 c-style: "K&R" -*- */
3 /*
4 libgpiv - library for Particle Image Velocimetry
6 Copyright (C) 2002, 2003, 2004, 2005, 2006, 2007, 2008
7 Gerber van der Graaf
9 This file is part of libgpiv.
10 Libgpiv is free software; you can redistribute it and/or modify
11 it under the terms of the GNU General Public License as published by
12 the Free Software Foundation; either version 2, or (at your option)
13 any later version.
15 This program is distributed in the hope that it will be useful,
16 but WITHOUT ANY WARRANTY; without even the implied warranty of
17 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
18 GNU General Public License for more details.
20 You should have received a copy of the GNU General Public License
21 along with this program; if not, write to the Free Software Foundation,
22 Inc., 59 Temple Place - Suite 330, Boston, MA 02111-1307, USA.
26 -----------------------------------------------------------------------
27 FILENAME: valid_proc.c
28 LIBRARY: libgpiv
29 EXTERNAL FUNCTIONS: gpiv_valid_residu
30 gpiv_valid_peaklck
31 gpiv_valid_residu_stats
32 gpiv_valid_errvec
33 gpiv_valid_gradient
34 gpiv_valid_threshold
35 gpiv_cumhisto_eqdatbin_gnuplot
37 LAST MODIFICATION DATE: $Id: valid.c,v 1.25 2008-09-25 13:19:53 gerber Exp $
38 --------------------------------------------------------------------- */
39 #include <gpiv.h>
41 #define SNR_ERR 99.0
42 #define MIN_VECTORS4MEDIAN 3 /* Minimum number of valid vectors needed
43 to calculate median */
44 #define RESIDU_EPSI 0.1
47 * Local functions
50 static int
51 compare_float (const void * a,
52 const void * b)
53 /*-----------------------------------------------------------------------------
54 * DESCRIPTION:
55 * Compares two float numbers
57 * PROTOTYPE LOCATATION:
58 * valid.h
60 * INPUTS:
61 * a: first float number
62 * b: second float number
64 * OUTPUTS:
66 * RETURNS:
67 * int: -1, 0 or +1 for a < b, a = b or a > b respectively
68 *---------------------------------------------------------------------------*/
70 float *la = (float *) a, *lb = (float *) b;
72 if (*la > *lb)
73 return 1;
74 else if (*la < *lb)
75 return -1;
76 else
77 return 0;
81 static float
82 median_residu (const guint i,
83 const guint j,
84 const guint neighbors,
85 const gboolean incl_point,
86 const gboolean norm,
87 const GpivPivData *data,
88 guint *i_median_x,
89 guint *j_median_x,
90 guint *i_median_y,
91 guint *j_median_y
93 /*-----------------------------------------------------------------------------
94 * DESCRIPTION:
95 * If incl_point. calculates residu value of point i, j from a NxN data
96 * array out of data[i][j] and the indices of the median velocity
97 * If not incl_point, returns the median residu from NxN data array,
98 * calculated from dx(i) - dx(median), i,= 1, ..,8
100 * PROTOTYPE LOCATATION:
101 * valid.h
103 * INPUTS:
104 * i: horizontal index of data to be investigated
105 * j: vertical index of data to be investigated
106 * neighbors: number of neighboring PIV data (typically 3x3, 5x5)
107 * data: piv dataset
108 * incl_point: flag to include data point under investigation
109 * norm: flag to return residu for normalizing
111 * OUTPUTS:
112 * i_median_x: i-index of median for x displacement
113 * j_median_x: j-index of median for x displacement
114 * i_median_y: i-index of median for y displacement
115 * j_median_y: j-index of median for y displacement
117 * RETURNS:
118 * median residu
119 *---------------------------------------------------------------------------*/
121 int k, l, m = 0;
122 gint vlength = 0, median_index;
123 const gint N = neighbors - 2;
124 float *dx_m, *dy_m, *dx_org, *dy_org;
125 float *r_m;
126 float residu;
129 * Obtain length of array
131 #pragma parallel for shared(data, incl_point, vlength) private(k, l)
132 for (k = i - N; k <= i + N; k++) {
133 if ((k >= 0) && (k < data->ny)) {
134 for (l = j - N; l <= j + N; l++) {
135 if ((l >= 0) && (l < data->nx)) {
136 if ((incl_point
137 && data->peak_no[k][l] != -1)
138 || (!incl_point
139 && (k != i || l != j)
140 && data->peak_no[k][l] != -1) ) {
141 vlength++;
148 if (data->peak_no[i][j] != -1 && vlength >= MIN_VECTORS4MEDIAN ) {
149 dx_m = gpiv_vector(vlength);
150 dy_m = gpiv_vector(vlength);
151 dx_org = gpiv_vector(vlength);
152 dy_org = gpiv_vector(vlength);
153 if (norm) {
154 r_m = gpiv_vector(vlength);
157 * Write the absolute neighbouring velocities and its components to a
158 * 1-dimensional array
160 m = 0;
161 for (k = i - N; k <= i + N; k++) {
162 if ((k >= 0) && (k < data->ny)) {
163 for (l = j - N; l <= j + N; l++) {
164 if ((l >= 0) && (l < data->nx)) {
165 if ((incl_point
166 && data->peak_no[k][l] != -1)
167 || (!incl_point
168 && (k != i || l != j)
169 && data->peak_no[k][l] != -1) ) {
170 dx_m[m] = data->dx[k][l];
171 dx_org[m] = data->dx[k][l];
172 dy_m[m] = data->dy[k][l];
173 dy_org[m] = data->dy[k][l];
174 m++;
182 * Sorting dx_m and dy_m arrays and searching median index and value
184 /* #pragma omp sections */
185 /* { */
186 /* #pragma omp section */
187 /* { */
188 qsort(dx_m, vlength, sizeof(float), compare_float);
189 /* } */
190 /* #pragma omp section */
191 /* { */
192 qsort(dy_m, vlength, sizeof(float), compare_float);
193 /* } */
194 /* } */
196 if (incl_point) {
197 median_index = (int) (vlength - 1) / 2;
198 } else {
199 median_index = (int) (vlength) / 2;
202 if (norm) {
204 * Obtain all residus from surrounding data, sorting and picking the median
206 for (k = 0; k < vlength; k++) {
207 r_m[k] = sqrt((dx_m[k] - dx_m[median_index]) *
208 (dx_m[k] - dx_m[median_index]) +
209 (dy_m[k] - dy_m[median_index]) *
210 (dy_m[k] - dy_m[median_index]));
212 qsort(r_m, vlength, sizeof(float), compare_float);
213 residu = r_m[median_index];
214 } else {
216 * Obtain residu from difference between current displacement at (i,j)
217 * and median displacement from the surroundings.
219 residu = sqrt((data->dx[i][j] - dx_m[median_index]) *
220 (data->dx[i][j] - dx_m[median_index]) +
221 (data->dy[i][j] - dy_m[median_index]) *
222 (data->dy[i][j] - dy_m[median_index]));
226 * Search the indexes for the 2-dim arrays dx and dy belonging
227 * to mediaan_index
229 m = 0;
230 for (k = i - N; k <= i + N; k++) {
231 if ((k >= 0) && (k < data->ny)) {
232 for (l = j - N; l <= j + N; l++) {
233 if ((l >= 0) && (l < data->nx)) {
234 if ((incl_point
235 && data->peak_no[k][l] != -1)
236 || (!incl_point
237 && (k != i || l != j)
238 && data->peak_no[k][l] != -1) ) {
239 if (dx_org[m] == dx_m[median_index]) {
240 *i_median_x = k - i;
241 *j_median_x = l - j;
243 if (dy_org[m] == dy_m[median_index]) {
244 *i_median_y = k - i;
245 *j_median_y = l - j;
247 m++;
254 gpiv_free_vector (dx_m);
255 gpiv_free_vector (dx_org);
256 gpiv_free_vector (dy_m);
257 gpiv_free_vector (dy_org);
258 if (norm) {
259 gpiv_free_vector (r_m);
262 } else {
263 residu = SNR_ERR;
266 return residu;
271 gchar *
272 gpiv_valid_residu (GpivPivData *piv_data,
273 const GpivValidPar *valid_par,
274 const gboolean incl_point
276 /*-----------------------------------------------------------------------------
277 * DESCRIPTION:
278 * Calculates residu values (at the inner points) of a PIV data set and
279 * applies to snr member of out_data
281 * PROTOTYPE LOCATATION:
282 * valid.h
284 * INPUTS:
285 * piv_data: piv dataset
286 * valid_par: validation parameters
287 * incl_point: flag to include data point under investigation
289 * OUTPUTS:
290 * out_data: piv dataset containing residu values in snr
292 * RETURNS:
293 *---------------------------------------------------------------------------*/
295 gchar *err_msg = NULL;
296 gint i, j;
299 /* out_data = gpiv_cp_pivdata (piv_data); */
301 if (valid_par->residu_type == GPIV_VALID_RESIDUTYPE__SNR) {
302 /* Nothing to do here */
304 } else if (valid_par->residu_type == GPIV_VALID_RESIDUTYPE__MEDIAN) {
306 for (i = 1; i < piv_data->ny - 1; i++) {
307 for (j = 1; j < piv_data->nx - 1; j++) {
308 guint i_mx, j_mx, i_my, j_my;
310 piv_data->snr[i][j] =
311 median_residu (i, j, valid_par->neighbors,
312 incl_point, FALSE, piv_data,
313 &i_mx, &j_mx, &i_my, &j_my);
314 if (piv_data->snr[i][j] == SNR_ERR)
315 piv_data->peak_no[i][j] = -1;
319 } else if (valid_par->residu_type == GPIV_VALID_RESIDUTYPE__NORMMEDIAN) {
320 gfloat residu_from_median_dxdy = 1.111, residu_norm = 5.555;
322 for (i = 0; i < piv_data->ny; i++) {
323 for (j = 0; j < piv_data->nx; j++) {
324 guint i_mx, j_mx, i_my, j_my;
326 residu_norm =
327 median_residu (i, j, valid_par->neighbors,
328 FALSE, TRUE, piv_data,
329 &i_mx, &j_mx, &i_my, &j_my);
330 residu_from_median_dxdy =
331 median_residu (i, j, valid_par->neighbors,
332 FALSE, FALSE, piv_data,
333 &i_mx, &j_mx, &i_my, &j_my);
334 if (residu_from_median_dxdy == SNR_ERR) {
335 piv_data->snr[i][j] = SNR_ERR;
336 piv_data->peak_no[i][j] = -1;
337 } else if (residu_norm + RESIDU_EPSI != 0.0) {
338 piv_data->snr[i][j] =
339 residu_from_median_dxdy /
340 (residu_norm + RESIDU_EPSI);
341 } else {
342 piv_data->snr[i][j] = SNR_ERR;
343 piv_data->peak_no[i][j] = -1;
349 } else {
350 gpiv_warning("gpiv valid residu: should no arrive here");
354 return err_msg;
359 static void
360 interr_reg_nhpeak(const guint index_y,
361 const guint index_x,
362 const GpivImage *image,
363 const GpivPivPar *piv_par,
364 GpivPivData *piv_data,
365 GpivFt *ft,
366 int tid
368 /*-----------------------------------------------------------------------------
369 * DESCRIPTION:
370 * Interrogates the image(pair) at a single region at the next higher
371 * correlation peak from a previous analysis
373 * PROTOTYPE LOCATATION:
374 * valid.h
376 * INPUTS:
377 * index_y: vertical array index (row)
378 * index_x: horizontal array index (column)
379 * image_par: image parameters
380 * image: image
381 * piv_par: PIV evaluation parameters
383 * OUTPUTS:
384 * piv_data: piv dataset containing particle image displacements
386 * RETURNS:
387 *---------------------------------------------------------------------------*/
389 char c_line[GPIV_MAX_LINES_C][GPIV_MAX_CHARS];
390 int nc_lines = 0;
392 int int_size_0;
394 GpivPivPar *lo_piv_par = NULL;
395 int sweep = 1, sweep_last = 1;
396 int return_val;
397 int cmpr = 1, cmpr_fact = 1;
400 * Checking for memory allocation of input variables
402 g_assert (image->frame1[0] != NULL);
403 g_assert (image->frame2[0] != NULL);
406 * Local (actualized) parameters
408 lo_piv_par = gpiv_piv_cp_parameters (piv_par);
409 lo_piv_par->peak = piv_data->peak_no[index_y][index_x] + 1;
412 * Memory allocation of interrogation area's and packed interrogation area
413 * arrays. Define weight kernel values
415 int_size_0 = ft->size;
419 * Interrogate at a single point
421 gpiv_piv_interrogate_ia (index_y, index_x, image, lo_piv_par, sweep,
422 sweep_last, piv_data, ft);
427 static gboolean
428 subst_vector (GpivPivData *piv_data,
429 const GpivImage *image,
430 const GpivPivPar *piv_par,
431 const GpivValidPar *valid_par,
432 GpivFt *ft
434 /*-----------------------------------------------------------------------------
435 * DESCRIPTION:
436 * Substitutes data with residu values higher than threshold
438 * PROTOTYPE LOCATATION:
439 * valid.h
441 * INPUTS:
442 * piv_data: piv dataset containing particle image displacements
443 * image: image struct containing parameters and image frames
444 * piv_par: PIV evaluation parameters
445 * valid_par: validation parameters
447 * OUTPUTS:
448 * l_data: piv dataset containing particle image displacements
450 * RETURNS: flag set TRUE if any data of the input set has been substituted
451 *---------------------------------------------------------------------------*/
453 gboolean *outlier_found = NULL;
454 int max_nr_thr = -1;
455 int nr_thr = -1;
456 int tid = -1;
458 gint i, j;
459 guint k, l;
460 const gint N = valid_par->neighbors - 2;
461 char line[GPIV_MAX_CHARS], command[GPIV_MAX_CHARS];
462 guint peak_no;
463 gfloat row = 0, col = 0;
464 gfloat dx, dy, snr;
465 FILE *fp;
467 guint count = 0;
468 gfloat dx_sum = 0, dy_sum = 0;
470 #ifdef ENABLE_OMP
471 max_nr_thr = omp_get_max_threads();
472 #else
473 max_nr_thr = 1;
474 #endif
475 outlier_found = malloc(max_nr_thr * sizeof(gboolean));
476 for (i = 0; i < max_nr_thr; i++) {
477 outlier_found[i] = FALSE;
482 * Create a local PIV data set, so input data will not be modified halfway
483 * during substitution. Else, this might affect PIV displacements
484 * calculated from (eventually previously substituted) surrounding values.
486 GpivPivData *l_data = NULL;
489 l_data = gpiv_cp_pivdata (piv_data);
491 if (valid_par->subst_type ==
492 GPIV_VALID_SUBSTYPE__NONE) {
494 * no substitution, only sets peak_no to 0
497 #pragma omp parallel default(none) \
498 private (nr_thr, tid, j) \
499 shared(i, piv_data, valid_par, l_data, outlier_found)
501 #ifdef ENABLE_OMP
502 tid = omp_get_thread_num();
503 #else
504 tid=0; /* init -1, for single thread set to 0 to address 1st array element */
505 #endif
506 #pragma omp for /* schedule(static) */
507 for (i = 0; i < piv_data->ny; i++) {
508 for (j = 0; j < piv_data->nx; j++) {
510 if (piv_data->peak_no[i][j] != -1
511 && piv_data->snr[i][j] > valid_par->residu_max) {
512 outlier_found[tid] = TRUE;
514 l_data->dx[i][j] = piv_data->dx[i][j];
515 l_data->dy[i][j] = piv_data->dy[i][j];
516 l_data->snr[i][j] = piv_data->snr[i][j];
517 if (l_data->snr[i][j] == SNR_ERR) {
518 l_data->peak_no[i][j] = -1;
519 } else {
520 l_data->peak_no[i][j] = 0;
525 } /* end of parallel region */
528 } else if (valid_par->subst_type ==
529 GPIV_VALID_SUBSTYPE__L_MEAN) {
531 #pragma omp parallel default(none) \
532 private (nr_thr, tid, j, k, l, count, dx_sum, dy_sum) \
533 shared (i, piv_data, valid_par, l_data, outlier_found)
535 #ifdef ENABLE_OMP
536 tid = omp_get_thread_num();
537 #else
538 tid=0; /* init -1, for single thread set to 0 to address 1st array element */
539 #endif
540 #pragma omp for /* schedule(static) */
541 for (i = 0; i < piv_data->ny; i++) {
542 for (j = 0; j < piv_data->nx; j++) {
543 count = 0;
544 dx_sum = 0.0;
545 dy_sum = 0.0;
547 if (piv_data->peak_no[i][j] != -1
548 && piv_data->snr[i][j] > valid_par->residu_max) {
549 outlier_found[tid] = TRUE;
551 for (k = i - N; k <= i + N; k++) {
552 if ((k >= 0) && (k < piv_data->ny)) {
553 for (l = j - N; l <= j + N; l++) {
554 if ((l >= 0) && (l < piv_data->nx)) {
556 * Exclude the point under investigation for calculating the mean
558 if ((k != 0) || (l != 0)) {
559 dx_sum += piv_data->dx[k][l];
560 dy_sum += piv_data->dy[k][l];
561 count++;
567 l_data->dx[i][j] = dx_sum / count;
568 l_data->dy[i][j] = dy_sum / count;
569 l_data->snr[i][j] = piv_data->snr[i][j];
570 l_data->peak_no[i][j] = 0;
574 } /* end of parallel region */
576 } else if (valid_par->subst_type ==
577 GPIV_VALID_SUBSTYPE__MEDIAN) {
579 * Substitution with median particle displacements
581 #pragma omp parallel default(none) /* "for" NOT missing! */ \
582 private (tid,j) \
583 shared (i,piv_data,valid_par,outlier_found,l_data)
585 #ifdef ENABLE_OMP
586 tid = omp_get_thread_num();
587 #else
588 tid=0; /* init -1, for single thread set to 0 to address 1st array element */
589 #endif
590 guint i_mx, j_mx, i_my, j_my;
591 #pragma omp for /* schedule(static) */
592 for (i = 0; i < piv_data->ny; i++) {
593 for (j = 0; j < piv_data->nx; j++) {
595 if (piv_data->peak_no[i][j] != -1
596 && piv_data->snr[i][j] > valid_par->residu_max) {
597 outlier_found[tid] = TRUE;
599 l_data->snr[i][j] =
600 median_residu (i, j, valid_par->neighbors,
601 FALSE, FALSE, piv_data,
602 &i_mx, &j_mx, &i_my, &j_my);
603 l_data->dx[i][j] = piv_data->dx[i + i_mx][j + j_mx];
604 l_data->dy[i][j] = piv_data->dy[i + i_my][j + j_my];
605 if (l_data->snr[i][j] == SNR_ERR) {
606 l_data->peak_no[i][j] = -1;
607 } else {
608 l_data->peak_no[i][j] = 0;
614 } /* end of parallel region */
616 } else if (valid_par->subst_type ==
617 GPIV_VALID_SUBSTYPE__COR_PEAK) {
619 * substitution with image analyse of next higher correlation
620 * peak (or higher)
622 #pragma omp parallel default(none) /* "for" NOT missing! */ \
623 private (nr_thr, tid, j) \
624 shared (i, piv_data, valid_par, l_data, image ,piv_par, \
625 outlier_found, stdout, \
628 #ifdef ENABLE_OMP
629 tid = omp_get_thread_num();
630 #else
631 tid=0; /* init -1, for single thread set to 0 to address 1st array element */
632 #endif
633 #pragma omp for /* schedule(static) */
634 for (i = 0; i < piv_data->ny; i++) {
635 for (j = 0; j < piv_data->nx; j++) {
636 if (piv_data->peak_no[i][j] != -1
637 && piv_data->snr[i][j] > valid_par->residu_max) {
638 outlier_found[tid] = TRUE;
640 interr_reg_nhpeak (i, j, image, piv_par, l_data, ft,
641 tid);
645 } /* end of parallel region */
650 * Copy back to original values and de-allocate local PIV data set
652 gpiv_ovwrt_pivdata (l_data, piv_data);
653 gpiv_free_pivdata (l_data);
654 /* return outlier_found; */
655 /* check outlier_found[] thread-array here on TRUE value */
656 for (i=0; i < max_nr_thr; i++) {
657 if (outlier_found[i] == TRUE) {
658 free(outlier_found);
659 return TRUE; /* function returns with TRUE on first outlier_found[]==TRUE */
662 free(outlier_found);
663 return FALSE; /* no outlier_found[] element was TRUE */
667 static void
668 cumhisto_eqdatbin (const GpivPivData *data,
669 GpivBinData *klass
671 /*-----------------------------------------------------------------------------
672 * DESCRIPTION:
673 * Calculating cumulative histogram from GpivPivData (NOT from GpivScalarData!)
674 * with an equal number of date per bin of klass
676 * PROTOTYPE LOCATATION:
677 * valid.h
679 * INPUTS:
680 * data: piv dataset
682 * OUTPUTS:
683 * klass: histogram dataset
685 * RETURNS:
686 *---------------------------------------------------------------------------*/
688 gint i, j, k, nresidus = 0;
689 gint nx = data->nx, ny = data->ny, **peak_no = data->peak_no;
690 gfloat **snr = data->snr;
691 gfloat delta, *residu;
692 gfloat *bound = klass->bound, *centre = klass->centre;
693 gdouble fract, yval;
694 gint *count = klass->count, nbins = klass->nbins;
695 gint total_ndata = 0, nresidus_bin = 0;
697 g_assert (data->point_x != NULL);
698 g_assert (data->point_y != NULL);
699 g_assert (data->dx != NULL);
700 g_assert (data->dy != NULL);
701 g_assert (data->snr != NULL);
702 g_assert (data->peak_no != NULL);
704 g_assert (klass->count != NULL);
705 g_assert (klass->bound != NULL);
706 g_assert (klass->centre != NULL);
709 for (i = 0, k = 0; i < ny; i++) {
710 for (j = 0; j < nx; j++) {
711 if (peak_no[i][j] != -1 && snr[i][j] != 0) k++;
714 nresidus = k;
715 nresidus_bin = nresidus / nbins;
717 residu = gpiv_vector(nresidus);
718 for (i = 0, k = 0; i < ny; i++) {
719 for (j = 0; j < nx; j++) {
720 if (peak_no[i][j] != -1 && snr[i][j] != 0) {
721 residu[k] = snr[i][j];
722 k++;
729 * sorting snr data
731 qsort(residu, nresidus, sizeof(float), compare_float);
734 * find lower boundaries of bins
737 for (i = 0; i < nbins; i++) {
738 for (j = 0; j < nresidus_bin; j++) {
739 if (j == 0) {
740 klass->bound[i] = log (1.0/(1.0 - (double) i / (double) nbins));
741 /* klass->bound[i] = (double) i / (double) nbins; */
742 klass->centre[i] = residu[(i+1) * nresidus_bin];
747 klass->min = bound[0];
748 klass->max = bound[nbins - 1];
749 gpiv_free_vector(residu);
754 GpivBinData *
755 gpiv_valid_peaklck (const GpivPivData *piv_data,
756 const guint nbins
758 /*-----------------------------------------------------------------------------
759 * DESCRIPTION:
760 * Calculating histogram of sub-pixel displacements to check on peaklocking
762 * PROTOTYPE LOCATATION:
763 * valid.h
765 * INPUTS:
767 * OUTPUTS:
769 * RETURNS:
770 *---------------------------------------------------------------------------*/
772 GpivBinData *klass = NULL;
774 gchar *err_msg = NULL;
775 gint i, j, k;
776 gfloat delta, fract;
777 gfloat *bound, *centre;
778 gint *count;
781 if ((err_msg = gpiv_check_alloc_pivdata (piv_data)) != NULL) {
782 gpiv_warning ("%s", err_msg);
783 return NULL;
786 if ((klass = gpiv_alloc_bindata (nbins)) == NULL) {
787 gpiv_warning ("gpiv_valid_peaklck: failing gpiv_alloc_bindata");
788 return NULL;
791 count = klass->count;
792 bound = klass->bound;
793 centre = klass->centre;
794 delta = 1. / nbins;
796 for (i = 0; i < nbins; i++) {
797 centre[i] = (float) i *delta;
798 bound[i] = -delta / 2.0 + (float) i *delta;
799 count[i] = 0;
803 * Subdividing fractional particle displacements in bins
805 #pragma omp parallel for shared(piv_data, bound, count) \
806 private(i, j, k, fract)
807 for (i = 0; i < piv_data->ny; i++) {
808 for (j = 0; j < piv_data->nx; j++) {
809 fract = fabs(piv_data->dx[i][j] - (int) piv_data->dx[i][j]);
810 for (k = 0; k < nbins; k++) {
811 if ((fract >= bound[k]) && (fract < bound[k + 1])) {
812 count[k] = count[k] + 1;
815 fract = fabs(piv_data->dy[i][j] - (int) piv_data->dy[i][j]);
816 for (k = 0; k < nbins; k++) {
817 if ((fract >= bound[k]) && (fract < bound[k + 1])) {
818 count[k] = count[k] + 1;
825 return klass;
830 GpivBinData *
831 gpiv_valid_residu_stats (const GpivPivData *piv_data,
832 const guint nbins,
833 GpivLinRegData *linreg
835 /*-----------------------------------------------------------------------------
836 * DESCRIPTION:
837 * Calculates cumulative histogram of residus
839 * PROTOTYPE LOCATATION:
840 * valid.h
842 * INPUTS:
843 * linreg: linear regression data structure
844 * klass: histogram data
846 * OUTPUTS:
847 * piv_data: piv dataset containing residu values in snr
849 * RETURNS:
850 * GpivBinData containing histogram
851 *---------------------------------------------------------------------------*/
853 GpivBinData *klass = NULL;
854 gchar *err_msg = NULL;
855 int i, return_val;
856 double *x, *y;
859 if ((klass = gpiv_alloc_bindata (nbins)) == NULL) {
860 gpiv_warning ("gpiv_valid_residu_stats: failing gpiv_alloc_bindata(%d)", nbins);
863 cumhisto_eqdatbin (piv_data, klass);
864 x = gpiv_dvector (klass->nbins);
865 y = gpiv_dvector (klass->nbins);
867 for (i = 0; i < klass->nbins; i++) {
868 x[i] = (double) klass->bound[i];
869 y[i] = (double) klass->centre[i];
872 if (return_val =
873 gsl_fit_linear (x, 1, y, 1, klass->nbins, &linreg->c0,
874 &linreg->c1, &linreg->cov00, &linreg->cov01,
875 &linreg->cov11, &linreg->sumsq) == 1) {
876 err_msg = "gpiv_valid_residu_stats: error from gsl_fit_linear";
877 g_warning("%s", err_msg);
878 return NULL;
882 gpiv_free_dvector(x);
883 gpiv_free_dvector(y);
884 return klass;
889 gchar *
890 gpiv_valid_errvec (GpivPivData *piv_data,
891 const GpivImage *image,
892 const GpivPivPar *piv_par,
893 const GpivValidPar *valid_par,
894 GpivFt *ft,
895 const gboolean interrogate_valid
897 /*-----------------------------------------------------------------------------
898 * DESCRIPTION:
899 * Searches the erroneous vectors in a PIV data set and
900 * substitutes with new values, if possible
902 * PROTOTYPE LOCATATION:
903 * valid.h
905 * INPUTS:
906 * image: image struct
907 * piv_par: PIV evaluation parameters
908 * valid_par: validation parameters
909 * piv_data: PIV dataset containing particle image displacements
910 * interrogate_valid: validation during (iterative) interrogation process
912 * RETURNS:
913 * GpivPivData on success or NULL on failure
914 *---------------------------------------------------------------------------*/
916 GpivPivData *l_data = NULL;
917 gchar *err_msg = NULL;
919 gint i, j;
920 gboolean outlier_found = TRUE;
921 gint count = 0;
922 enum SubstitutionType l_subst_type = valid_par->subst_type;
923 gint l_peak = piv_par->peak;
924 GpivPivPar *l_piv_par = NULL;
925 GpivValidPar *l_valid_par = NULL;
929 * Checking input pivdata
931 #ifdef DEBUG
932 g_message ("gpiv_valid_errvec:: piv_data nx = %d ny = %d",
933 piv_data->nx, piv_data->ny);
934 #endif
936 if (piv_data->point_x == NULL
937 || piv_data->point_y == NULL
938 || piv_data->dx == NULL
939 || piv_data->dy == NULL
940 || piv_data->snr == NULL
941 || piv_data->snr == NULL
942 || piv_data->peak_no == NULL) {
943 err_msg = "gpiv_valid_errvec: piv_data->* == NULL";
944 gpiv_warning ("%s", err_msg);
945 return err_msg;
949 l_piv_par = gpiv_piv_cp_parameters (piv_par);
950 l_valid_par = gpiv_valid_cp_parameters (valid_par);
953 while (outlier_found && count < GPIV_VALID_MAX_SWEEP) {
956 * Calculates and substitutes snr data with residu values higher than threshold
958 l_data = gpiv_cp_pivdata (piv_data);
959 gpiv_valid_residu (l_data, l_valid_par, TRUE);
960 if (l_valid_par->subst_type != GPIV_VALID_SUBSTYPE__NONE
961 && interrogate_valid) {
963 * Test data with different types of substitutions if used during
964 * (iterative) image interrogation
966 if (count <= 1) {
967 l_valid_par->subst_type = GPIV_VALID_SUBSTYPE__COR_PEAK;
968 l_piv_par->peak = count + 2;
969 } else {
970 l_valid_par->subst_type = l_subst_type;
971 l_piv_par->peak = l_peak;
975 outlier_found = subst_vector (l_data, image, l_piv_par, l_valid_par,
976 ft);
978 if (l_valid_par->subst_type == GPIV_VALID_SUBSTYPE__NONE) {
979 outlier_found = FALSE;
983 * piv_data are updated with corrected values
984 * l_data will be made free for an eventually next loop
985 * BUGFIX: possible memleak? Comment: when outside if {} causes crashing.
987 if (outlier_found) {
988 count++;
989 gpiv_ovwrt_pivdata (l_data, piv_data);
991 gpiv_free_pivdata (l_data);
995 return err_msg;
1000 void
1001 gpiv_valid_gradient (const GpivPivPar *piv_par,
1002 GpivPivData *piv_data
1004 /*-----------------------------------------------------------------------------
1005 * DESCRIPTION:
1006 * Searches vectors in a PIV data set that exceed the maximum gradient
1007 * (dU x dt/int_size > GPIV_GRADIENT_THRESHOLD)
1009 * PROTOTYPE LOCATATION:
1010 * valid.h
1012 * INPUTS:
1013 * piv_par: PIV evaluation parameters
1014 * valid_par: validation parameters
1016 * OUTPUTS:
1017 * piv_data: piv dataset containing peak_no = -1 for exceeded maxima
1019 * RETURNS:
1020 *---------------------------------------------------------------------------*/
1022 int i, j, diff_order = 1;
1023 double grad_dx, delta_dx, grad_dy, delta_dy;
1025 g_return_if_fail (piv_data->point_x != NULL);
1026 g_return_if_fail (piv_data->point_y != NULL);
1027 g_return_if_fail (piv_data->dx != NULL);
1028 g_return_if_fail (piv_data->dy != NULL);
1029 g_return_if_fail (piv_data->snr != NULL);
1030 g_return_if_fail (piv_data->peak_no != NULL);
1032 /* BUGFIX: test op patch. Gerber */
1033 for (i = diff_order; i < piv_data->ny - diff_order; i++) {
1034 for (j = diff_order; j < piv_data->nx - diff_order; j++) {
1036 if(piv_data->peak_no[i][j] != -1
1037 /* && piv_data->peak_no[i-1][j] != -1 && */
1038 /* piv_data->peak_no[i][j-1] != -1 && */
1039 /* piv_data->peak_no[i][j+1] != -1 */
1041 grad_dx = (piv_data->dx[i+1][j] - piv_data->dx[i-1][j]) /
1042 (2 * piv_par->int_shift);
1043 delta_dx = fabs(grad_dx) * piv_par->int_size_f;
1045 piv_data->snr[i][j] = delta_dx;
1046 grad_dy = (piv_data->dy[i][j+1] - piv_data->dy[i][j-1]) /
1047 (2 * piv_par->int_shift);
1048 delta_dy = fabs(grad_dy) * piv_par->int_size_f;
1049 if (delta_dx > GPIV_GRADIENT_THRESHOLD ||
1050 delta_dy > GPIV_GRADIENT_THRESHOLD)
1051 piv_data->peak_no[i][j] = -1;
1059 * exclude all data near the boundaries of the dataset
1061 for (i=0; i < diff_order; i++) {
1062 for (j=0; j < piv_data->nx; j++) {
1063 piv_data->peak_no[i][j] = 0;
1067 for (i=0; i < piv_data->ny; i++) {
1068 for (j=0; j < diff_order; j++) {
1069 piv_data->peak_no[i][j] = 0;
1073 for (i=piv_data->ny - diff_order; i < piv_data->ny; i++) {
1074 for (j=0; j < piv_data->nx; j++) {
1075 piv_data->peak_no[i][j] = 0;
1079 for (i=0; i < piv_data->ny; i++) {
1080 for (j=piv_data->nx - diff_order; j < piv_data->nx; j++) {
1081 piv_data->peak_no[i][j] = 0;
1092 void
1093 gpiv_cumhisto_eqdatbin_gnuplot (const gchar *fname_out,
1094 const gchar *title,
1095 const gchar *GNUPLOT_DISPLAY_COLOR,
1096 const gint GNUPLOT_DISPLAY_SIZE,
1097 const GpivLinRegData *linreg
1099 /*-----------------------------------------------------------------------------
1100 * DESCRIPTION:
1101 * Plots data on screen with gnuplot
1103 * PROTOTYPE LOCATATION:
1104 * valid.h
1106 * INPUTS:
1108 * OUTPUTS:
1110 * RETURNS:
1111 *---------------------------------------------------------------------------*/
1114 FILE *fp_cmd;
1115 gchar *fname_cmd="/tmp/gpiv_gnuplot.cmd";
1116 gchar *function_name="gpiv_histo_gnuplot";
1117 gchar command[GPIV_MAX_CHARS];
1119 if((fp_cmd=fopen(fname_cmd,"w"))==NULL) {
1120 fprintf (stderr,"\n%s:%s error: Failure opening %s for output\n",
1121 LIBNAME, function_name, fname_cmd);
1122 exit(1);
1125 fprintf (fp_cmd,"\nset xlabel \"-ln(1-i/nbins)\"");
1126 fprintf (fp_cmd,"\nset ylabel \"residu (pixels)\"");
1127 fprintf (fp_cmd,"\nplot \"%s\" title \"%s\" with boxes, %f + %f * x", /* with boxes */
1128 fname_out, title, linreg->c0, linreg->c1);
1129 fprintf (fp_cmd,"\npause -1 \"Hit return to exit\"");
1130 fprintf (fp_cmd,"\nquit");
1132 fclose (fp_cmd);
1135 snprintf(command, GPIV_MAX_CHARS, "gnuplot -bg %s -geometry %dx%d %s",
1136 GNUPLOT_DISPLAY_COLOR, GNUPLOT_DISPLAY_SIZE,
1137 GNUPLOT_DISPLAY_SIZE, fname_cmd);
1139 if (system (command) != 0) {
1140 fprintf (stderr,"\n%s:%s could not exec shell command\n",
1141 LIBNAME, function_name);
1142 exit(1);
1149 gfloat
1150 gpiv_valid_threshold (const GpivPivPar *piv_par,
1151 const GpivValidPar *valid_par,
1152 const GpivLinRegData *linreg
1154 /*-----------------------------------------------------------------------------
1155 * DESCRIPTION:
1156 * Calculates threshold value (residu_max) from residus.
1157 * Will need the int_size_f from the GpivPivPar struct!
1159 * PROTOTYPE LOCATATION:
1160 * valid.h
1162 * INPUTS:
1164 * OUTPUTS:
1166 * RETURNS:
1167 *---------------------------------------------------------------------------*/
1169 gfloat search_radius = (float) piv_par->int_size_f / 4.0;
1170 gfloat residu_max = - linreg->c1 * log((1.0 - valid_par->data_yield) /
1171 valid_par->data_yield *
1172 linreg->c1 / search_radius);
1173 return residu_max;
1176 #undef SNR_ERR
1177 #undef MIN_VECTORS4MEDIAN
1178 #undef RESIDU_EPSI