SLIP  1.4
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
correlation.hpp
Go to the documentation of this file.
1 /*
2  * Copyright(c):
3  * Signal Image and Communications (SIC) Department
4  * http://www.sic.sp2mi.univ-poitiers.fr/
5  * - University of Poitiers, France http://www.univ-poitiers.fr
6  * - XLIM Institute UMR CNRS 7252 http://www.xlim.fr/
7  *
8  * and
9  *
10  * D2 Fluid, Thermic and Combustion
11  * - University of Poitiers, France http://www.univ-poitiers.fr
12  * - PPRIME Institute - UPR CNRS 3346 http://www.pprime.fr
13  * - ISAE-ENSMA http://www.ensma.fr
14  *
15  * Contributor(s):
16  * The SLIP team,
17  * Benoit Tremblais <tremblais_AT_sic.univ-poitiers.fr>,
18  * Laurent David <laurent.david_AT_lea.univ-poitiers.fr>,
19  * Ludovic Chatellier <ludovic.chatellier_AT_univ-poitiers.fr>,
20  * Lionel Thomas <lionel.thomas_AT_univ-poitiers.fr>,
21  * Denis Arrivault <arrivault_AT_sic.univ-poitiers.fr>,
22  * Julien Dombre <julien.dombre_AT_univ-poitiers.fr>.
23  *
24  * Description:
25  * The Simple Library of Image Processing (SLIP) is a new image processing
26  * library. It is written in the C++ language following as much as possible
27  * the ISO/ANSI C++ standard. It is consequently compatible with any system
28  * satisfying the ANSI C++ complience. It works on different Unix , Linux ,
29  * Mircrosoft Windows and Mac OS X plateforms. SLIP is a research library that
30  * was created by the Signal, Image and Communications (SIC) departement of
31  * the XLIM, UMR 7252 CNRS Institute in collaboration with the Fluids, Thermic
32  * and Combustion departement of the P', UPR 3346 CNRS Institute of the
33  * University of Poitiers.
34  *
35  * The SLIP Library source code has been registered to the APP (French Agency
36  * for the Protection of Programs) by the University of Poitiers and CNRS,
37  * under registration number IDDN.FR.001.300034.000.S.P.2010.000.21000.
38 
39  * http://www.sic.sp2mi.univ-poitiers.fr/slip/
40  *
41  * This software is governed by the CeCILL-C license under French law and
42  * abiding by the rules of distribution of free software. You can use,
43  * modify and/ or redistribute the software under the terms of the CeCILL-C
44  * license as circulated by CEA, CNRS and INRIA at the following URL
45  * http://www.cecill.info.
46  * As a counterpart to the access to the source code and rights to copy,
47  * modify and redistribute granted by the license, users are provided only
48  * with a limited warranty and the software's author, the holder of the
49  * economic rights, and the successive licensors have only limited
50  * liability.
51  *
52  * In this respect, the user's attention is drawn to the risks associated
53  * with loading, using, modifying and/or developing or reproducing the
54  * software by the user in light of its specific status of free software,
55  * that may mean that it is complicated to manipulate, and that also
56  * therefore means that it is reserved for developers and experienced
57  * professionals having in-depth computer knowledge. Users are therefore
58  * encouraged to load and test the software's suitability as regards their
59  * requirements in conditions enabling the security of their systems and/or
60  * data to be ensured and, more generally, to use and operate it in the
61  * same conditions as regards security.
62  *
63  * The fact that you are presently reading this means that you have had
64  * knowledge of the CeCILL-C license and that you accept its terms.
65  */
66 
67 
75 #ifndef SLIP_CORRELATION_HPP
76 #define SLIP_CORRELATION_HPP
77 
78 
79 #include <algorithm>
80 #include <numeric>
81 #include <cassert>
82 #include <cmath>
83 #include <complex>
84 #include "statistics.hpp"
85 #include "FFT.hpp"
86 #include "Array.hpp"
87 #include "error.hpp"
88 #include "Array3d.hpp"
89 #include "convolution.hpp"
90 #include "linear_algebra.hpp"
92 #include "arithmetic_op.hpp"
93 #include "complex_cast.hpp"
94 
95 
96 namespace slip
97 {
98 
99 
100 
122  {
123  CC,
127  };
128 
129 
131  /* @{ */
132 
178  template<typename T, typename InputIterator1, typename InputIterator2>
179  inline
180  T std_crosscorrelation(InputIterator1 first, InputIterator1 last,
181  InputIterator2 first2)
182  {
183  return std::inner_product(first,last,first2,T());
184  }
185 
218  template<typename T,
219  typename InputIterator1,
220  typename InputIterator2,
221  typename MaskIterator>
222  inline
223  T std_crosscorrelation_mask(InputIterator1 first1,
224  InputIterator1 last1,
225  MaskIterator mask_first,
226  InputIterator2 first2,
227  typename std::iterator_traits<MaskIterator>::value_type value=typename
228  std::iterator_traits<MaskIterator>::value_type(1))
229 {
230  assert(first1 != last1);
231  T init = T(0);
232  for (; first1 != last1; ++first1,++first2,++mask_first)
233  {
234  if(*mask_first == value)
235  {
236  init +=(*first1)*(*first2);
237  }
238  }
239  return init;
240 }
241 
273 template<typename T,
274  typename InputIterator1,
275  typename InputIterator2,
276  typename Predicate>
277 
278  T std_crosscorrelation_if(InputIterator1 first1,
279  InputIterator1 last1,
280  InputIterator2 first2,
281  Predicate pred)
282  {
283  assert(first1 != last1);
284  T init = T(0);
285  for (; first1 != last1; ++first1, ++first2)
286  {
287  if(pred(*first1))
288  {
289  init +=(*first1)*(*first2);
290  }
291  }
292  return init;
293  }
294 
337  template<typename T, typename InputIterator1, typename InputIterator2>
338  inline
339  T cc(InputIterator1 first, InputIterator1 last,
340  InputIterator2 first2)
341  {
342  return slip::std_crosscorrelation<T>(first,last,first2);
343  }
344 
377 template<typename T,
378  typename InputIterator1,
379  typename InputIterator2,
380  typename MaskIterator>
381  inline
382  T cc_mask(InputIterator1 first1,
383  InputIterator1 last1,
384  MaskIterator mask_first,
385  InputIterator2 first2,
386 
387  typename std::iterator_traits<MaskIterator>::value_type value=typename
388  std::iterator_traits<MaskIterator>::value_type(1))
389 {
390  return slip::std_crosscorrelation_mask<T>(first1,last1,mask_first,first2,value);
391 }
392 
424 template<typename T,
425  typename InputIterator1,
426  typename InputIterator2,
427  typename Predicate>
428 
429  T cc_if(InputIterator1 first1,
430  InputIterator1 last1,
431  InputIterator2 first2,
432  Predicate pred)
433 {
434  return slip::std_crosscorrelation_if<T>(first1,last1,first2,pred);
435 }
436 
437 
438 
439 
486  template<typename Real, typename InputIterator1, typename InputIterator2>
487  inline
488  Real centered_crosscorrelation(InputIterator1 first, InputIterator1 last,
489  InputIterator2 first2,
490  typename std::iterator_traits<InputIterator1>::value_type mean1 = typename std::iterator_traits<InputIterator1>::value_type(),
491  typename std::iterator_traits<InputIterator1>::value_type mean2 = typename std::iterator_traits<InputIterator2>::value_type())
492  {
493 
494  Real __init1 = Real();
495 
496  for (; first != last; ++first, ++first2)
497  {
498  __init1 = __init1 + (*first - mean1) * (*first2 - mean2);
499  }
500 
501  return __init1;
502  }
503 
537 template<typename Real,
538  typename InputIterator1,
539  typename InputIterator2,
540  typename MaskIterator>
541  inline
542  Real centered_crosscorrelation_mask(InputIterator1 first1,InputIterator1 last1,
543  MaskIterator mask_first,
544  InputIterator2 first2,
545 
546  typename std::iterator_traits<InputIterator1>::value_type mean1 = typename std::iterator_traits<InputIterator1>::value_type(),
547  typename std::iterator_traits<InputIterator1>::value_type mean2 = typename std::iterator_traits<InputIterator2>::value_type(),
548  typename std::iterator_traits<MaskIterator>::value_type value=typename
549  std::iterator_traits<MaskIterator>::value_type(1))
550 {
551  Real __init1 = Real();
552  for (; first1 != last1; ++first1, ++first2,++mask_first)
553  {
554  if(*mask_first == value)
555  {
556  __init1 = __init1 + (*first1 - mean1) * (*first2 - mean2);
557  }
558  }
559  return __init1;
560 }
561 
596 template<typename Real,
597  typename InputIterator1,
598  typename InputIterator2,
599  typename Predicate>
600 inline
601 Real centered_crosscorrelation_if(InputIterator1 first1,InputIterator1 last1,
602  InputIterator2 first2,
603  Predicate pred,
604  typename std::iterator_traits<InputIterator1>::value_type mean1 = typename std::iterator_traits<InputIterator1>::value_type(),
605  typename std::iterator_traits<InputIterator1>::value_type mean2 = typename std::iterator_traits<InputIterator2>::value_type())
606 {
607 
608 
609  Real __init1 = Real();
610  for (; first1 != last1; ++first1, ++first2)
611  {
612  if(pred(*first1))
613  {
614  __init1 = __init1 + (*first1 - mean1) * (*first2 - mean2);
615  }
616  }
617  return __init1;
618 }
619 
665  template<typename Real,
666  typename InputIterator1,
667  typename InputIterator2>
668  inline
669  Real ccc(InputIterator1 first, InputIterator1 last,
670  InputIterator2 first2,
671  typename std::iterator_traits<InputIterator1>::value_type mean1 = typename std::iterator_traits<InputIterator1>::value_type(),
672  typename std::iterator_traits<InputIterator1>::value_type mean2 = typename std::iterator_traits<InputIterator2>::value_type())
673  {
674  return slip::centered_crosscorrelation<Real>(first,last,first2,mean1,mean2);
675  }
676 
710 template<typename Real,
711  typename InputIterator1,
712  typename InputIterator2,
713  typename MaskIterator>
714  inline
715  Real ccc_mask(InputIterator1 first1,InputIterator1 last1,
716  MaskIterator mask_first,
717  InputIterator2 first2,
718 
719  typename std::iterator_traits<InputIterator1>::value_type mean1 = typename std::iterator_traits<InputIterator1>::value_type(),
720  typename std::iterator_traits<InputIterator1>::value_type mean2 = typename std::iterator_traits<InputIterator2>::value_type(),
721  typename std::iterator_traits<MaskIterator>::value_type value=typename
722  std::iterator_traits<MaskIterator>::value_type(1))
723 
724 {
725  return slip::centered_crosscorrelation_mask<Real>(first1,last1,mask_first,first2,mean1,mean2,value);
726 }
727 
763 template<typename Real,
764  typename InputIterator1,
765  typename InputIterator2,
766  typename Predicate>
767  inline
768  Real ccc_if(InputIterator1 first1,InputIterator1 last1,
769  InputIterator2 first2,Predicate pred,
770  typename std::iterator_traits<InputIterator1>::value_type mean1 = typename std::iterator_traits<InputIterator1>::value_type(),
771  typename std::iterator_traits<InputIterator1>::value_type mean2 = typename std::iterator_traits<InputIterator2>::value_type())
772 {
773  return slip::centered_crosscorrelation_if<Real>(first1,last1,first2,pred,mean1,mean2);
774 }
775 
776 
824  template<typename Real,
825  typename InputIterator1,
826  typename InputIterator2>
827  inline
828  Real pseudo_normalized_crosscorrelation(InputIterator1 first, InputIterator1 last,
829  InputIterator2 first2,
830  typename std::iterator_traits<InputIterator1>::value_type mean1 = typename std::iterator_traits<InputIterator1>::value_type(),
831  typename std::iterator_traits<InputIterator1>::value_type mean2 = typename std::iterator_traits<InputIterator2>::value_type())
832  {
833 
834  Real __init1 = Real();
835  Real __init2 = Real();
836  Real __cc = Real();
837  for (; first != last; ++first, ++first2)
838  {
839  __init1 = __init1 + (*first - mean1) * (*first - mean1);
840  __init2 = __init2 + (*first2 - mean2) * (*first2 - mean2);
841  __cc = __cc + *first * *first2;
842  }
843 
844  return __cc / std::sqrt(__init1 * __init2);
845  }
881 template<typename Real, typename InputIterator1, typename InputIterator2,typename MaskIterator >
882  inline
883  Real pseudo_normalized_crosscorrelation_mask(InputIterator1 first1,
884  InputIterator1 last1,
885  MaskIterator mask_first,
886  InputIterator2 first2,
887  typename std::iterator_traits<InputIterator1>::value_type mean1 = typename std::iterator_traits<InputIterator1>::value_type(),
888  typename std::iterator_traits<InputIterator1>::value_type mean2 = typename std::iterator_traits<InputIterator2>::value_type(),typename std::iterator_traits<MaskIterator>::value_type value=typename
889  std::iterator_traits<MaskIterator>::value_type(1))
890  {
891  Real __init1 = Real();
892  Real __init2 = Real();
893  Real __cc = Real();
894  for (; first1 != last1; ++first1, ++first2,++mask_first)
895  {
896  if(*mask_first == value)
897  {
898 
899  Real __diff1 = (*first1 - mean1);
900  Real __diff2 = (*first2 - mean2);
901  __init1 = __init1 + (__diff1 * __diff1);
902  __init2 = __init2 + (__diff2 * __diff2);
903  __cc = __cc + (*first1 * *first2);
904  }
905  }
906  return __cc / std::sqrt(__init1 * __init2);
907  }
908 
909 
945 template<typename Real, typename InputIterator1, typename InputIterator2,typename Predicate>
946  inline
947  Real pseudo_normalized_crosscorrelation_if(InputIterator1 first1,InputIterator1 last1,
948  InputIterator2 first2,Predicate pred,
949  typename std::iterator_traits<InputIterator1>::value_type mean1 = typename std::iterator_traits<InputIterator1>::value_type(),
950  typename std::iterator_traits<InputIterator1>::value_type mean2 = typename std::iterator_traits<InputIterator2>::value_type())
951 
952  {
953 
954  Real __init1 = Real();
955  Real __init2 = Real();
956  Real __cc = Real();
957 
958  for (; first1 != last1; ++first1, ++first2)
959  {
960  if(pred(*first1))
961  {
962  Real __diff1 = (*first1 - mean1);
963  Real __diff2 = (*first2 - mean2);
964  __init1 = __init1 + (__diff1 * __diff1);
965  __init2 = __init2 + (__diff2 * __diff2);
966  __cc = __cc + (*first1 * *first2);
967  }
968  }
969  return __cc / std::sqrt(__init1 * __init2);
970  }
971 
972 
1020  template<typename Real, typename InputIterator1, typename InputIterator2>
1021  inline
1022  Real pncc(InputIterator1 first, InputIterator1 last,
1023  InputIterator2 first2,
1024  typename std::iterator_traits<InputIterator1>::value_type mean1 = typename std::iterator_traits<InputIterator1>::value_type(),
1025  typename std::iterator_traits<InputIterator1>::value_type mean2 = typename std::iterator_traits<InputIterator2>::value_type())
1026  {
1027 
1028  return slip::pseudo_normalized_crosscorrelation<Real>(first,last,first2,mean1,mean2);
1029  }
1030 
1064 template<typename Real,
1065  typename InputIterator1,
1066  typename InputIterator2,
1067  typename MaskIterator >
1068  inline
1069  Real pncc_mask(InputIterator1 first1,
1070  InputIterator1 last1,
1071  MaskIterator mask_first,
1072  InputIterator2 first2,
1073  typename std::iterator_traits<InputIterator1>::value_type mean1 = typename std::iterator_traits<InputIterator1>::value_type(),
1074  typename std::iterator_traits<InputIterator1>::value_type mean2 = typename std::iterator_traits<InputIterator2>::value_type(),typename std::iterator_traits<MaskIterator>::value_type value=typename
1075  std::iterator_traits<MaskIterator>::value_type(1))
1076  {
1077  return slip::pseudo_normalized_crosscorrelation_mask<Real>(first1,last1,mask_first,first2,mean1,mean2,value);
1078  }
1079 
1080 
1115 template<typename Real,
1116  typename InputIterator1,
1117  typename InputIterator2,
1118  typename Predicate>
1119 inline
1120 Real pncc_if(InputIterator1 first1,InputIterator1 last1,
1121  InputIterator2 first2,Predicate pred,
1122  typename std::iterator_traits<InputIterator1>::value_type mean1 = typename std::iterator_traits<InputIterator1>::value_type(),
1123  typename std::iterator_traits<InputIterator1>::value_type mean2 = typename std::iterator_traits<InputIterator2>::value_type())
1124 
1125  {
1126  return slip::pseudo_normalized_crosscorrelation_if<Real>(first1,last1,first2,pred,mean1,mean2);
1127  }
1175  template<typename Real,
1176  typename InputIterator1,
1177  typename InputIterator2>
1178  inline
1179  Real normalized_crosscorrelation(InputIterator1 first,
1180  InputIterator1 last,
1181  InputIterator2 first2,
1182  typename std::iterator_traits<InputIterator1>::value_type mean1 = typename std::iterator_traits<InputIterator1>::value_type(),
1183  typename std::iterator_traits<InputIterator1>::value_type mean2 = typename std::iterator_traits<InputIterator2>::value_type())
1184  {
1185 
1186  //computation of the normalization term
1187  Real __init0 = Real();
1188  Real __init1 = Real();
1189  Real __init2 = Real();
1190 
1191  for (; first != last; ++first, ++first2)
1192  {
1193  Real __diff1 = (*first - mean1);
1194  Real __diff2 = (*first2 - mean2);
1195  __init0 = __init0 + __diff1 * __diff2;
1196  __init1 = __init1 + __diff1 * __diff1;
1197  __init2 = __init2 + __diff2 * __diff2;
1198  }
1199 
1200  return __init0 / std::sqrt(__init1 * __init2);
1201  }
1202 
1236 template<typename Real, typename InputIterator1, typename InputIterator2,typename MaskIterator >
1237  inline
1238  Real normalized_crosscorrelation_mask(InputIterator1 first1,
1239  InputIterator1 last1,
1240  MaskIterator mask_first,
1241  InputIterator2 first2,
1242 
1243  typename std::iterator_traits<InputIterator1>::value_type mean1 = typename std::iterator_traits<InputIterator1>::value_type(),
1244  typename std::iterator_traits<InputIterator1>::value_type mean2 = typename std::iterator_traits<InputIterator2>::value_type(),typename std::iterator_traits<MaskIterator>::value_type value=typename
1245  std::iterator_traits<MaskIterator>::value_type(1))
1246  {
1247 
1248  //computation of the normalization term
1249  Real __init0 = Real();
1250  Real __init1 = Real();
1251  Real __init2 = Real();
1252 
1253  for (; first1 != last1; ++first1, ++first2,++mask_first)
1254  {
1255  if(*mask_first == value)
1256  {
1257  Real __diff1 = (*first1 - mean1);
1258  Real __diff2 = (*first2 - mean2);
1259  __init0 = __init0 + __diff1 * __diff2;
1260  __init1 = __init1 + __diff1 * __diff1;
1261  __init2 = __init2 + __diff2 * __diff2;
1262  }
1263  }
1264  return __init0 / std::sqrt(__init1 * __init2);
1265  }
1266 
1300 template<typename Real,
1301  typename InputIterator1,
1302  typename InputIterator2,
1303  typename Predicate>
1304  inline
1305  Real normalized_crosscorrelation_if(InputIterator1 first1,InputIterator1 last1,
1306  InputIterator2 first2,Predicate pred,
1307  typename std::iterator_traits<InputIterator1>::value_type mean1 = typename std::iterator_traits<InputIterator1>::value_type(),
1308  typename std::iterator_traits<InputIterator1>::value_type mean2 = typename std::iterator_traits<InputIterator2>::value_type())
1309 
1310  {
1311 
1312  //computation of the normalization term
1313  Real __init0 = Real();
1314  Real __init1 = Real();
1315  Real __init2 = Real();
1316 
1317  for (; first1 != last1; ++first1, ++first2)
1318  {
1319  if(pred(*first1))
1320  {
1321  Real __diff1 = (*first1 - mean1);
1322  Real __diff2 = (*first2 - mean2);
1323  __init0 = __init0 + __diff1 * __diff2;
1324  __init1 = __init1 + __diff1 * __diff1;
1325  __init2 = __init2 + __diff2 * __diff2;
1326  }
1327  }
1328  return __init0 / std::sqrt(__init1 * __init2);
1329  }
1330 
1376  template<typename Real, typename InputIterator1, typename InputIterator2>
1377  inline
1378  Real ncc(InputIterator1 first,
1379  InputIterator1 last,
1380  InputIterator2 first2,
1381  typename std::iterator_traits<InputIterator1>::value_type mean1 = typename std::iterator_traits<InputIterator1>::value_type(),
1382  typename std::iterator_traits<InputIterator1>::value_type mean2 = typename std::iterator_traits<InputIterator2>::value_type())
1383  {
1384  return slip::normalized_crosscorrelation<Real>(first,last,first2,mean1,mean2);
1385  }
1386 
1424 template<typename Real, typename InputIterator1, typename InputIterator2,typename MaskIterator >
1425  inline
1426  Real ncc_mask(InputIterator1 first1,
1427  InputIterator1 last1,
1428  MaskIterator mask_first,
1429  InputIterator2 first2,
1430 
1431  typename std::iterator_traits<InputIterator1>::value_type mean1 = typename std::iterator_traits<InputIterator1>::value_type(),
1432  typename std::iterator_traits<InputIterator1>::value_type mean2 = typename std::iterator_traits<InputIterator2>::value_type(),typename std::iterator_traits<MaskIterator>::value_type value=typename
1433  std::iterator_traits<MaskIterator>::value_type(1))
1434 
1435 {
1436  return slip::normalized_crosscorrelation_mask<Real>(first1,last1,mask_first,first2,mean1,mean2,value);
1437 }
1438 
1471 template<typename Real, typename InputIterator1, typename InputIterator2,typename Predicate>
1472  inline
1473  Real ncc_if(InputIterator1 first1,InputIterator1 last1,
1474  InputIterator2 first2,Predicate pred,
1475  typename std::iterator_traits<InputIterator1>::value_type mean1 = typename std::iterator_traits<InputIterator1>::value_type(),
1476  typename std::iterator_traits<InputIterator1>::value_type mean2 = typename std::iterator_traits<InputIterator2>::value_type())
1477 
1478  {
1479  return slip::centered_crosscorrelation_if<Real>(first1,last1,first2,pred,mean1,mean2);
1480 
1481  }
1482 
1483  /* @} */
1484 
1486  /* @{ */
1516  template<typename RandomAccessIterator1, typename RandomAccessIterator2>
1517  inline
1518  void autocorrelation_full(RandomAccessIterator1 first,
1519  RandomAccessIterator1 last,
1520  RandomAccessIterator2 auto_first,
1521  RandomAccessIterator2 auto_last)
1522 
1523  {
1524 
1525  //computation of the means
1526  typedef typename std::iterator_traits<RandomAccessIterator2>::value_type value_type;
1527  typedef typename std::iterator_traits<RandomAccessIterator2>::difference_type _Difference;
1528 
1529  assert(last - first > 0);
1530  assert((auto_last - auto_first) == (last-first)*2 - 1);
1531  _Difference med = _Difference(last - first) - 1;
1532  for(_Difference i = 0; i <= med; ++i)
1533  {
1534  value_type __init0 = value_type();
1535  RandomAccessIterator1 first_tmp = first;
1536  RandomAccessIterator2 first_tmp2 = first + i;
1537  for (; first_tmp2 != last; ++first_tmp, ++first_tmp2)
1538  {
1539  value_type __diff1 = slip::conj(*first_tmp);
1540  value_type __diff2 = *first_tmp2;
1541  __init0 = __init0 + __diff1 * __diff2;
1542  }
1543  value_type val = __init0;
1544 
1545  *(auto_first + med + i) = val;
1546  *(auto_first + med - i) = val;
1547  }
1548  }
1549 
1579  template<typename RandomAccessIterator1, typename RandomAccessIterator2>
1580  inline
1581  void autocorrelation(RandomAccessIterator1 first,
1582  RandomAccessIterator1 last,
1583  RandomAccessIterator2 auto_first,
1584  RandomAccessIterator2 auto_last)
1585 
1586  {
1587 
1588  typedef typename std::iterator_traits<RandomAccessIterator2>::value_type value_type;
1589  typedef typename std::iterator_traits<RandomAccessIterator2>::difference_type _Difference;
1590 
1591  assert(last - first > 0);
1592  assert((auto_last - auto_first) == (last-first));
1593  _Difference med = _Difference(last - first);
1594  for(_Difference i = 0; i < med; ++i)
1595  {
1596  value_type __init0 = value_type();
1597  RandomAccessIterator1 first_tmp = first;
1598  RandomAccessIterator2 first_tmp2 = first + i;
1599  for (; first_tmp2 != last; ++first_tmp, ++first_tmp2)
1600  {
1601  value_type __diff1 = slip::conj(*first_tmp);
1602  value_type __diff2 = *first_tmp2;
1603  __init0 = __init0 + __diff1 * __diff2;
1604  }
1605  value_type val = __init0;
1606 
1607  *(auto_first + i) = val;
1608  }
1609  }
1610 
1640  template<typename RandomAccessIterator1, typename RandomAccessIterator2>
1641  inline
1642  void biased_autocorrelation(RandomAccessIterator1 first,
1643  RandomAccessIterator1 last,
1644  RandomAccessIterator2 auto_first,
1645  RandomAccessIterator2 auto_last)
1646 
1647  {
1648  typedef typename std::iterator_traits<RandomAccessIterator2>::value_type value_type;
1649  typedef typename std::iterator_traits<RandomAccessIterator2>::difference_type _Difference;
1650  slip::autocorrelation(first,last,auto_first,auto_last);
1651 
1652  _Difference N = (auto_last - auto_first);
1653  value_type one_over_N = value_type(1)/value_type(N);
1654  std::transform(auto_first,auto_last,auto_first,std::bind2nd(std::multiplies<value_type>(),one_over_N));
1655 
1656  }
1657 
1687  template<typename RandomAccessIterator1, typename RandomAccessIterator2>
1688  inline
1689  void unbiased_autocorrelation(RandomAccessIterator1 first,
1690  RandomAccessIterator1 last,
1691  RandomAccessIterator2 auto_first,
1692  RandomAccessIterator2 auto_last)
1693 
1694  {
1695 
1696  typedef typename std::iterator_traits<RandomAccessIterator2>::value_type value_type;
1697  typedef typename std::iterator_traits<RandomAccessIterator2>::difference_type _Difference;
1698 
1699  assert(last - first > 0);
1700  assert((auto_last - auto_first) == (last-first));
1701  _Difference med = _Difference(last - first);
1702  for(_Difference i = 0; i < med; ++i)
1703  {
1704  value_type __init0 = value_type();
1705  RandomAccessIterator1 first_tmp = first;
1706  RandomAccessIterator2 first_tmp2 = first + i;
1707  for (; first_tmp2 != last; ++first_tmp, ++first_tmp2)
1708  {
1709  value_type __diff1 = slip::conj(*first_tmp);
1710  value_type __diff2 = *first_tmp2;
1711  __init0 = __init0 + __diff1 * __diff2;
1712  }
1713  value_type val = __init0;
1714 
1715  *(auto_first + i) = val / value_type(med-i);
1716  }
1717  }
1718 
1748  template<typename RandomAccessIterator1, typename RandomAccessIterator2>
1749  inline
1750  void autocovariance(RandomAccessIterator1 first,
1751  RandomAccessIterator1 last,
1752  RandomAccessIterator2 auto_first,
1753  RandomAccessIterator2 auto_last)
1754 
1755  {
1756 
1757  typedef typename std::iterator_traits<RandomAccessIterator2>::value_type value_type;
1758  typedef typename std::iterator_traits<RandomAccessIterator2>::difference_type _Difference;
1759 
1760  assert(last - first > 0);
1761  assert((auto_last - auto_first) == (last-first));
1762  //coputes the mean
1763  value_type mean = slip::mean<value_type>(first,last);
1764  _Difference med = _Difference(last - first);
1765  for(_Difference i = 0; i < med; ++i)
1766  {
1767  value_type __init0 = value_type();
1768  RandomAccessIterator1 first_tmp = first;
1769  RandomAccessIterator2 first_tmp2 = first + i;
1770  for (; first_tmp2 != last; ++first_tmp, ++first_tmp2)
1771  {
1772  value_type __diff1 = (slip::conj(*first_tmp) - mean);
1773  value_type __diff2 = (*first_tmp2 - mean);
1774  __init0 = __init0 + __diff1 * __diff2;
1775  }
1776  value_type val = __init0;
1777 
1778  *(auto_first + i) = val;
1779  }
1780  }
1781 
1811  template<typename RandomAccessIterator1, typename RandomAccessIterator2>
1812  inline
1813  void biased_autocovariance(RandomAccessIterator1 first,
1814  RandomAccessIterator1 last,
1815  RandomAccessIterator2 auto_first,
1816  RandomAccessIterator2 auto_last)
1817 
1818  {
1819  typedef typename std::iterator_traits<RandomAccessIterator2>::value_type value_type;
1820  typedef typename std::iterator_traits<RandomAccessIterator2>::difference_type _Difference;
1821  slip::autocovariance(first,last,auto_first,auto_last);
1822 
1823  _Difference N = (auto_last - auto_first);
1824  value_type one_over_N = value_type(1)/value_type(N);
1825  std::transform(auto_first,auto_last,auto_first,std::bind2nd(std::multiplies<value_type>(),one_over_N));
1826  typedef typename std::iterator_traits<RandomAccessIterator2>::value_type value_type;
1827 
1828  }
1858  template<typename RandomAccessIterator1, typename RandomAccessIterator2>
1859  inline
1860  void unbiased_autocovariance(RandomAccessIterator1 first,
1861  RandomAccessIterator1 last,
1862  RandomAccessIterator2 auto_first,
1863  RandomAccessIterator2 auto_last)
1864 
1865  {
1866  typedef typename std::iterator_traits<RandomAccessIterator2>::value_type value_type;
1867  typedef typename std::iterator_traits<RandomAccessIterator2>::difference_type _Difference;
1868 
1869  assert(last - first > 0);
1870  assert((auto_last - auto_first) == (last-first));
1871  //coputes the mean
1872  value_type mean = slip::mean<value_type>(first,last);
1873  _Difference med = _Difference(last - first);
1874  for(_Difference i = 0; i < med; ++i)
1875  {
1876  value_type __init0 = value_type();
1877  RandomAccessIterator1 first_tmp = first;
1878  RandomAccessIterator2 first_tmp2 = first + i;
1879  for (; first_tmp2 != last; ++first_tmp, ++first_tmp2)
1880  {
1881  value_type __diff1 = (slip::conj(*first_tmp) - mean);
1882  value_type __diff2 = (*first_tmp2 - mean);
1883  __init0 = __init0 + __diff1 * __diff2;
1884  }
1885  value_type val = __init0;
1886 
1887  *(auto_first + i) = val / value_type(med-i);
1888  }
1889  }
1890 
1912  template <typename RandomAccessIterator, typename Matrix>
1913  inline
1914  void biased_autocorrelation_matrix(RandomAccessIterator first,
1915  RandomAccessIterator last,
1916  Matrix& A)
1917  {
1918  typedef typename std::iterator_traits<RandomAccessIterator>::value_type value_type;
1919 
1920  assert(A.rows() == A.cols());
1921  assert(int(last -first) == int(A.rows()));
1922  slip::Array<value_type> Auto(A.rows());
1923  slip::biased_autocorrelation(first,last,
1924  Auto.begin(),Auto.end());
1925  slip::toeplitz(Auto.begin(),Auto.end(),A);
1926 
1927  }
1928 
1929 
1951  template <typename RandomAccessIterator, typename Matrix>
1952  inline
1953  void biased_autocovariance_matrix(RandomAccessIterator first,
1954  RandomAccessIterator last,
1955  Matrix& A)
1956  {
1957  typedef typename std::iterator_traits<RandomAccessIterator>::value_type value_type;
1958 
1959  // _Difference size = last -first;
1960  assert(A.rows() == A.cols());
1961  assert(int(last -first) == int(A.rows()));
1962  slip::Array<value_type> Auto(A.rows());
1963  slip::biased_autocovariance(first,last,
1964  Auto.begin(),Auto.end());
1965  slip::toeplitz(Auto.begin(),Auto.end(),A);
1966 
1967  }
1968 
1969  /* @} */
1970 
1972  /* @{ */
1973 
2010  template<typename Real, typename InputIterator2d1, typename InputIterator2d2, typename OutputIterator2d>
2011  inline
2012  void crosscorrelation2d(InputIterator2d1 in1_upper_left, InputIterator2d1 in1_bottom_right,
2013  InputIterator2d2 in2_upper_left, InputIterator2d2 in2_bottom_right,
2014  OutputIterator2d out_upper_left, OutputIterator2d out_bottom_right,
2016  {
2017  typename InputIterator2d1::difference_type size2din1 = in1_bottom_right - in1_upper_left;
2018  typename InputIterator2d2::difference_type size2din2 = in2_bottom_right - in2_upper_left;
2019  typename OutputIterator2d::difference_type size2dout = out_bottom_right - out_upper_left;
2020  assert(size2din1 == size2dout);
2021  typename OutputIterator2d::difference_type size2d(size2din1[0] + size2din2[0] - 1, size2din1[1] + size2din2[1] - 1);
2022  //Input containers resized
2023  slip::Array2d<Real> I(size2d[0],size2d[1],Real());
2024  slip::Box2d<int> box1(size2din2[0]/2,size2din2[1]/2,size2din2[0]/2 + size2din1[0]-1,size2din2[1]/2 + size2din1[1]-1);
2025  std::copy(in1_upper_left,in1_bottom_right,I.upper_left(box1));
2026  if(t==NCC)
2027  {
2028  slip::Box2d<int> box2;
2030  typename InputIterator2d2::value_type mean2 = slip::mean<typename InputIterator2d2::value_type>(in2_upper_left,in2_bottom_right);
2031  for(size_t i=0; i<(size_t)size2dout[0]; i++)
2032  for(size_t j=0; j<(size_t)size2dout[1]; j++)
2033  {
2034  box2.set_coord(i,j,size2din2[0]-1+i,size2din2[1]-1+j);
2035  typename InputIterator2d1::value_type mean1 = slip::mean<typename InputIterator2d1::value_type>(I.upper_left(box2),I.bottom_right(box2));
2036  d.set_coord(i,j);
2037  *(out_upper_left + d) = normalized_crosscorrelation<Real>(in2_upper_left,in2_bottom_right,I.upper_left(box2),mean2,mean1);
2038  }
2039  }
2040  else if(t==PNCC)
2041  {
2042  slip::Box2d<int> box2;
2044  typename InputIterator2d2::value_type mean2 = slip::mean<typename InputIterator2d2::value_type>(in2_upper_left,in2_bottom_right);
2045  for(size_t i=0; i<(size_t)size2dout[0]; i++)
2046  for(size_t j=0; j<(size_t)size2dout[1]; j++)
2047  {
2048  box2.set_coord(i,j,size2din2[0]-1+i,size2din2[1]-1+j);
2049  typename InputIterator2d1::value_type mean1 = slip::mean<typename InputIterator2d1::value_type>(I.upper_left(box2),I.bottom_right(box2));
2050  d.set_coord(i,j);
2051  *(out_upper_left + d) = pseudo_normalized_crosscorrelation<Real>(in2_upper_left,in2_bottom_right,I.upper_left(box2),mean2,mean1);
2052  }
2053  }
2054  else if(t==CCC)
2055  {
2056  slip::Box2d<int> box2;
2058  typename InputIterator2d2::value_type mean2 = slip::mean<typename InputIterator2d2::value_type>(in2_upper_left,in2_bottom_right);
2059  for(size_t i=0; i<(size_t)size2dout[0]; i++)
2060  for(size_t j=0; j<(size_t)size2dout[1]; j++)
2061  {
2062  box2.set_coord(i,j,size2din2[0]-1+i,size2din2[1]-1+j);
2063  typename InputIterator2d1::value_type mean1 = slip::mean<typename InputIterator2d1::value_type>(I.upper_left(box2),I.bottom_right(box2));
2064  d.set_coord(i,j);
2065  *(out_upper_left + d) = ccc<Real>(in2_upper_left,in2_bottom_right,I.upper_left(box2),mean2,mean1);
2066  }
2067  }
2068  else if(t==CC)
2069  {
2070  slip::Box2d<int> box2;
2072  for(size_t i=0; i<(size_t)size2dout[0]; i++)
2073  for(size_t j=0; j<(size_t)size2dout[1]; j++)
2074  {
2075  box2.set_coord(i,j,size2din2[0]-1+i,size2din2[1]-1+j);
2076  d.set_coord(i,j);
2077  *(out_upper_left + d) = cc<Real>(in2_upper_left,in2_bottom_right,I.upper_left(box2));
2078  }
2079  }
2080  else
2081  {
2082  std::cerr<< slip::CROSSCORRELATION_BAD_ALGO << std::endl;
2083  slip::Box2d<int> box2;
2085  for(size_t i=0; i<(size_t)size2dout[0]; i++)
2086  for(size_t j=0; j<(size_t)size2dout[1]; j++)
2087  {
2088  box2.set_coord(i,j,size2din2[0]-1+i,size2din2[1]-1+j);
2089  d.set_coord(i,j);
2090  *(out_upper_left + d) = cc<Real>(in2_upper_left,in2_bottom_right,I.upper_left(box2));
2091  }
2092  }
2093  }
2094 
2115  template<typename Real, typename InputIterator2d1, typename InputIterator2d2, typename OutputIterator2d>
2116  inline
2117  void std_crosscorrelation2d(InputIterator2d1 in1_upper_left, InputIterator2d1 in1_bottom_right,
2118  InputIterator2d2 in2_upper_left, InputIterator2d2 in2_bottom_right,
2119  OutputIterator2d out_upper_left,
2120  OutputIterator2d out_bottom_right)
2121  {
2122  slip::crosscorrelation2d<Real>(in1_upper_left,in1_bottom_right,in2_upper_left,in2_bottom_right,out_upper_left,out_bottom_right,CC);
2123  }
2124 
2145  template<typename Real, typename InputIterator2d1, typename InputIterator2d2, typename OutputIterator2d>
2146  inline
2147  void centered_crosscorrelation2d(InputIterator2d1 in1_upper_left, InputIterator2d1 in1_bottom_right,
2148  InputIterator2d2 in2_upper_left, InputIterator2d2 in2_bottom_right,
2149  OutputIterator2d out_upper_left,
2150  OutputIterator2d out_bottom_right)
2151  {
2152  slip::crosscorrelation2d<Real>(in1_upper_left,in1_bottom_right,in2_upper_left,in2_bottom_right,out_upper_left,out_bottom_right,CCC);
2153  }
2154 
2155 
2176  template<typename Real, typename InputIterator2d1, typename InputIterator2d2, typename OutputIterator2d>
2177  inline
2178  void pseudo_normalized_crosscorrelation2d(InputIterator2d1 in1_upper_left, InputIterator2d1 in1_bottom_right,
2179  InputIterator2d2 in2_upper_left, InputIterator2d2 in2_bottom_right,
2180  OutputIterator2d out_upper_left,
2181  OutputIterator2d out_bottom_right)
2182  {
2183  slip::crosscorrelation2d<Real>(in1_upper_left,in1_bottom_right,in2_upper_left,in2_bottom_right,out_upper_left,out_bottom_right,PNCC);
2184  }
2185 
2205  template<typename Real, typename InputIterator2d1, typename InputIterator2d2, typename OutputIterator2d>
2206  inline
2207  void normalized_crosscorrelation2d(InputIterator2d1 in1_upper_left, InputIterator2d1 in1_bottom_right,
2208  InputIterator2d2 in2_upper_left, InputIterator2d2 in2_bottom_right,
2209  OutputIterator2d out_upper_left,
2210  OutputIterator2d out_bottom_right)
2211  {
2212  slip::crosscorrelation2d<Real>(in1_upper_left,in1_bottom_right,in2_upper_left,in2_bottom_right,out_upper_left,out_bottom_right,NCC);
2213  }
2214 /* @} */
2215 
2217  /* @{ */
2218 
2246 template<typename InputIterator1,
2247  typename InputIterator2,
2248  typename OutputIterator>
2249  inline
2250  void fft_crosscorrelation(InputIterator1 first, InputIterator1 last,
2251  InputIterator2 first2,InputIterator2 last2,
2252  OutputIterator result_first,
2253  OutputIterator result_last)
2254  {
2255  assert((last - first) != 0);
2256  assert((last2 - first2) != 0);
2257  assert((result_last-result_first) >= (last - first) + (last2 - first2) -1);
2258 
2259  std::size_t size= static_cast<std::size_t>((last-first) +(last2-first2)-1);
2260 
2261  typedef typename std::iterator_traits<InputIterator1>::value_type value_type;
2262  typedef typename slip::lin_alg_traits<value_type>::value_type Real;
2263 
2264 
2265  slip::Array<value_type> signal(size,first,last);
2266  slip::Array<value_type> kernel(size,first2,last2);
2267 
2268  slip::Array<std::complex<Real> > fft1(size,std::complex<Real>());
2269  slip::Array<std::complex<Real> > fft2(size,std::complex<Real>());
2270 
2271  slip::real_fft(signal.begin(), signal.end(), fft1.begin());
2272  slip::real_fft(kernel.begin(),kernel.end(), fft2.begin());
2273 
2274  for (std::size_t i = 0; i < size; ++i)
2275  {
2276  fft1[i] *= std::conj(fft2[i]);
2277  }
2278  slip::ifft(fft1.begin(), fft1.end(), fft2.begin());
2279  slip::fftshift(fft2.begin(),fft2.end());
2280  std::transform(fft2.begin(),fft2.end(),result_first,
2282  }
2283 
2318 template<typename InputIterator1,
2319  typename InputIterator2,
2320  typename OutputIterator>
2321  inline
2322  void fft_crosscorrelation_same(InputIterator1 first, InputIterator1 last,
2323  InputIterator2 first2,InputIterator2 last2,
2324  OutputIterator result_first,
2325  OutputIterator result_last)
2326  {
2327  assert((last - first) == (last2 - first2));
2328  assert((last - first) == (result_last - result_first));
2329  assert((last - first) == (last2 - first2));
2330 
2331  std::size_t size= static_cast<std::size_t>((last-first) +(last2-first2)-1);
2332 
2333  typedef typename std::iterator_traits<InputIterator1>::value_type value_type;
2334  typedef typename slip::lin_alg_traits<value_type>::value_type Real;
2335 
2336 
2337  slip::Array<value_type> signal(size,first,last);
2338  slip::Array<value_type> kernel(size,first2,last2);
2339 
2340  slip::Array<std::complex<Real> > fft1(size,std::complex<Real>());
2341  slip::Array<std::complex<Real> > fft2(size,std::complex<Real>());
2342 
2343  slip::real_fft(signal.begin(), signal.end(), fft1.begin());
2344  slip::real_fft(kernel.begin(),kernel.end(), fft2.begin());
2345 
2346  for (std::size_t i = 0; i < size; ++i)
2347  {
2348  fft1[i] *= std::conj(fft2[i]);
2349  }
2350  slip::ifft(fft1.begin(), fft1.end(), fft2.begin());
2351  slip::fftshift(fft2.begin(),fft2.end());
2352  std::transform(fft2.begin()+((last2-first2)/2),
2353  fft2.begin()+((last2-first2)/2) + (last-first),result_first,
2354  slip::un_real<std::complex<Real>,Real>());
2355  }
2356 
2357 
2358 
2393  template<typename RandomAccessIterator1,
2394  typename RandomAccessIterator2,
2395  typename RandomAccessIterator3>
2396  inline
2397  void fft_circular_crosscorrelation(RandomAccessIterator1 first , RandomAccessIterator1 last,
2398  RandomAccessIterator2 first2, RandomAccessIterator2 last2,
2399  RandomAccessIterator3 result_first, RandomAccessIterator3 result_last )
2400  {
2401  // assert((last - first) == (last2 - first2));
2402  // assert((last - first) == (result_last - result_first));
2403  // assert((last - first) == (last2 - first2));
2404 
2405  // typedef typename std::iterator_traits<RandomAccessIterator1>::value_type value_type;
2406  // typedef typename slip::lin_alg_traits<value_type>::value_type T;
2407 
2408 
2409  // std::size_t N = static_cast<std::size_t>(last-first);
2410  // slip::Array<std::complex<T> > fft1(N);
2411  // slip::Array<std::complex<T> > fft2(N);
2412  // slip::real_fft(first,last,fft1.begin());
2413  // slip::real_fft(first2,last2,fft2.begin());
2414  // typename slip::Array2d<std::complex<T> >::iterator it_fft1 = fft1.begin();
2415  // typename slip::Array2d<std::complex<T> >::iterator it_fft1_last = fft1.end();
2416 
2417  // typename slip::Array2d<std::complex<T> >::iterator it_fft2 = fft2.begin();
2418  // for(; it_fft1 != it_fft1_last; ++it_fft1, ++it_fft2)
2419  // {
2420  // *it_fft1 *= std::conj(*it_fft2);
2421  // }
2422  // slip::ifft(fft1.begin(),fft1.end(),fft2.begin());
2423  // slip::fftshift(fft2.begin(),fft2.end());
2424  // std::transform(fft2.begin(),fft2.end(),result_first,
2425  // slip::un_real<std::complex<T>,T>());
2426 
2427 
2428  slip::fft_circular_convolution(first,last,
2429  std::reverse_iterator<RandomAccessIterator2>(last2),
2430  std::reverse_iterator<RandomAccessIterator2>(first2),
2431  result_first,result_last);
2432  }
2469 template<typename InputIterator1,
2470  typename InputIterator2,
2471  typename OutputIterator>
2472  inline
2473  void complex_fft_crosscorrelation(InputIterator1 first, InputIterator1 last,
2474  InputIterator2 first2,InputIterator2 last2,
2475  OutputIterator result_first,
2476  OutputIterator result_last)
2477  {
2478  assert((last - first) != 0);
2479  assert((last2 - first2) != 0);
2480  assert((result_last-result_first) >= (last - first) + (last2 - first2) -1);
2481 
2482 
2483  std::size_t size= static_cast<std::size_t>((last-first) +(last2-first2)-1);
2484 
2485  typedef typename std::iterator_traits<InputIterator1>::value_type value_type;
2486  typedef typename slip::lin_alg_traits<value_type>::value_type Real;
2487 
2488 
2489  slip::Array<value_type> signal(size,first,last);
2490  slip::Array<value_type> kernel(size,first2,last2);
2491 
2492  slip::Array<std::complex<Real> > fft1(size,std::complex<Real>());
2493  slip::Array<std::complex<Real> > fft2(size,std::complex<Real>());
2494 
2495  slip::complex_fft(signal.begin(), signal.end(), fft1.begin());
2496  slip::complex_fft(kernel.begin(),kernel.end(), fft2.begin());
2497 
2498  for (std::size_t i = 0; i < size; ++i)
2499  {
2500  fft1[i] *= std::conj(fft2[i]);
2501  }
2502  slip::ifft(fft1.begin(), fft1.end(), result_first);
2503  slip::fftshift(result_first,result_first+size);
2504  }
2505 
2546  template<typename InputIterator2d1, typename InputIterator2d2, typename OutputIterator2d>
2547  inline
2548  void fft_circular_crosscorrelation2d(InputIterator2d1 in1_upper_left,
2549  InputIterator2d1 in1_bottom_right,
2550  InputIterator2d2 in2_upper_left,
2551  InputIterator2d2 in2_bottom_right,
2552  OutputIterator2d out_upper_left,
2553  OutputIterator2d out_bottom_right)
2554  {
2555 
2556  assert((in1_bottom_right - in1_upper_left) == (out_bottom_right - out_upper_left));
2557 
2558  assert((in2_bottom_right - in2_upper_left) == (out_bottom_right - out_upper_left));
2559  typename InputIterator2d1::difference_type size2din1 = in1_bottom_right - in1_upper_left;
2560  typename InputIterator2d2::difference_type size2din2 = in2_bottom_right - in2_upper_left;
2561  typename OutputIterator2d::difference_type size2dout = out_bottom_right - out_upper_left;
2563 
2564  slip::Array2d<std::complex<Real> > fft1(size2din1[0],size2din1[1]);
2565  slip::Array2d<std::complex<Real> > fft2(size2din1[0],size2din1[1]);
2566 
2567  slip::Array2d<Real> Kernel(size2din1[0],size2din1[1]);
2568  typedef std::reverse_iterator<InputIterator2d2> reverse_iterator;
2569  std::copy(reverse_iterator(in2_bottom_right- slip::DPoint2d<int>(1,0)),
2570  reverse_iterator(in2_upper_left),
2571  Kernel.upper_left());
2572 
2573  slip::real_fft2d(in1_upper_left, in1_bottom_right, fft1.upper_left());
2574  slip::real_fft2d(Kernel.upper_left(),Kernel.bottom_right(),
2575  fft2.upper_left());
2576  slip::multiplies(fft1.begin(),fft1.end(),fft2.begin(),fft1.begin());
2577 
2578 
2579 
2580  slip::ifft2d(fft1.upper_left(), fft1.bottom_right(), fft2.upper_left());
2581  slip::fftshift2d(fft2.upper_left(),fft2.bottom_right());
2582  std::transform(fft2.upper_left(),fft2.bottom_right(),out_upper_left,slip::un_real<std::complex<Real>,Real>());
2583 
2584  }
2585 
2625  template<typename InputIterator2d1, typename InputIterator2d2, typename OutputIterator2d>
2626  inline
2627  void fft_crosscorrelation2d(InputIterator2d1 in1_upper_left,
2628  InputIterator2d1 in1_bottom_right,
2629  InputIterator2d2 in2_upper_left,
2630  InputIterator2d2 in2_bottom_right,
2631  OutputIterator2d out_upper_left,
2632  OutputIterator2d out_bottom_right)
2633  {
2634  assert(((in1_bottom_right - in1_upper_left)[0]+(in2_bottom_right - in2_upper_left)[0] - 1) == (out_bottom_right - out_upper_left)[0]);
2635  assert(((in1_bottom_right - in1_upper_left)[1]+(in2_bottom_right - in2_upper_left)[1] - 1) == (out_bottom_right - out_upper_left)[1]);
2636 
2637  typename InputIterator2d1::difference_type size2din1 = in1_bottom_right - in1_upper_left;
2638  typename InputIterator2d2::difference_type size2din2 = in2_bottom_right - in2_upper_left;
2639  typename OutputIterator2d::difference_type size2dout = out_bottom_right - out_upper_left;
2640 
2641  typedef typename std::iterator_traits<InputIterator2d1>::value_type value_type;
2642  typedef typename slip::lin_alg_traits<value_type>::value_type Real;
2643 
2644  std::size_t rows = size2din1[0] + size2din2[0] - 1;
2645  std::size_t cols = size2din1[1] + size2din2[1] - 1;
2646  //zero-padd Signal and Kernel data
2647  slip::Array2d<Real> Signal_zero_padded(rows,cols,Real());
2648  slip::Array2d<Real> Kernel_zero_padded(rows,cols,Real());
2649  slip::Box2d<int> signal_box(0,0,
2650  static_cast<int>(size2din1[0]-1),
2651  static_cast<int>(size2din1[1]-1));
2652  slip::Box2d<int> kernel_box(0,0,
2653  static_cast<int>(size2din2[0]-1),
2654  static_cast<int>(size2din2[1]-1));
2655  std::copy(in1_upper_left,in1_bottom_right,Signal_zero_padded.upper_left(signal_box));
2656  //rotate the kernel by 180°
2657  typedef std::reverse_iterator<InputIterator2d2> reverse_iterator;
2658  std::copy(reverse_iterator(in2_bottom_right- slip::DPoint2d<int>(1,0)),
2659  reverse_iterator(in2_upper_left),
2660  Kernel_zero_padded.upper_left(kernel_box));
2661 
2662 
2663  slip::Array2d<std::complex<Real> > fft1(rows,cols,std::complex<Real>());
2664  slip::Array2d<std::complex<Real> > fft2(rows,cols,std::complex<Real>());
2665 
2666  //computes the two fft
2667  slip::real_fft2d(Signal_zero_padded,fft1);
2668  slip::real_fft2d(Kernel_zero_padded,fft2);
2669  //computes fft1*fft2
2670  slip::multiplies(fft1.begin(),fft1.end(),fft2.begin(),fft1.begin());
2671  //ifft2d
2672  slip::ifft2d(fft1.upper_left(), fft1.bottom_right(), fft2.upper_left());
2673 
2674  //the result should be real
2675  std::transform(fft2.begin(),fft2.end(),out_upper_left,slip::un_real<std::complex<Real>,Real>());
2676 
2677 
2678  }
2679 
2680 
2681 
2714 template<typename InputIterator2d1,
2715  typename InputIterator2d2,
2716  typename OutputIterator2d>
2717 inline
2718 void fft_crosscorrelation2d_same(InputIterator2d1 in1_upper_left,
2719  InputIterator2d1 in1_bottom_right,
2720  InputIterator2d2 in2_upper_left,
2721  InputIterator2d2 in2_bottom_right,
2722  OutputIterator2d out_upper_left,
2723  OutputIterator2d out_bottom_right)
2724  {
2725  assert((in1_bottom_right - in1_upper_left)[0] == (out_bottom_right - out_upper_left)[0]);
2726  assert((in1_bottom_right - in1_upper_left)[1] == (out_bottom_right - out_upper_left)[1]);
2727 
2728  typename InputIterator2d1::difference_type size2din1 = in1_bottom_right - in1_upper_left;
2729  typename InputIterator2d2::difference_type size2din2 = in2_bottom_right - in2_upper_left;
2730  typename OutputIterator2d::difference_type size2dout = out_bottom_right - out_upper_left;
2731 
2732  typedef typename std::iterator_traits<InputIterator2d1>::value_type value_type;
2733  typedef typename slip::lin_alg_traits<value_type>::value_type Real;
2734 
2735  std::size_t rows = size2din1[0] + size2din2[0] - 1;
2736  std::size_t cols = size2din1[1] + size2din2[1] - 1;
2737  //zero-padd Signal and Kernel data
2738  slip::Array2d<Real> Signal_zero_padded(rows,cols,Real());
2739  slip::Array2d<Real> Kernel_zero_padded(rows,cols,Real());
2740  slip::Box2d<int> signal_box(0,0,
2741  static_cast<int>(size2din1[0]-1),
2742  static_cast<int>(size2din1[1]-1));
2743  slip::Box2d<int> kernel_box(0,0,
2744  static_cast<int>(size2din2[0]-1),
2745  static_cast<int>(size2din2[1]-1));
2746  std::copy(in1_upper_left,in1_bottom_right,Signal_zero_padded.upper_left(signal_box));
2747  //rotate the kernel by 180°
2748  typedef std::reverse_iterator<InputIterator2d2> reverse_iterator;
2749  std::copy(reverse_iterator(in2_bottom_right- slip::DPoint2d<int>(1,0)),
2750  reverse_iterator(in2_upper_left),
2751  Kernel_zero_padded.upper_left(kernel_box));
2752 
2753  slip::Array2d<std::complex<Real> > fft1(rows,cols,std::complex<Real>());
2754  slip::Array2d<std::complex<Real> > fft2(rows,cols,std::complex<Real>());
2755 
2756  //computes the two fft
2757  slip::real_fft2d(Signal_zero_padded,fft1);
2758  slip::real_fft2d(Kernel_zero_padded,fft2);
2759  //computes fft1*fft2
2760  slip::multiplies(fft1.begin(),fft1.end(),fft2.begin(),fft1.begin());
2761  //ifft2d
2762  slip::ifft2d(fft1.upper_left(), fft1.bottom_right(), fft2.upper_left());
2763 
2764  //the result should be real
2765 
2766  slip::Box2d<int> out_box(size2din2[0]/2,size2din2[1]/2,
2767  static_cast<int>(size2din2[0]/2 + size2din1[0] - 1),
2768  static_cast<int>(size2din2[1]/2 + size2din1[1] - 1));
2769 
2770  std::transform(fft2.upper_left(out_box),
2771  fft2.bottom_right(out_box),
2772  out_upper_left,slip::un_real<std::complex<Real>,Real>());
2773 
2774  }
2775 
2776 
2809  template<typename InputIterator3d1,
2810  typename InputIterator3d2,
2811  typename OutputIterator3d>
2812  inline
2813  void fft_circular_crosscorrelation3d(InputIterator3d1 in1_front_upper_left,
2814  InputIterator3d1 in1_back_bottom_right,
2815  InputIterator3d2 in2_front_upper_left,
2816  InputIterator3d2 in2_back_bottom_right,
2817  OutputIterator3d out_front_upper_left,
2818  OutputIterator3d out_back_bottom_right)
2819  {
2820  assert((in1_back_bottom_right - in1_front_upper_left) == (out_back_bottom_right - out_front_upper_left));
2821  assert((in2_back_bottom_right - in2_front_upper_left) == (out_back_bottom_right - out_front_upper_left));
2822 
2823  typename InputIterator3d1::difference_type size3din1 = in1_back_bottom_right - in1_front_upper_left;
2824  typename InputIterator3d2::difference_type size3din2 = in2_back_bottom_right - in2_front_upper_left;
2825  typename OutputIterator3d::difference_type size3dout = out_back_bottom_right - out_front_upper_left;
2827 
2828  slip::Array3d<std::complex<Real> > fft1(size3din1[0],size3din1[1],size3din1[2]);
2829  slip::Array3d<std::complex<Real> > fft2(size3din1[0],size3din1[1],size3din1[2]);
2830 
2831  slip::real_fft3d(in1_front_upper_left, in1_back_bottom_right, fft1.front_upper_left());
2832  slip::real_fft3d(in2_front_upper_left, in2_back_bottom_right, fft2.front_upper_left());
2833 
2834  for (std::size_t k=0;k<(size_t)size3din1[0];k++)
2835  for (std::size_t i=0;i<(size_t)size3din1[1];i++)
2836  for (std::size_t j=0;j<(size_t)size3din1[2];j++)
2837  {
2838  fft1[k][i][j] *= std::conj(fft2[k][i][j]);
2839  }
2842  std::transform(fft2.front_upper_left(),fft2.back_bottom_right(),out_front_upper_left,slip::un_real<std::complex<Real>,Real>());
2843  }
2844 
2845 
2883  template<typename InputIterator3d1, typename InputIterator3d2, typename OutputIterator3d>
2884  inline
2885  void
2886  fft_crosscorrelation3d(InputIterator3d1 in1_front_upper_left,
2887  InputIterator3d1 in1_back_bottom_right,
2888  InputIterator3d2 in2_front_upper_left,
2889  InputIterator3d2 in2_back_bottom_right,
2890  OutputIterator3d out_front_upper_left,
2891  OutputIterator3d out_back_bottom_right)
2892  {
2893  assert(((in1_back_bottom_right - in1_front_upper_left)[0]+(in2_back_bottom_right - in2_front_upper_left)[0] - 1) == (out_back_bottom_right - out_front_upper_left)[0]);
2894  assert(((in1_back_bottom_right - in1_front_upper_left)[1]+(in2_back_bottom_right - in2_front_upper_left)[1] - 1) == (out_back_bottom_right - out_front_upper_left)[1]);
2895 assert(((in1_back_bottom_right - in1_front_upper_left)[2]+(in2_back_bottom_right - in2_front_upper_left)[2] - 1) == (out_back_bottom_right - out_front_upper_left)[2]);
2896 
2897  typename InputIterator3d1::difference_type size3din1 = in1_back_bottom_right - in1_front_upper_left;
2898  typename InputIterator3d2::difference_type size3din2 = in2_back_bottom_right - in2_front_upper_left;
2899  typename OutputIterator3d::difference_type size3dout = out_back_bottom_right - out_front_upper_left;
2900 
2902 
2903  std::size_t slices = size3din1[0] + size3din2[0] - 1;
2904  std::size_t rows = size3din1[1] + size3din2[1] - 1;
2905  std::size_t cols = size3din1[2] + size3din2[2] - 1;
2906 
2907  //zero-padd Signal and Kernel data
2908  slip::Array3d<Real> Signal_zero_padded(slices,rows,cols,Real());
2909  slip::Array3d<Real> Kernel_zero_padded(slices,rows,cols,Real());
2910  slip::Box3d<int> signal_box(0,0,0,
2911  static_cast<int>(size3din1[0]-1),
2912  static_cast<int>(size3din1[1]-1),
2913  static_cast<int>(size3din1[2]-1));
2914  slip::Box3d<int> kernel_box(0,0,0,
2915  static_cast<int>(size3din2[0]-1),
2916  static_cast<int>(size3din2[1]-1),
2917  static_cast<int>(size3din2[2]-1));
2918  std::copy(in1_front_upper_left,
2919  in1_back_bottom_right,
2920  Signal_zero_padded.front_upper_left(signal_box));
2921  //rotate the kernel
2922  typedef std::reverse_iterator<InputIterator3d2> reverse_iterator;
2923  std::copy(reverse_iterator(in2_back_bottom_right- slip::DPoint3d<int>(1,1,0)),
2924  reverse_iterator(in2_front_upper_left),
2925  Kernel_zero_padded.front_upper_left(kernel_box));
2926 
2927 
2928  slip::Array3d<std::complex<Real> > fft1(slices,rows,cols,std::complex<Real>());
2929  slip::Array3d<std::complex<Real> > fft2(slices,rows,cols,std::complex<Real>());
2930  //computes the two fft
2931  slip::real_fft3d(Signal_zero_padded.front_upper_left(),
2932  Signal_zero_padded.back_bottom_right(),
2933  fft1.front_upper_left());
2934  slip::real_fft3d(Kernel_zero_padded.front_upper_left(),
2935  Kernel_zero_padded.back_bottom_right(),
2936  fft2.front_upper_left());
2937  //computes fft1*fft2
2938  slip::multiplies(fft1.begin(),fft1.end(),fft2.begin(),fft1.begin());
2939 
2940  //ifft3d
2941  slip::ifft3d(fft1.front_upper_left(), fft1.back_bottom_right(),
2942  fft2.front_upper_left());
2943  //fftshift3d
2945  //the result should be real
2946  std::transform(fft2.front_upper_left(),fft2.back_bottom_right(),
2947  out_front_upper_left,
2949  }
2950 
2951 
2988  template<typename InputIterator3d1,
2989  typename InputIterator3d2,
2990  typename OutputIterator3d>
2991  inline
2992  void
2993  fft_crosscorrelation3d_same(InputIterator3d1 in1_front_upper_left,
2994  InputIterator3d1 in1_back_bottom_right,
2995  InputIterator3d2 in2_front_upper_left,
2996  InputIterator3d2 in2_back_bottom_right,
2997  OutputIterator3d out_front_upper_left,
2998  OutputIterator3d out_back_bottom_right)
2999  {
3000  assert( (in1_back_bottom_right - in1_front_upper_left)
3001  == (out_back_bottom_right - out_front_upper_left));
3002 
3003  typename InputIterator3d1::difference_type size3din1 = in1_back_bottom_right - in1_front_upper_left;
3004  typename InputIterator3d2::difference_type size3din2 = in2_back_bottom_right - in2_front_upper_left;
3005  typename OutputIterator3d::difference_type size3dout = out_back_bottom_right - out_front_upper_left;
3006 
3008 
3009  std::size_t slices = size3din1[0] + size3din2[0] - 1;
3010  std::size_t rows = size3din1[1] + size3din2[1] - 1;
3011  std::size_t cols = size3din1[2] + size3din2[2] - 1;
3012 
3013  //zero-padd Signal and Kernel data
3014  slip::Array3d<Real> Signal_zero_padded(slices,rows,cols,Real());
3015  slip::Array3d<Real> Kernel_zero_padded(slices,rows,cols,Real());
3016  slip::Box3d<int> signal_box(0,0,0,
3017  static_cast<int>(size3din1[0]-1),
3018  static_cast<int>(size3din1[1]-1),
3019  static_cast<int>(size3din1[2]-1));
3020  slip::Box3d<int> kernel_box(0,0,0,
3021  static_cast<int>(size3din2[0]-1),
3022  static_cast<int>(size3din2[1]-1),
3023  static_cast<int>(size3din2[2]-1));
3024  std::copy(in1_front_upper_left,
3025  in1_back_bottom_right,
3026  Signal_zero_padded.front_upper_left(signal_box));
3027  //rotate the kernel
3028  typedef std::reverse_iterator<InputIterator3d2> reverse_iterator;
3029  std::copy(reverse_iterator(in2_back_bottom_right- slip::DPoint3d<int>(1,1,0)),
3030  reverse_iterator(in2_front_upper_left),
3031  Kernel_zero_padded.front_upper_left(kernel_box));
3032 
3033 
3034  slip::Array3d<std::complex<Real> > fft1(slices,rows,cols,std::complex<Real>());
3035  slip::Array3d<std::complex<Real> > fft2(slices,rows,cols,std::complex<Real>());
3036  //computes the two fft
3037  slip::real_fft3d(Signal_zero_padded.front_upper_left(),
3038  Signal_zero_padded.back_bottom_right(),
3039  fft1.front_upper_left());
3040  slip::real_fft3d(Kernel_zero_padded.front_upper_left(),
3041  Kernel_zero_padded.back_bottom_right(),
3042  fft2.front_upper_left());
3043 
3044  //computes fft1*fft2
3045  slip::multiplies(fft1.begin(),fft1.end(),fft2.begin(),fft1.begin());
3046 
3047  //ifft3d
3048  slip::ifft3d(fft1.front_upper_left(), fft1.back_bottom_right(),
3049  fft2.front_upper_left());
3050  //fftshift3d
3052  //the result should be real
3053  slip::Box3d<int> out_box(size3din2[0]/2,size3din2[1]/2,size3din2[2]/2,
3054  static_cast<int>(size3din2[0]/2 + size3din1[0] - 1),
3055  static_cast<int>(size3din2[1]/2 + size3din1[1] - 1),
3056  static_cast<int>(size3din2[2]/2 + size3din1[2] - 1));
3057 
3058  std::transform(fft2.front_upper_left(out_box),
3059  fft2.back_bottom_right(out_box),
3060  out_front_upper_left,slip::un_real<std::complex<Real>,Real>());
3061 
3062  }
3063 
3064  /* @} */
3065 
3066 }//slip::
3067 
3068 #endif //SLIP_CORRELATION_HPP
void pseudo_normalized_crosscorrelation2d(InputIterator2d1 in1_upper_left, InputIterator2d1 in1_bottom_right, InputIterator2d2 in2_upper_left, InputIterator2d2 in2_bottom_right, OutputIterator2d out_upper_left, OutputIterator2d out_bottom_right)
Computes the pseudo normalized crosscorrelation between two Images.
T cc_mask(InputIterator1 first1, InputIterator1 last1, MaskIterator mask_first, InputIterator2 first2, typename std::iterator_traits< MaskIterator >::value_type value=typename std::iterator_traits< MaskIterator >::value_type(1))
Computes the standard crosscorrelation between two sequences according to a mask sequence ...
void set_coord(const CoordType &x11, const CoordType &x12, const CoordType &x21, const CoordType &x22)
Mutator of the coordinates of the Box2d.
Definition: Box2d.hpp:351
void fft_circular_crosscorrelation2d(InputIterator2d1 in1_upper_left, InputIterator2d1 in1_bottom_right, InputIterator2d2 in2_upper_left, InputIterator2d2 in2_bottom_right, OutputIterator2d out_upper_left, OutputIterator2d out_bottom_right)
Computes the circular crosscorrelation between two 2D sequences using fft2d.
iterator2d upper_left()
Returns a read/write iterator2d that points to the first element of the Array2d. It points to the upp...
Definition: Array2d.hpp:2744
void fft_circular_crosscorrelation3d(InputIterator3d1 in1_front_upper_left, InputIterator3d1 in1_back_bottom_right, InputIterator3d2 in2_front_upper_left, InputIterator3d2 in2_back_bottom_right, OutputIterator3d out_front_upper_left, OutputIterator3d out_back_bottom_right)
Computes the standard crosscorrelation between two 3D sequences using fft3d.
void complex_fft_crosscorrelation(InputIterator1 first, InputIterator1 last, InputIterator2 first2, InputIterator2 last2, OutputIterator result_first, OutputIterator result_last)
Computes the standard complex crosscorrelation between two sequences using fft.
void fftshift(ForwardIterator first, ForwardIterator last)
Performs a shift of a container, for use with the fft and ifft functions, in order to move the freque...
Definition: FFT.hpp:166
void fftshift2d(ForwardIterator2D upper_left, ForwardIterator2D bottom_right)
Performs a shift of a 2d container, for use with the fft2d and ifft2d functions, in order to move the...
Definition: FFT.hpp:191
iterator3d front_upper_left()
Returns a read/write iterator3d that points to the first element of the Array3d. It points to the fro...
Definition: Array3d.hpp:4882
T cc_if(InputIterator1 first1, InputIterator1 last1, InputIterator2 first2, Predicate pred)
Computes the standard crosscorrelation between two sequences according to a Predicate. .
void fft_circular_convolution(RandomAccessIterator1 first, RandomAccessIterator1 last, RandomAccessIterator2 first2, RandomAccessIterator2 last2, RandomAccessIterator3 conv_first, RandomAccessIterator3 conv_last)
Computes the FFT circular convolution.
Real pseudo_normalized_crosscorrelation(InputIterator1 first, InputIterator1 last, InputIterator2 first2, typename std::iterator_traits< InputIterator1 >::value_type mean1=typename std::iterator_traits< InputIterator1 >::value_type(), typename std::iterator_traits< InputIterator1 >::value_type mean2=typename std::iterator_traits< InputIterator2 >::value_type())
Computes the pseudo normalized crosscorrelation between two sequences: .
void fft_circular_crosscorrelation(RandomAccessIterator1 first, RandomAccessIterator1 last, RandomAccessIterator2 first2, RandomAccessIterator2 last2, RandomAccessIterator3 result_first, RandomAccessIterator3 result_last)
Computes the FFT circular crosscorrelation.
Real ncc_if(InputIterator1 first1, InputIterator1 last1, InputIterator2 first2, Predicate pred, typename std::iterator_traits< InputIterator1 >::value_type mean1=typename std::iterator_traits< InputIterator1 >::value_type(), typename std::iterator_traits< InputIterator1 >::value_type mean2=typename std::iterator_traits< InputIterator2 >::value_type())
Computes the standard normalized crosscorrelation between two sequences according to a Predicate: ...
Difference of Point3D class, specialization of DPoint<CoordType,DIM> with DIM = 3.
T std_crosscorrelation_if(InputIterator1 first1, InputIterator1 last1, InputIterator2 first2, Predicate pred)
Computes the standard crosscorrelation between two sequences according to a Predicate. .
void autocorrelation_full(RandomAccessIterator1 first, RandomAccessIterator1 last, RandomAccessIterator2 auto_first, RandomAccessIterator2 auto_last)
Computes the autocorrelations from lag -(last-first) to (last-first) of a range [first,last).
void centered_crosscorrelation2d(InputIterator2d1 in1_upper_left, InputIterator2d1 in1_bottom_right, InputIterator2d2 in2_upper_left, InputIterator2d2 in2_bottom_right, OutputIterator2d out_upper_left, OutputIterator2d out_bottom_right)
Computes the centered crosscorrelation between two Images.
This is a three-dimensional dynamic and generic container. This container statisfies the Bidirectionn...
Definition: Array3d.hpp:109
Real centered_crosscorrelation(InputIterator1 first, InputIterator1 last, InputIterator2 first2, typename std::iterator_traits< InputIterator1 >::value_type mean1=typename std::iterator_traits< InputIterator1 >::value_type(), typename std::iterator_traits< InputIterator1 >::value_type mean2=typename std::iterator_traits< InputIterator2 >::value_type())
Computes the centered crosscorrelation between two sequences: .
const_iterator begin() const
Returns a read-only (constant) iterator that points to the first element in the Array3d. Iteration is done in ordinary element order.
Definition: Array3d.hpp:3644
void biased_autocorrelation(RandomAccessIterator1 first, RandomAccessIterator1 last, RandomAccessIterator2 auto_first, RandomAccessIterator2 auto_last)
Computes the biased autocorrelations from lag -(last-first) to (last-first) of a range [first...
void fft_crosscorrelation2d(InputIterator2d1 in1_upper_left, InputIterator2d1 in1_bottom_right, InputIterator2d2 in2_upper_left, InputIterator2d2 in2_bottom_right, OutputIterator2d out_upper_left, OutputIterator2d out_bottom_right)
Computes the standard crosscorrelation between two 2D sequences using fft2d.
iterator begin()
Returns a read/write iterator that points to the first element in the Array. Iteration is done in ord...
Definition: Array.hpp:1077
iterator2d bottom_right()
Returns a read/write iterator2d that points to the past the end element of the Array2d. It points to past the end element of the bottom right element of the Array2d.
Definition: Array2d.hpp:2759
size_type cols() const
Returns the number of columns (second dimension size) in the Matrix.
Definition: Matrix.hpp:3512
T cc(InputIterator1 first, InputIterator1 last, InputIterator2 first2)
Computes the standard crosscorrelation between two sequences: . Multiplies successive elements from t...
Real functor. Return the real part of x.
Definition: macros.hpp:98
void biased_autocovariance(RandomAccessIterator1 first, RandomAccessIterator1 last, RandomAccessIterator2 auto_first, RandomAccessIterator2 auto_last)
Computes the biased autocovariance from lag -(last-first) to 0 of a range [first,last).
Real ccc_mask(InputIterator1 first1, InputIterator1 last1, MaskIterator mask_first, InputIterator2 first2, typename std::iterator_traits< InputIterator1 >::value_type mean1=typename std::iterator_traits< InputIterator1 >::value_type(), typename std::iterator_traits< InputIterator1 >::value_type mean2=typename std::iterator_traits< InputIterator2 >::value_type(), typename std::iterator_traits< MaskIterator >::value_type value=typename std::iterator_traits< MaskIterator >::value_type(1))
Computes the centered crosscorrelation between two sequences according to a mask sequence ...
Real pseudo_normalized_crosscorrelation_if(InputIterator1 first1, InputIterator1 last1, InputIterator2 first2, Predicate pred, typename std::iterator_traits< InputIterator1 >::value_type mean1=typename std::iterator_traits< InputIterator1 >::value_type(), typename std::iterator_traits< InputIterator1 >::value_type mean2=typename std::iterator_traits< InputIterator2 >::value_type())
Computes the standard pseudo normalized crosscorrelation between two sequences according to a Predica...
void autocorrelation(RandomAccessIterator1 first, RandomAccessIterator1 last, RandomAccessIterator2 auto_first, RandomAccessIterator2 auto_last)
Computes the autocorrelations from lag -(last-first) to (last-first) of a range [first,last).
void real_fft2d(InputBidirectionalIterator2d in_upper_left, InputBidirectionalIterator2d in_bottom_right, OutputBidirectionalIterator2d out_upper_left)
Computes the real fft2d of a container with 2d iterators.
Definition: FFT.hpp:2553
Real ncc(InputIterator1 first, InputIterator1 last, InputIterator2 first2, typename std::iterator_traits< InputIterator1 >::value_type mean1=typename std::iterator_traits< InputIterator1 >::value_type(), typename std::iterator_traits< InputIterator1 >::value_type mean2=typename std::iterator_traits< InputIterator2 >::value_type())
Computes the standard normalized crosscorrelation between two sequences: .
Provides common linear algebra algorithms.
Real ccc_if(InputIterator1 first1, InputIterator1 last1, InputIterator2 first2, Predicate pred, typename std::iterator_traits< InputIterator1 >::value_type mean1=typename std::iterator_traits< InputIterator1 >::value_type(), typename std::iterator_traits< InputIterator1 >::value_type mean2=typename std::iterator_traits< InputIterator2 >::value_type())
Computes the standard crosscorrelation between two sequences according to a Predicate: ...
void unbiased_autocovariance(RandomAccessIterator1 first, RandomAccessIterator1 last, RandomAccessIterator2 auto_first, RandomAccessIterator2 auto_last)
Computes the unbiased autocovariance from lag -(last-first) to 0 of a range [first,last).
T std_crosscorrelation_mask(InputIterator1 first1, InputIterator1 last1, MaskIterator mask_first, InputIterator2 first2, typename std::iterator_traits< MaskIterator >::value_type value=typename std::iterator_traits< MaskIterator >::value_type(1))
Computes the standard crosscorrelation between two sequences according to a mask sequence ...
void fft_crosscorrelation3d_same(InputIterator3d1 in1_front_upper_left, InputIterator3d1 in1_back_bottom_right, InputIterator3d2 in2_front_upper_left, InputIterator3d2 in2_back_bottom_right, OutputIterator3d out_front_upper_left, OutputIterator3d out_back_bottom_right)
Computes the standard crosscorrelation between two 3D sequences using fft3d.
Provides a class to manipulate 1d dynamic and generic arrays.
void ifft3d(InputBidirectionalIterator3d in_front_upper_left, InputBidirectionalIterator3d in_back_bottom_right, OutputBidirectionalIterator3d out_front_upper_left)
Computes the ifft3d of a container with 3d iterators.
Definition: FFT.hpp:2852
T std_crosscorrelation(InputIterator1 first, InputIterator1 last, InputIterator2 first2)
Computes the standard crosscorrelation between two sequences: . Multiplies successive elements from t...
void biased_autocovariance_matrix(RandomAccessIterator first, RandomAccessIterator last, Matrix &A)
Constructs the biased autocorrelation matrix from given a range.
Real normalized_crosscorrelation(InputIterator1 first, InputIterator1 last, InputIterator2 first2, typename std::iterator_traits< InputIterator1 >::value_type mean1=typename std::iterator_traits< InputIterator1 >::value_type(), typename std::iterator_traits< InputIterator1 >::value_type mean2=typename std::iterator_traits< InputIterator2 >::value_type())
Computes the standard normalized crosscorrelation between two sequences: .
void ifft(InputIter in_begin, InputIter in_end, OutputIter out_begin)
Computes the inverse fft of a container.
Definition: FFT.hpp:2507
void fft_crosscorrelation_same(InputIterator1 first, InputIterator1 last, InputIterator2 first2, InputIterator2 last2, OutputIterator result_first, OutputIterator result_last)
Computes the FFT central part of the crosscorrelation.
Real ccc(InputIterator1 first, InputIterator1 last, InputIterator2 first2, typename std::iterator_traits< InputIterator1 >::value_type mean1=typename std::iterator_traits< InputIterator1 >::value_type(), typename std::iterator_traits< InputIterator1 >::value_type mean2=typename std::iterator_traits< InputIterator2 >::value_type())
Computes the centered crosscorrelation between two sequences: .
void autocovariance(RandomAccessIterator1 first, RandomAccessIterator1 last, RandomAccessIterator2 auto_first, RandomAccessIterator2 auto_last)
Computes the autocovariance from lag -(last-first) to 0 of a range [first,last).
Real pncc_mask(InputIterator1 first1, InputIterator1 last1, MaskIterator mask_first, InputIterator2 first2, typename std::iterator_traits< InputIterator1 >::value_type mean1=typename std::iterator_traits< InputIterator1 >::value_type(), typename std::iterator_traits< InputIterator1 >::value_type mean2=typename std::iterator_traits< InputIterator2 >::value_type(), typename std::iterator_traits< MaskIterator >::value_type value=typename std::iterator_traits< MaskIterator >::value_type(1))
Computes the pseudo normalized crosscorrelation between two sequences according to a mask sequence : ...
iterator end()
Returns a read/write iterator that points one past the last element in the Array. Iteration is done i...
Definition: Array.hpp:1084
Provides some statistics algorithms.
void unbiased_autocorrelation(RandomAccessIterator1 first, RandomAccessIterator1 last, RandomAccessIterator2 auto_first, RandomAccessIterator2 auto_last)
Computes the unbiased autocorrelations from lag -(last-first) to (last-first) of a range [first...
Real centered_crosscorrelation_mask(InputIterator1 first1, InputIterator1 last1, MaskIterator mask_first, InputIterator2 first2, typename std::iterator_traits< InputIterator1 >::value_type mean1=typename std::iterator_traits< InputIterator1 >::value_type(), typename std::iterator_traits< InputIterator1 >::value_type mean2=typename std::iterator_traits< InputIterator2 >::value_type(), typename std::iterator_traits< MaskIterator >::value_type value=typename std::iterator_traits< MaskIterator >::value_type(1))
Computes the centered crosscorrelation between two sequences according to a mask sequence ...
void multiplies(InputIterator1 __first1, InputIterator1 __last1, InputIterator2 __first2, OutputIterator __result)
Computes the pointwise product of two ranges.
void copy(_II first, _II last, _OI output_first)
Copy algorithm optimized for slip iterators.
Definition: copy_ext.hpp:177
Difference of Point2D class, specialization of DPoint<CoordType,DIM> with DIM = 2.
Definition: Array2d.hpp:129
Provides some macros which are used for using complex as real.
Value_T mean(InputIterator first, InputIterator last)
Computes the mean value of a range .
Definition: statistics.hpp:314
Real pseudo_normalized_crosscorrelation_mask(InputIterator1 first1, InputIterator1 last1, MaskIterator mask_first, InputIterator2 first2, typename std::iterator_traits< InputIterator1 >::value_type mean1=typename std::iterator_traits< InputIterator1 >::value_type(), typename std::iterator_traits< InputIterator1 >::value_type mean2=typename std::iterator_traits< InputIterator2 >::value_type(), typename std::iterator_traits< MaskIterator >::value_type value=typename std::iterator_traits< MaskIterator >::value_type(1))
Computes the pseudo normalized crosscorrelation between two sequences according to a mask sequence : ...
Real normalized_crosscorrelation_mask(InputIterator1 first1, InputIterator1 last1, MaskIterator mask_first, InputIterator2 first2, typename std::iterator_traits< InputIterator1 >::value_type mean1=typename std::iterator_traits< InputIterator1 >::value_type(), typename std::iterator_traits< InputIterator1 >::value_type mean2=typename std::iterator_traits< InputIterator2 >::value_type(), typename std::iterator_traits< MaskIterator >::value_type value=typename std::iterator_traits< MaskIterator >::value_type(1))
Computes the standard normalized crosscorrelation between two sequences according to a mask sequence ...
void set_coord(const CoordType &dx1, const CoordType &dx2)
Accessor/Mutator of the coordinates of DPoint2d.
Definition: DPoint2d.hpp:253
iterator end()
Returns a read/write iterator that points one past the last element in the Array2d. Iteration is done in ordinary element order.
Definition: Array2d.hpp:2352
HyperVolume< T > sqrt(const HyperVolume< T > &M)
void complex_fft(InputIter in_begin, InputIter in_end, OutputIter out_begin)
Computes the complex fft of a container.
Definition: FFT.hpp:2470
Provides some FFT algorithms.
Numerical matrix class. This container statisfies the BidirectionnalContainer concepts of the STL...
Provides some algorithms to computes arithmetical operations on ranges.
Provides SLIP error messages.
void fft_crosscorrelation3d(InputIterator3d1 in1_front_upper_left, InputIterator3d1 in1_back_bottom_right, InputIterator3d2 in2_front_upper_left, InputIterator3d2 in2_back_bottom_right, OutputIterator3d out_front_upper_left, OutputIterator3d out_back_bottom_right)
Computes the standard crosscorrelation between two 3D sequences using fft3d.
iterator3d back_bottom_right()
Returns a read/write iterator3d that points to the past the end element of the Array3d. It points to past the end element of the back bottom right element of the Array3d.
Definition: Array3d.hpp:4897
size_type rows() const
Returns the number of rows (first dimension size) in the Matrix.
Definition: Matrix.hpp:3499
Real pncc(InputIterator1 first, InputIterator1 last, InputIterator2 first2, typename std::iterator_traits< InputIterator1 >::value_type mean1=typename std::iterator_traits< InputIterator1 >::value_type(), typename std::iterator_traits< InputIterator1 >::value_type mean2=typename std::iterator_traits< InputIterator2 >::value_type())
Computes the pseudo normalized crosscorrelation between two sequences: .
std::iterator_traits< RandomAccessIterator1 >::value_type inner_product(RandomAccessIterator1 first1, RandomAccessIterator1 last1, RandomAccessIterator2 first2, typename std::iterator_traits< RandomAccessIterator1 >::value_type init=typename std::iterator_traits< RandomAccessIterator1 >::value_type())
Computes the inner_product of two ranges X and Y: .
Provides some convolution algorithms.
Real pncc_if(InputIterator1 first1, InputIterator1 last1, InputIterator2 first2, Predicate pred, typename std::iterator_traits< InputIterator1 >::value_type mean1=typename std::iterator_traits< InputIterator1 >::value_type(), typename std::iterator_traits< InputIterator1 >::value_type mean2=typename std::iterator_traits< InputIterator2 >::value_type())
Computes the standard pseudo normalized crosscorrelation between two sequences according to a Predica...
void ifft2d(InputBidirectionalIterator2d in_upper_left, InputBidirectionalIterator2d in_bottom_right, OutputBidirectionalIterator2d out_upper_left)
Computes the ifft2d of a container with 2d iterators.
Definition: FFT.hpp:2629
void fftshift3d(ForwardIterator3D front_upper_left, ForwardIterator3D back_bottom_right)
Performs a shift of a 2d container, for use with the fft2d and ifft2d functions, in order to move the...
Definition: FFT.hpp:230
Provides a class to manipulate 3d dynamic and generic arrays.
void normalized_crosscorrelation2d(InputIterator2d1 in1_upper_left, InputIterator2d1 in1_bottom_right, InputIterator2d2 in2_upper_left, InputIterator2d2 in2_bottom_right, OutputIterator2d out_upper_left, OutputIterator2d out_bottom_right)
Computes the normalized crosscorrelation between two Images.
void real_fft(InputIter in_begin, InputIter in_end, OutputIter out_begin)
Computes the real fft of a container.
Definition: FFT.hpp:2430
void std_crosscorrelation2d(InputIterator2d1 in1_upper_left, InputIterator2d1 in1_bottom_right, InputIterator2d2 in2_upper_left, InputIterator2d2 in2_bottom_right, OutputIterator2d out_upper_left, OutputIterator2d out_bottom_right)
Computes the standard crosscorrelation between two Images.
void fft_crosscorrelation2d_same(InputIterator2d1 in1_upper_left, InputIterator2d1 in1_bottom_right, InputIterator2d2 in2_upper_left, InputIterator2d2 in2_bottom_right, OutputIterator2d out_upper_left, OutputIterator2d out_bottom_right)
Computes the standard crosscorrelation between two 2D sequences using fft2d.
void biased_autocorrelation_matrix(RandomAccessIterator first, RandomAccessIterator last, Matrix &A)
Constructs the biased autocorrelation matrix from given a range.
Real centered_crosscorrelation_if(InputIterator1 first1, InputIterator1 last1, InputIterator2 first2, Predicate pred, typename std::iterator_traits< InputIterator1 >::value_type mean1=typename std::iterator_traits< InputIterator1 >::value_type(), typename std::iterator_traits< InputIterator1 >::value_type mean2=typename std::iterator_traits< InputIterator2 >::value_type())
Computes the standard crosscorrelation between two sequences according to a Predicate: ...
Real normalized_crosscorrelation_if(InputIterator1 first1, InputIterator1 last1, InputIterator2 first2, Predicate pred, typename std::iterator_traits< InputIterator1 >::value_type mean1=typename std::iterator_traits< InputIterator1 >::value_type(), typename std::iterator_traits< InputIterator1 >::value_type mean2=typename std::iterator_traits< InputIterator2 >::value_type())
Computes the standard normalized crosscorrelation between two sequences according to a Predicate: ** ...
iterator begin()
Returns a read/write iterator that points to the first element in the Array2d. Iteration is done in o...
Definition: Array2d.hpp:2345
void crosscorrelation2d(InputIterator2d1 in1_upper_left, InputIterator2d1 in1_bottom_right, InputIterator2d2 in2_upper_left, InputIterator2d2 in2_bottom_right, OutputIterator2d out_upper_left, OutputIterator2d out_bottom_right, CROSSCORRELATION_TYPE t)
Computes the crosscorrelation between two Images.
This is a linear (one-dimensional) dynamic template container. This container statisfies the RandomAc...
Definition: Array.hpp:94
const std::string CROSSCORRELATION_BAD_ALGO
Definition: error.hpp:83
void real_fft3d(InputBidirectionalIterator3d in_front_upper_left, InputBidirectionalIterator3d in_back_bottom_right, OutputBidirectionalIterator3d out_front_upper_left)
Computes the real fft3d of a container with 3d iterators.
Definition: FFT.hpp:2804
CROSSCORRELATION_TYPE
Choose between different types of crosscorrelation.
Real ncc_mask(InputIterator1 first1, InputIterator1 last1, MaskIterator mask_first, InputIterator2 first2, typename std::iterator_traits< InputIterator1 >::value_type mean1=typename std::iterator_traits< InputIterator1 >::value_type(), typename std::iterator_traits< InputIterator1 >::value_type mean2=typename std::iterator_traits< InputIterator2 >::value_type(), typename std::iterator_traits< MaskIterator >::value_type value=typename std::iterator_traits< MaskIterator >::value_type(1))
Computes the standard normalized crosscorrelation between two sequences according to a mask sequence ...
void fft_crosscorrelation(InputIterator1 first, InputIterator1 last, InputIterator2 first2, InputIterator2 last2, OutputIterator result_first, OutputIterator result_last)
Computes the standard crosscorrelation between two sequences using fft.
void toeplitz(RandomAccessIterator first, RandomAccessIterator last, Matrix &A)
Constructs the Toeplitz matrix from given a range.
int conj(const int arg)
This is a two-dimensional dynamic and generic container. This container statisfies the Bidirectionnal...
Definition: Array2d.hpp:135