libcm is a C development framework with an emphasis on audio signal processing applications.
Você não pode selecionar mais de 25 tópicos Os tópicos devem começar com uma letra ou um número, podem incluir traços ('-') e podem ter até 35 caracteres.

cmVectOpsTemplateHdr.h 44KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717
  1. //( { label:misc desc:"Miscellaneous vector operations." kw:[vop] }
  2. // Compute the cummulative sum of sbp[dn]. Equivalent to Matlab cumsum().
  3. VECT_OP_TYPE* VECT_OP_FUNC(CumSum)(VECT_OP_TYPE* dbp, unsigned dn, const VECT_OP_TYPE* sbp );
  4. // Returns true if all values in each vector are equal.
  5. bool VECT_OP_FUNC(Equal)( const VECT_OP_TYPE* s0p, const VECT_OP_TYPE* s1p, unsigned sn );
  6. // Same as Matlab linspace() v[i] = i * (limit-1)/n
  7. VECT_OP_TYPE* VECT_OP_FUNC(LinSpace)( VECT_OP_TYPE* dbp, unsigned dn, VECT_OP_TYPE base, VECT_OP_TYPE limit );
  8. //======================================================================================================================
  9. //)
  10. //( { label:Print desc:"Vector printing functions." kw:[vop] }
  11. // Setting fieldWidth or decPltCnt to to negative values result in fieldWidth == 10 or decPlCnt == 4
  12. //
  13. void VECT_OP_FUNC(Printf)( cmRpt_t* rpt, unsigned rn, unsigned cn, const VECT_OP_TYPE* dbp, int fieldWidth, int decPlCnt, const char* fmt, unsigned flags );
  14. void VECT_OP_FUNC(Print)( cmRpt_t* rpt, unsigned rn, unsigned cn, const VECT_OP_TYPE* dbp );
  15. void VECT_OP_FUNC(PrintE)( cmRpt_t* rpt, unsigned rn, unsigned cn, const VECT_OP_TYPE* dbp );
  16. void VECT_OP_FUNC(PrintLf)( const char* label, cmRpt_t* rpt, unsigned rn, unsigned cn, const VECT_OP_TYPE* dbp, unsigned fieldWidth, unsigned decPlCnt, const char* fmt );
  17. void VECT_OP_FUNC(PrintL)( const char* label, cmRpt_t* rpt, unsigned rn, unsigned cn, const VECT_OP_TYPE* dbp );
  18. void VECT_OP_FUNC(PrintLE)( const char* label, cmRpt_t* rpt, unsigned rn, unsigned cn, const VECT_OP_TYPE* dbp );
  19. //======================================================================================================================
  20. //)
  21. //( { label:Normalization desc:"Normalization and standardization functions." kw:[vop] }
  22. // Normalize the vector of proabilities by dividing through by the sum.
  23. // This leaves the relative proportions of each value unchanged while producing a total probability of 1.0.
  24. //
  25. VECT_OP_TYPE* VECT_OP_FUNC(NormalizeProbabilityVV)(VECT_OP_TYPE* dbp, unsigned dn, const VECT_OP_TYPE* sbp);
  26. VECT_OP_TYPE* VECT_OP_FUNC(NormalizeProbability)(VECT_OP_TYPE* dbp, unsigned dn);
  27. VECT_OP_TYPE* VECT_OP_FUNC(NormalizeProbabilityN)(VECT_OP_TYPE* dbp, unsigned dn, unsigned stride);
  28. //
  29. // Standardize the columns of the matrix by subtracting the mean and dividing by the standard deviation.
  30. // uV[dcn] returns the mean of the data and is optional.
  31. // sdV[dcn] return the standard deviation of the data and is optional.
  32. VECT_OP_TYPE* VECT_OP_FUNC(StandardizeRows)( VECT_OP_TYPE* dbp, unsigned drn, unsigned dcn, VECT_OP_TYPE* uV, VECT_OP_TYPE* sdV );
  33. VECT_OP_TYPE* VECT_OP_FUNC(StandardizeCols)( VECT_OP_TYPE* dbp, unsigned drn, unsigned dcn, VECT_OP_TYPE* uV, VECT_OP_TYPE* sdV );
  34. //
  35. // Normalize by dividing through by the max. value.
  36. // dp[] ./= max(dp). Returns the index of the max value.
  37. unsigned VECT_OP_FUNC(NormToMax)( VECT_OP_TYPE* dp, unsigned dn );
  38. //
  39. // Normalize by dividing through by the max. absolute value.
  40. // db[] .*= fact / abs(max(dp));
  41. unsigned VECT_OP_FUNC(NormToAbsMax)( VECT_OP_TYPE* dp, unsigned dn, VECT_OP_TYPE fact );
  42. //======================================================================================================================
  43. //)
  44. //( { label:"Mean and variance" desc:"Compute mean and variance." kw:[vop] }
  45. VECT_OP_TYPE VECT_OP_FUNC(Mean)( const VECT_OP_TYPE* sp, unsigned sn );
  46. VECT_OP_TYPE VECT_OP_FUNC(MeanN)( const VECT_OP_TYPE* sp, unsigned sn, unsigned stride );
  47. //
  48. // Take the mean of each column/row of a matrix.
  49. // Set 'dim' to 0 to return mean of columns else return mean of rows.
  50. VECT_OP_TYPE* VECT_OP_FUNC(MeanM)( VECT_OP_TYPE* dp, const VECT_OP_TYPE* sp, unsigned srn, unsigned scn, unsigned dim );
  51. //
  52. // Take the mean of the first 'cnt' element of each column/row of a matrix.
  53. // Set 'dim' to 0 to return mean of columns else return mean of rows.
  54. // If 'cnt' is greater than the number of elements in the column/row then 'cnt' is
  55. // reduced to the number of elements in the column/row.
  56. VECT_OP_TYPE* VECT_OP_FUNC(MeanM2)( VECT_OP_TYPE* dp, const VECT_OP_TYPE* sp, unsigned srn, unsigned scn, unsigned dim, unsigned cnt );
  57. //
  58. // Find the mean of the data points returned by srcFuncPtr(argPtr,i) and return it in dp[dim].
  59. // 'dim' is both the size of dp[] and the length of each data point returned by srcFuncPtr().
  60. // srcFuncPtr() will be called 'cnt' times but it may return NULL on some calls if the associated
  61. // data point should not be included in the mean calculation.
  62. VECT_OP_TYPE* VECT_OP_FUNC(Mean2)( VECT_OP_TYPE* dp, const VECT_OP_TYPE* (*srcFuncPtr)(void* arg, unsigned idx ), unsigned dim, unsigned cnt, void* argPtr );
  63. //
  64. // avgPtr is optional - set to NULL to compute the average
  65. VECT_OP_TYPE VECT_OP_FUNC(Variance)( const VECT_OP_TYPE* sp, unsigned sn, const VECT_OP_TYPE* avgPtr );
  66. VECT_OP_TYPE VECT_OP_FUNC(VarianceN)(const VECT_OP_TYPE* sp, unsigned sn, unsigned stride, const VECT_OP_TYPE* avgPtr );
  67. //
  68. // Set dim=0 to return variance of columns otherwise return variance or rows.
  69. VECT_OP_TYPE* VECT_OP_FUNC(VarianceM)(VECT_OP_TYPE* dp, const VECT_OP_TYPE* sp, unsigned srn, unsigned scn, const VECT_OP_TYPE* avgPtr, unsigned dim );
  70. //======================================================================================================================
  71. //)
  72. //( { label:"Covariance" desc:"Matrix covariance" kw:[vop] }
  73. // Calculate the sample covariance matrix from a set of Gaussian distributed multidimensional data.
  74. // sp[dn,scn] is the data set.
  75. // dn is the dimensionality of the data.
  76. // scn is the count of data points
  77. // up[dn] is an optional mean vector. If up == NULL then the mean of the data is calculated internally.
  78. // selIdxV[scn] can be used to select a subset of datapoints to process.
  79. // If selIdxV[] is non-NULL then only columns where selIdxV[i]==selKey will be processed.
  80. //
  81. // dp[dn,dn] = covar( sp[dn,scn], u[dn] )
  82. void VECT_OP_FUNC(GaussCovariance)(VECT_OP_TYPE* dp, unsigned dn, const VECT_OP_TYPE* sp, unsigned scn, const VECT_OP_TYPE* up, const unsigned* selIdxV, unsigned selKey );
  83. // Calculate the sample covariance matrix.
  84. // dp[ dn*dn ] - output matrix
  85. // dn - dimensionality of the data
  86. // srcFuncPtr - User defined function which is called to return a pointer to a data vector at index 'idx'.
  87. // The returned data vector must contain 'dn' elements. The function should return NULL
  88. // if the data point associated with 'idx' should not be included in the covariance calculation.
  89. // sn - count of data vectors
  90. // userPtr - User arg. passed to srcFuncPtr.
  91. // uV[ dn ] - mean of the data set (optional)
  92. // Note that this function computes the covariance matrix in 2 serial passes (1 if the mean vector is given)
  93. // through the 'sn' data points.
  94. // The result of this function are identical to the octave cov() function.
  95. void VECT_OP_FUNC(GaussCovariance2)(VECT_OP_TYPE* dp, unsigned dn, const VECT_OP_TYPE* (*srcFuncPtr)(void* userPtr, unsigned idx), unsigned sn, void* userPtr, const VECT_OP_TYPE* uV, const unsigned* selIdxV, unsigned selKey );
  96. //======================================================================================================================
  97. //)
  98. //( { label:"Float point normal" desc:"Evaluate the 'normalness of floating point values." kw:[vop] }
  99. // Returns true if all values are 'normal' according the the C macro 'isnormal'.
  100. // This function will return false if any of the values are zero.
  101. bool VECT_OP_FUNC(IsNormal)( const VECT_OP_TYPE* sp, unsigned sn );
  102. // Returns true if all values are 'normal' or zero according the the C macro 'isnormal'.
  103. // This function accepts zeros as normal.
  104. bool VECT_OP_FUNC(IsNormalZ)(const VECT_OP_TYPE* sp, unsigned sn );
  105. // Set dp[dn] to the indexes of the non-normal numbers in sp[dn].
  106. // Returns the count of indexes stored in dp[].
  107. unsigned VECT_OP_FUNC(FindNonNormal)( unsigned* dp, unsigned dn, const VECT_OP_TYPE* sp );
  108. unsigned VECT_OP_FUNC(FindNonNormalZ)( unsigned* dp, unsigned dn, const VECT_OP_TYPE* sp );
  109. //======================================================================================================================
  110. //)
  111. //( { label:"Measure" desc:"Measure features of a vector." kw:[vop] }
  112. // Successive call to to ZeroCrossCount should preserve the value pointed to by delaySmpPtr.
  113. unsigned VECT_OP_FUNC(ZeroCrossCount)( const VECT_OP_TYPE* sp, unsigned n, VECT_OP_TYPE* delaySmpPtr);
  114. // Calculuate the sum of the squares of all elements in bp[bn].
  115. VECT_OP_TYPE VECT_OP_FUNC(SquaredSum)( const VECT_OP_TYPE* bp, unsigned bn );
  116. // sn must be <= wndSmpCnt. If sn < wndSmpCnt then sp[sn] is treated as a
  117. // a partially filled buffer padded with wndSmpCnt-sn zeros.
  118. // rms = sqrt( sum(sp[1:sn] .* sp[1:sn]) / wndSmpCnt )
  119. VECT_OP_TYPE VECT_OP_FUNC(RMS)( const VECT_OP_TYPE* sp, unsigned sn, unsigned wndSmpCnt );
  120. // This function handles the case were sn is not an integer multiple of
  121. // wndSmpCnt or hopSmpCnt. In this case the function computes zero
  122. // padded RMS values for windows which go past the end of sp[sn].
  123. VECT_OP_TYPE* VECT_OP_FUNC(RmsV)( VECT_OP_TYPE* dp, unsigned dn, const VECT_OP_TYPE* sp, unsigned sn, unsigned wndSmpCnt, unsigned hopSmpCnt );
  124. // Return the magnitude (Euclidean Norm) of a vector.
  125. VECT_OP_TYPE VECT_OP_FUNC(EuclidNorm)( const VECT_OP_TYPE* sp, unsigned sn );
  126. VECT_OP_TYPE VECT_OP_FUNC(AlphaNorm)(const VECT_OP_TYPE* sp, unsigned sn, VECT_OP_TYPE alpha );
  127. //======================================================================================================================
  128. //)
  129. //( { label:"Distance" desc:"Calculate various vector distances." kw:[vop] }
  130. // Return the Itakura-Saito distance between a modelled power spectrum (up) and another power spectrum (sp).
  131. VECT_OP_TYPE VECT_OP_FUNC(ItakuraDistance)( const VECT_OP_TYPE* up, const VECT_OP_TYPE* sp, unsigned sn );
  132. // Return the cosine distance between two vectors.
  133. VECT_OP_TYPE VECT_OP_FUNC(CosineDistance)( const VECT_OP_TYPE* s0P, const VECT_OP_TYPE* s1p, unsigned sn );
  134. // Return the Euclidean distance between two vectors
  135. VECT_OP_TYPE VECT_OP_FUNC(EuclidDistance)( const VECT_OP_TYPE* s0p, const VECT_OP_TYPE* s1p, unsigned sn );
  136. // Return the Manhattan distance between two vectors
  137. VECT_OP_TYPE VECT_OP_FUNC(L1Distance)( const VECT_OP_TYPE* s0p, const VECT_OP_TYPE* s1p, unsigned sn );
  138. // Return the Mahalanobis distance between a vector and the mean of the distribution.
  139. // The mean vector could be replaced with another vector drawn from the same distribution in which
  140. // case the returned value would reflect the distance between the two vectors.
  141. // 'sn' is the dimensionality of the data.
  142. // up[D] and invCovM[sn,sn] are the mean and inverse of the covariance matrix of the distribution from
  143. // which sp[D] is drawn.
  144. VECT_OP_TYPE VECT_OP_FUNC(MahalanobisDistance)( const VECT_OP_TYPE* sp, unsigned sn, const VECT_OP_TYPE* up, const VECT_OP_TYPE* invCovM );
  145. // Return the KL distance between two probability distributions up[sn] and sp[sn].
  146. // Since up[] and sp[] are probability distributions they must sum to 1.0.
  147. VECT_OP_TYPE VECT_OP_FUNC(KL_Distance)( const VECT_OP_TYPE* up, const VECT_OP_TYPE* sp, unsigned sn );
  148. // Return the KL distance between a prototype vector up[sn] and another vector sp[sn].
  149. // This function first normalizes the two vectors to sum to 1.0 before calling
  150. // VECT_OP_FUNC(KL_Distance)(up,sp,sn);
  151. VECT_OP_TYPE VECT_OP_FUNC(KL_Distance2)( const VECT_OP_TYPE* up, const VECT_OP_TYPE* sp, unsigned sn );
  152. // Measure the Euclidean distance between a vector and all the columns in a matrix.
  153. // If dv[scn] is no NULL then return the Euclidean distance from sv[scn] to each column of sm[srn,scn].
  154. // The function returns the index of the closest data point (column) in sm[].
  155. unsigned VECT_OP_FUNC(EuclidDistanceVM)( VECT_OP_TYPE* dv, const VECT_OP_TYPE* sv, const VECT_OP_TYPE* sm, unsigned srn, unsigned scn );
  156. // Measure the distance between each column in s0M[ rn, s0cn ] and
  157. // each column in s1M[rn, s1cn ]. If dM is non-NULL store the
  158. // result in dM[s1cn, s0cn]. The difference between s0M[:,0] and s1M[:,0]
  159. // is stored in dM[0,0], the diff. between s0M[:,1] and s1M[:,1] is stored
  160. // in dM[1,0], etc. If mvV[s0cn] is non-NULL then minV[i] is set with
  161. // the distance from s0M[:,i] to the nearest column in s1M[]. If miV[s0cn]
  162. // is non-NULL then it is set with the column index of s1M[] which is
  163. // closest to s0M[:,i]. In other words mvV[i] gives the distance to column
  164. // miV[i] from column s0M[:,i].
  165. // In those cases where the distane from a prototype (centroid) to the data point
  166. // is not the same as from the data point to the centroid then s1M[] is considered
  167. // to hold the prototypes and s0M[] is considered to hold the data points.
  168. // The distance function returns the distance from a prototype 'cV[dimN]' to
  169. // an datapoint dV[dimN]. 'dimN' is the dimensionality of the data vector
  170. // and is threfore equal to 'rn'.
  171. void VECT_OP_FUNC(DistVMM)(
  172. VECT_OP_TYPE* dM, // dM[s1cn,s0cn] return distance mtx (optional)
  173. VECT_OP_TYPE* mvV, // mvV[s0cn] distance to closest data point in s0M[]. (optional)
  174. unsigned* miV, // miV[s0cn] column index into s1M[] of closest data point to s0M[:,i]. (optional)
  175. unsigned rn, // dimensionality of the data and the row count for s0M[] and s1M[]
  176. const VECT_OP_TYPE* s0M, // s0M[rn,s0cn] contains one data point per column
  177. unsigned s0cn, // count of data points (count of columns in s0M[]
  178. const VECT_OP_TYPE* s1M, // s1M[rn,s1cn] contains one prototype per column
  179. unsigned s1cn, // count of prototypes (count of columns in s1m[]
  180. VECT_OP_TYPE (*distFunc)( void* userPtr, const VECT_OP_TYPE* cV, const VECT_OP_TYPE* dV, unsigned dimN ),
  181. void* userPtr );
  182. //======================================================================================================================
  183. //)
  184. //( { label:"Select columns" desc:"Select columns based on distance." kw:[vop] }
  185. // Select 'selIdxN' columns from sM[srn,scn].
  186. // dM[srn,selIdxN] receives copies of the selected columns.
  187. // selIdxV[selIdxN] receives the column indexes of the selected columns.
  188. // Both dM[] and selIdxV[] are optional.
  189. // In each case the first selected point is chosen at random.
  190. // SelectRandom() then selects the following selIdxN-1 points at random.
  191. // SelectMaxDist() selects the next selIdxN-1 points by selecting
  192. // the point whose combined distance to the previously selected points
  193. // is greatest. SelectMaxAvgDist() selectes the points whose combined
  194. // average distance is greatest relative the the previously selected
  195. // points.
  196. void VECT_OP_FUNC(SelectRandom)( VECT_OP_TYPE* dM, unsigned* selIdxV, unsigned selIdxN, const VECT_OP_TYPE* sM, unsigned srn, unsigned scn );
  197. void VECT_OP_FUNC(SelectMaxDist)( VECT_OP_TYPE* dM, unsigned* selIdxV, unsigned selIdxN, const VECT_OP_TYPE* sM, unsigned srn, unsigned scn, VECT_OP_TYPE (*distFunc)( void* userPtr, const VECT_OP_TYPE* s0V, const VECT_OP_TYPE* s1V, unsigned sn ), void* distUserPtr );
  198. void VECT_OP_FUNC(SelectMaxAvgDist)( VECT_OP_TYPE* dM, unsigned* selIdxV, unsigned selIdxN, const VECT_OP_TYPE* sM, unsigned srn, unsigned scn, VECT_OP_TYPE (*distFunc)( void* userPtr, const VECT_OP_TYPE* s0V, const VECT_OP_TYPE* s1V, unsigned sn ), void* distUserPtr );
  199. //======================================================================================================================
  200. //)
  201. //( { label:"Matrix multiplication" desc:"Various matrix multiplication operations." kw:[vop] }
  202. // Return the sum of the products (dot product)
  203. VECT_OP_TYPE VECT_OP_FUNC(MultSumVV)( const VECT_OP_TYPE* s0p, const VECT_OP_TYPE* s1p, unsigned sn );
  204. VECT_OP_TYPE VECT_OP_FUNC(MultSumVS)( const VECT_OP_TYPE* s0p, unsigned sn, VECT_OP_TYPE s );
  205. // Number of elements in the dest vector is expected to be the same
  206. // as the number of source matrix rows.
  207. // mcn gives the number of columns in the source matrix which is
  208. // expected to match the number of elements in the source vector.
  209. // dbp[dn,1] = mp[dn,mcn] * vp[mcn,1]
  210. VECT_OP_TYPE* VECT_OP_FUNC(MultVMV)( VECT_OP_TYPE* dbp, unsigned dn, const VECT_OP_TYPE* mp, unsigned mcn, const VECT_OP_TYPE* vp );
  211. // Multiply a row vector with a matrix to produce a row vector.
  212. // dbp[1,dn] = v[1,vn] * m[vn,dn]
  213. VECT_OP_TYPE* VECT_OP_FUNC(MultVVM)( VECT_OP_TYPE* dbp, unsigned dn, const VECT_OP_TYPE* vp, unsigned vn, const VECT_OP_TYPE* mp );
  214. // Same as MultVMtV() except M is transposed as part of the multiply.
  215. // mrn gives the number of rows in m[] and number of elements in vp[]
  216. // dpb[dn] = mp[mrn,dn] * vp[mrn]
  217. VECT_OP_TYPE* VECT_OP_FUNC(MultVMtV)( VECT_OP_TYPE* dbp, unsigned dn, const VECT_OP_TYPE* mp, unsigned mrn, const VECT_OP_TYPE* vp );
  218. // Same as MultVMV() but where the matrix is diagonal.
  219. VECT_OP_TYPE* VECT_OP_FUNC(MultDiagVMV)( VECT_OP_TYPE* dbp, unsigned dn, const VECT_OP_TYPE* mp, unsigned mcn, const VECT_OP_TYPE* vp );
  220. // Generalized matrix multiply.
  221. // If transposition is selected for M0 or M1 then the given dimension represent the size of the matrix 'after' the transposion.
  222. // d[drn,dcn] = alpha * op(m0[drn,m0cn_m1rn]) * op(m1[m0cn_m1rn,dcn]) + beta * d[drn,dcn]
  223. /// See enum { kTranpsoseM0Fl=0x01, kTransposeM1Fl=0x02 } in cmVectOps for flags.
  224. VECT_OP_TYPE* VECT_OP_FUNC(MultMMM1)(VECT_OP_TYPE* dbp, unsigned drn, unsigned dcn, VECT_OP_TYPE alpha, const VECT_OP_TYPE* m0, const VECT_OP_TYPE* m1, unsigned m0cn_m1rn, VECT_OP_TYPE beta, unsigned flags );
  225. // Same a VECT_OP_FUNC(MultMMM1) except allows the operation on a sub-matrix by providing the physical (memory) row count rather than the logical (matrix) row count.
  226. VECT_OP_TYPE* VECT_OP_FUNC(MultMMM2)(VECT_OP_TYPE* dbp, unsigned drn, unsigned dcn, VECT_OP_TYPE alpha, const VECT_OP_TYPE* m0, const VECT_OP_TYPE* m1, unsigned m0cn_m1rn, VECT_OP_TYPE beta, unsigned flags, unsigned dprn, unsigned m0prn, unsigned m1prn );
  227. // d[drn,dcn] = m0[drn,m0cn] * m1[m1rn,dcn]
  228. VECT_OP_TYPE* VECT_OP_FUNC(MultMMM)( VECT_OP_TYPE* dbp, unsigned drn, unsigned dcn, const VECT_OP_TYPE* m0, const VECT_OP_TYPE* m1, unsigned m0cn_m1rn );
  229. // same as MultMMM() except second source matrix is transposed prior to the multiply
  230. VECT_OP_TYPE* VECT_OP_FUNC(MultMMMt)(VECT_OP_TYPE* dbp, unsigned drn, unsigned dcn, const VECT_OP_TYPE* m0, const VECT_OP_TYPE* m1, unsigned m0cn_m1rn );
  231. //======================================================================================================================
  232. //)
  233. //( { label:"Linear algebra" desc:"Miscellaneous linear algebra operations. Determinant, Inversion, Cholesky decompostion. Linear system solver." kw:[vop] }
  234. // Initialize dbp[dn,dn] as a square symetric positive definite matrix using values
  235. // from a random uniform distribution. This is useful for initializing random
  236. // covariance matrices as used by multivariate Gaussian distributions
  237. // If t is non-NULL it must point to a block of scratch memory of t[dn,dn].
  238. // If t is NULL then scratch memory is internally allocated and deallocated.
  239. VECT_OP_TYPE* VECT_OP_FUNC(RandSymPosDef)( VECT_OP_TYPE* dbp, unsigned dn, VECT_OP_TYPE* t );
  240. // Compute the determinant of any square matrix.
  241. VECT_OP_TYPE VECT_OP_FUNC(DetM)( const VECT_OP_TYPE* sp, unsigned srn );
  242. // Compute the determinant of a diagonal matrix.
  243. VECT_OP_TYPE VECT_OP_FUNC(DetDiagM)( const VECT_OP_TYPE* sp, unsigned srn);
  244. // Compute the log determinant of any square matrix.
  245. VECT_OP_TYPE VECT_OP_FUNC(LogDetM)( const VECT_OP_TYPE* sp, unsigned srn );
  246. // Compute the log determinant of a diagonal matrix.
  247. VECT_OP_TYPE VECT_OP_FUNC(LogDetDiagM)( const VECT_OP_TYPE* sp, unsigned srn);
  248. // Compute the inverse of a square matrix. Returns NULL if the matrix is not invertable.
  249. // 'drn' is the dimensionality of the data.
  250. VECT_OP_TYPE* VECT_OP_FUNC(InvM)( VECT_OP_TYPE* dp, unsigned drn );
  251. // Compute the inverse of a diagonal matrix. Returns NULL if the matrix is not invertable.
  252. VECT_OP_TYPE* VECT_OP_FUNC(InvDiagM)( VECT_OP_TYPE* dp, unsigned drn );
  253. // Solve a linear system of the form AX=B where A[an,an] is square.
  254. // Since A is square B must have 'an' rows.
  255. // Result is returned in B.
  256. // Returns a pointer to B on success or NULL on fail.
  257. // NOTE: Both A and B are overwritten by this operation.
  258. VECT_OP_TYPE* VECT_OP_FUNC(SolveLS)( VECT_OP_TYPE* A, unsigned an, VECT_OP_TYPE* B, unsigned bcn );
  259. // Perform a Cholesky decomposition of the square symetric matrix U[un,un].
  260. // The factorization has the form: A=U'TU.
  261. // If the factorization is successful A is set to U and a pointer to A is returned.
  262. // Note that the lower triangle of A is not overwritten. See CholZ().
  263. // If the factorization fails NULL is returned.
  264. VECT_OP_TYPE* VECT_OP_FUNC(Chol)(VECT_OP_TYPE* A, unsigned an );
  265. // Same as Chol() but sets the lower triangle of U to zero.
  266. // This is equivalent ot the Matlab version.
  267. VECT_OP_TYPE* VECT_OP_FUNC(CholZ)(VECT_OP_TYPE* U, unsigned un );
  268. // Calculate the best fit line: b0 + b1*x_i through the points x_i,y_i.
  269. // Set x to NULL if it uses sequential integers [0,1,2,3...]
  270. void VECT_OP_FUNC(Lsq1)(const VECT_OP_TYPE* x, const VECT_OP_TYPE* y, unsigned n, VECT_OP_TYPE* b0, VECT_OP_TYPE* b1 );
  271. //======================================================================================================================
  272. //)
  273. //( { label:"Stretch/Shrink" desc:"Stretch or shrink a vector by resampling." kw:[vop] }
  274. // Return the average value of the contents of sbp[] between two fractional indexes
  275. VECT_OP_TYPE VECT_OP_FUNC(FracAvg)( double bi, double ei, const VECT_OP_TYPE* sbp, unsigned sn );
  276. // Shrinking function - Decrease the size of sbp[] by averaging blocks of values into single values in dbp[]
  277. VECT_OP_TYPE* VECT_OP_FUNC(DownSampleAvg)( VECT_OP_TYPE* dbp, unsigned dn, const VECT_OP_TYPE* sbp, unsigned sn );
  278. // Stretching function - linear interpolate between points in sbp[] to fill dbp[] ... where dn > sn
  279. VECT_OP_TYPE* VECT_OP_FUNC(UpSampleInterp)( VECT_OP_TYPE* dbp, unsigned dn, const VECT_OP_TYPE* sbp, unsigned sn );
  280. // Stretch or shrink the sbp[] to fit into dbp[]
  281. VECT_OP_TYPE* VECT_OP_FUNC(FitToSize)( VECT_OP_TYPE* dbp, unsigned dn, const VECT_OP_TYPE* sbp, unsigned sn );
  282. // Stretch or shrink sV[] to fit into dV[] using a simple linear mapping.
  283. // When stretching (sn<dn) each source element is repeated dn/sn times
  284. // and the last fraction position is interpolated. When shrinking
  285. // (sn>dn) each dest value is formed by the average of sequential segments
  286. // of sn/dn source elements. Fractional values are used at the beginning
  287. // and end of each segment.
  288. VECT_OP_TYPE* VECT_OP_FUNC(LinearMap)(VECT_OP_TYPE* dV, unsigned dn, VECT_OP_TYPE* sV, unsigned sn );
  289. //======================================================================================================================
  290. //)
  291. //( { label:"Random number generation" desc:"Generate random numbers." kw:[vop] }
  292. // Generate a vector of uniformly distributed random numbers in the range minVal to maxVal.
  293. VECT_OP_TYPE* VECT_OP_FUNC(Random)( VECT_OP_TYPE* dbp, unsigned dn, VECT_OP_TYPE minVal, VECT_OP_TYPE maxVal );
  294. // Generate dn random numbers integers between 0 and wn-1 based on a the relative
  295. // weights in wp[wn]. Note thtat the weights do not have to sum to 1.0.
  296. unsigned* VECT_OP_FUNC(WeightedRandInt)( unsigned* dbp, unsigned dn, const VECT_OP_TYPE* wp, unsigned wn );
  297. // Generate a vector of normally distributed univariate random numbers
  298. VECT_OP_TYPE* VECT_OP_FUNC(RandomGauss)( VECT_OP_TYPE* dbp, unsigned dn, VECT_OP_TYPE mean, VECT_OP_TYPE var );
  299. // Generate a vector of normally distributed univariate random numbers where each value has been drawn from a
  300. // seperately parameterized Gaussian distribution. meanV[] and varV[] must both contain dn velues.
  301. VECT_OP_TYPE* VECT_OP_FUNC(RandomGaussV)( VECT_OP_TYPE* dbp, unsigned dn, const VECT_OP_TYPE* meanV, const VECT_OP_TYPE* varV );
  302. // Generate a matrix of multi-dimensional random values. Each column represents a single vector value and each row contains a dimension.
  303. // meanV[] and varV[] must both contain drn elements where each meanV[i],varV[i] pair parameterize one dimensions Gaussian distribution.
  304. VECT_OP_TYPE* VECT_OP_FUNC(RandomGaussM)( VECT_OP_TYPE* dbp, unsigned drn, unsigned dcn, const VECT_OP_TYPE* meanV, const VECT_OP_TYPE* varV );
  305. VECT_OP_TYPE* VECT_OP_FUNC(RandomGaussDiagM)( VECT_OP_TYPE* dbp, unsigned drn, unsigned dcn, const VECT_OP_TYPE* meanV, const VECT_OP_TYPE* diagCovarM );
  306. // Generate a matrix of multivariate random values drawn from a normal distribution.
  307. // The dimensionality of the values are 'drn'.
  308. // The count of returned values is 'dcn'.
  309. // meanV[drn] and covarM[drn,drn] parameterize the normal distribution.
  310. // The covariance matrix must be symetric and positive definite.
  311. // t[(drn*drn) ] points to scratch memory or is set to NULL if the function should
  312. // allocate the memory internally.
  313. // Based on octave function mvrnd.m.
  314. VECT_OP_TYPE* VECT_OP_FUNC(RandomGaussNonDiagM)( VECT_OP_TYPE* dbp, unsigned drn, unsigned dcn, const VECT_OP_TYPE* meanV, const VECT_OP_TYPE* covarM, VECT_OP_TYPE* t );
  315. // Same as RandomGaussNonDiagM() except requires the upper trianglular
  316. // Cholesky factor of the covar matrix in 'uM'.
  317. VECT_OP_TYPE* VECT_OP_FUNC(RandomGaussNonDiagM2)( VECT_OP_TYPE* dbp, unsigned drn, unsigned dcn, const VECT_OP_TYPE* meanV, const VECT_OP_TYPE* uM );
  318. // Generate a matrix of N*K multi-dimensional data points.
  319. // Where D is the dimensionality of the data. (D == drn).
  320. // K is the number of multi-dimensional PDF's (clusters).
  321. // N is the number of data points to generate per cluster.
  322. // dbp[ D, N*K ] contains the returned data point.
  323. // The first N columns is associated with the cluster 0,
  324. // the next N columns is associated with cluster 1, ...
  325. // meanM[ D, K ] and varM[D,K] parameterize the generating PDF.s for each cluster
  326. VECT_OP_TYPE* VECT_OP_FUNC(RandomGaussMM)( VECT_OP_TYPE* dbp, unsigned drn, unsigned dcn, const VECT_OP_TYPE* meanM, const VECT_OP_TYPE* varM, unsigned K );
  327. // Evaluate the univariate normal distribution defined by 'mean' and 'stdDev'.
  328. VECT_OP_TYPE* VECT_OP_FUNC(GaussPDF)( VECT_OP_TYPE* dbp, unsigned dn, const VECT_OP_TYPE* sbp, VECT_OP_TYPE mean, VECT_OP_TYPE stdDev );
  329. // Evaluate a multivariate normal distribution defined by meanV[D] and covarM[D,D]
  330. // at the data points held in the columns of xM[D,N]. Return the evaluation
  331. // results in the vector yV[N]. D is the dimensionality of the data. N is the number of
  332. // data points to evaluate and values to return in yV[N].
  333. // Set diagFl to true if covarM is diagonal.
  334. // The function fails and returns false if the covariance matrix is singular.
  335. bool VECT_OP_FUNC(MultVarGaussPDF)( VECT_OP_TYPE* yV, const VECT_OP_TYPE* xM, const VECT_OP_TYPE* meanV, const VECT_OP_TYPE* covarM, unsigned D, unsigned N, bool diagFl );
  336. // Same as multVarGaussPDF[] except takes the inverse covar mtx invCovarM[D,D]
  337. // and log determinant of covar mtx.
  338. // Always returns yV[].
  339. VECT_OP_TYPE* VECT_OP_FUNC(MultVarGaussPDF2)( VECT_OP_TYPE* yV, const VECT_OP_TYPE* xM, const VECT_OP_TYPE* meanV, const VECT_OP_TYPE* invCovarM, VECT_OP_TYPE logDet, unsigned D, unsigned N, bool diagFl );
  340. // Same as multVarGaussPDF[] except uses a function to obtain the data vectors.
  341. // srcFunc() can filter the data points by returning NULL if the data vector at frmIdx should
  342. // not be evaluated against the PDF. In this case yV[frmIdx] will be set to 0.
  343. VECT_OP_TYPE* VECT_OP_FUNC(MultVarGaussPDF3)(
  344. VECT_OP_TYPE* yV,
  345. const VECT_OP_TYPE* (*srcFunc)(void* funcDataPtr, unsigned frmIdx ),
  346. void* funcDataPtr,
  347. const VECT_OP_TYPE* meanV,
  348. const VECT_OP_TYPE* invCovarM,
  349. VECT_OP_TYPE logDet,
  350. unsigned D,
  351. unsigned N,
  352. bool diagFl );
  353. //======================================================================================================================
  354. //)
  355. //( { label:"Signal generators" desc:"Generate periodic signals." kw:[vop] }
  356. // The following functions all return the phase of the next value.
  357. unsigned VECT_OP_FUNC(SynthSine)( VECT_OP_TYPE* dbp, unsigned dn, unsigned phase, double srate, double hz );
  358. unsigned VECT_OP_FUNC(SynthCosine)( VECT_OP_TYPE* dbp, unsigned dn, unsigned phase, double srate, double hz );
  359. unsigned VECT_OP_FUNC(SynthSquare)( VECT_OP_TYPE* dbp, unsigned dn, unsigned phase, double srate, double hz, unsigned otCnt );
  360. unsigned VECT_OP_FUNC(SynthTriangle)( VECT_OP_TYPE* dbp, unsigned dn, unsigned phase, double srate, double hz, unsigned otCnt );
  361. unsigned VECT_OP_FUNC(SynthSawtooth)( VECT_OP_TYPE* dbp, unsigned dn, unsigned phase, double srate, double hz, unsigned otCnt );
  362. unsigned VECT_OP_FUNC(SynthPulseCos)( VECT_OP_TYPE* dbp, unsigned dn, unsigned phase, double srate, double hz, unsigned otCnt );
  363. unsigned VECT_OP_FUNC(SynthImpulse)( VECT_OP_TYPE* dbp, unsigned dn, unsigned phase, double srate, double hz );
  364. unsigned VECT_OP_FUNC(SynthPhasor)( VECT_OP_TYPE* dbp, unsigned dn, unsigned phase, double srate, double hz );
  365. // Return value should be passed back via delaySmp on the next call.
  366. VECT_OP_TYPE VECT_OP_FUNC(SynthPinkNoise)( VECT_OP_TYPE* dbp, unsigned dn, VECT_OP_TYPE delaySmp );
  367. //======================================================================================================================
  368. //)
  369. //( { label:"Exponential conversion" desc:"pow() and log() functions." kw:[vop] }
  370. // Raise dbp[] to the power 'expon'
  371. VECT_OP_TYPE* VECT_OP_FUNC(PowVS)( VECT_OP_TYPE* dbp, unsigned dn, VECT_OP_TYPE expon );
  372. VECT_OP_TYPE* VECT_OP_FUNC(PowVVS)( VECT_OP_TYPE* dbp, unsigned dn, const VECT_OP_TYPE* sp, VECT_OP_TYPE expon );
  373. // Take the natural log of all values in sbp[dn]. It is allowable for sbp point to the same array as dbp=.
  374. VECT_OP_TYPE* VECT_OP_FUNC(LogV)( VECT_OP_TYPE* dbp, unsigned dn, const VECT_OP_TYPE* sbp );
  375. //======================================================================================================================
  376. //)
  377. //( { label:"dB Conversions" desc:"Convert vectors between dB,linear and power representations." kw:[vop] }
  378. // Convert a magnitude (amplitude) spectrum to/from decibels.
  379. // It is allowable for dbp==sbp.
  380. VECT_OP_TYPE* VECT_OP_FUNC(AmplToDbVV)( VECT_OP_TYPE* dbp, unsigned dn, const VECT_OP_TYPE* sbp, VECT_OP_TYPE minDb );
  381. VECT_OP_TYPE* VECT_OP_FUNC(DbToAmplVV)( VECT_OP_TYPE* dbp, unsigned dn, const VECT_OP_TYPE* sbp);
  382. VECT_OP_TYPE* VECT_OP_FUNC(PowToDbVV)( VECT_OP_TYPE* dbp, unsigned dn, const VECT_OP_TYPE* sbp, VECT_OP_TYPE minDb );
  383. VECT_OP_TYPE* VECT_OP_FUNC(DbToPowVV)( VECT_OP_TYPE* dbp, unsigned dn, const VECT_OP_TYPE* sbp);
  384. VECT_OP_TYPE* VECT_OP_FUNC(LinearToDb)( VECT_OP_TYPE* dbp, unsigned dn, const VECT_OP_TYPE* sp, VECT_OP_TYPE mult );
  385. VECT_OP_TYPE* VECT_OP_FUNC(dBToLinear)( VECT_OP_TYPE* dbp, unsigned dn, const VECT_OP_TYPE* sp, VECT_OP_TYPE mult );
  386. VECT_OP_TYPE* VECT_OP_FUNC(AmplitudeToDb)( VECT_OP_TYPE* dbp, unsigned dn, const VECT_OP_TYPE* sp );
  387. VECT_OP_TYPE* VECT_OP_FUNC(PowerToDb)( VECT_OP_TYPE* dbp, unsigned dn, const VECT_OP_TYPE* sp );
  388. VECT_OP_TYPE* VECT_OP_FUNC(dBToAmplitude)( VECT_OP_TYPE* dbp, unsigned dn, const VECT_OP_TYPE* sp );
  389. VECT_OP_TYPE* VECT_OP_FUNC(dBToPower)( VECT_OP_TYPE* dbp, unsigned dn, const VECT_OP_TYPE* sp );
  390. //======================================================================================================================
  391. //)
  392. //( { label:"DSP Windows" desc:"DSP windowing functions." kw:[vop] }
  393. VECT_OP_TYPE VECT_OP_FUNC(KaiserBetaFromSidelobeReject)( double sidelobeRejectDb );
  394. VECT_OP_TYPE VECT_OP_FUNC(KaiserFreqResolutionFactor)( double sidelobeRejectDb );
  395. VECT_OP_TYPE* VECT_OP_FUNC(Kaiser)( VECT_OP_TYPE* dbp, unsigned dn, double beta );
  396. VECT_OP_TYPE* VECT_OP_FUNC(Gaussian)(VECT_OP_TYPE* dbp, unsigned dn, double mean, double variance );
  397. VECT_OP_TYPE* VECT_OP_FUNC(Hamming)( VECT_OP_TYPE* dbp, unsigned dn );
  398. VECT_OP_TYPE* VECT_OP_FUNC(Hann)( VECT_OP_TYPE* dbp, unsigned dn );
  399. VECT_OP_TYPE* VECT_OP_FUNC(Triangle)(VECT_OP_TYPE* dbp, unsigned dn );
  400. // The MATLAB equivalent Hamming and Hann windows.
  401. //VECT_OP_TYPE* VECT_OP_FUNC(HammingMatlab)(VECT_OP_TYPE* dbp, unsigned dn );
  402. VECT_OP_TYPE* VECT_OP_FUNC(HannMatlab)( VECT_OP_TYPE* dbp, unsigned dn );
  403. // Simulates the MATLAB GaussWin function. Set arg to 2.5 to simulate the default arg
  404. // as used by MATLAB.
  405. VECT_OP_TYPE* VECT_OP_FUNC(GaussWin)( VECT_OP_TYPE* dbp, unsigned dn, double arg );
  406. //======================================================================================================================
  407. //)
  408. //( { label:"DSP Filters" desc:"DSP filtering functions." kw:[vop] }
  409. // Direct form II algorithm based on the MATLAB implmentation
  410. // http://www.mathworks.com/access/helpdesk/help/techdoc/ref/filter.html#f83-1015962
  411. // The only difference between this function and the equivalent MATLAB filter() function
  412. // is that the first feedforward coeff is given as a seperate value. The first b coefficient
  413. // in this function is therefore the same as the second coefficient in the MATLAB function.
  414. // and the first a[] coefficient (which is generally set to 1.0) is skipped.
  415. // Example:
  416. // Matlab: b=[.5 .4 .3] a=[1 .2 .1]
  417. // Equiv: b0 = .5 b=[ .4 .3] a=[ .2 .1];
  418. //
  419. // y[yn] - output vector
  420. // x[xn] - input vector. xn must be <= yn. if xn < yn then the end of y[] is set to zero.
  421. // b0 - signal scale. This can also be seen as b[0] (which is not included in b[])
  422. // b[dn] - feedforward coeff's b[1..dn-1]
  423. // a[dn] - feedback coeff's a[1..dn-1]
  424. // d[dn+1] - delay registers - note that this array must be one element longer than the coeff arrays.
  425. //
  426. VECT_OP_TYPE* VECT_OP_FUNC(Filter)( VECT_OP_TYPE* y, unsigned yn, const VECT_OP_TYPE* x, unsigned xn, cmReal_t b0, const cmReal_t* b, const cmReal_t* a, cmReal_t* d, unsigned dn );
  427. struct cmFilter_str;
  428. //typedef cmRC_t (*VECT_OP_FUNC(FiltExecFunc_t))( struct acFilter_str* f, const VECT_OP_TYPE* x, unsigned xn, VECT_OP_TYPE* y, unsigned yn );
  429. VECT_OP_TYPE* VECT_OP_FUNC(FilterFilter)(struct cmFilter_str* f, cmRC_t (*func)( struct cmFilter_str* f, const VECT_OP_TYPE* x, unsigned xn, VECT_OP_TYPE* y, unsigned yn ), const cmReal_t bb[], unsigned bn, const cmReal_t aa[], unsigned an, const VECT_OP_TYPE* x, unsigned xn, VECT_OP_TYPE* y, unsigned yn );
  430. // Compute the coefficients of a low/high pass FIR filter
  431. // wndV[dn] gives the window function used to truncate the ideal low-pass impulse response.
  432. // Set wndV to NULL to use a unity window.
  433. // See enum { kHighPass_LPSincFl=0x01, kNormalize_LPSincFl=0x02 } in cmVectOps.h
  434. VECT_OP_TYPE* VECT_OP_FUNC(LP_Sinc)(VECT_OP_TYPE* dp, unsigned dn, const VECT_OP_TYPE* wndV, double srate, double fcHz, unsigned flags );
  435. //======================================================================================================================
  436. //)
  437. //( { label:"Spectral Masking" desc:"A collection of spectral masking functions." kw:[vop] }
  438. // Compute a set of filterCnt mel filter masks for wieghting magnitude spectra consisting of binCnt bins.
  439. // The spectrum is divided into bandCnt equal bands in the mel domain
  440. // Each row of the matrix contains the mask for a single filter band consisting of binCnt elements.
  441. // See enum{ kShiftMelFl=0x01, kNormalizeMelFl=0x02 } in cmVectOps.h
  442. // Set kShiftMelFl to shift the mel bands onto the nearest FFT bin.
  443. // Set kNormalizeMelFl to normalize the combined filters for unity gain.
  444. VECT_OP_TYPE* VECT_OP_FUNC(MelMask)( VECT_OP_TYPE* maskMtx, unsigned bandCnt, unsigned binCnt, double srate, unsigned flags );
  445. // Fill binIdxV[bandCnt] and cntV[bandCnt] with a bin to band map.
  446. // binIdx[] contains the first (minimum) bin index for a given band.
  447. // cntV[] contains the count of bins for each band.
  448. // bandCnt is the number of bark bands to return
  449. // The function returns the actual number of bands mapped which will always be <= 23.
  450. unsigned VECT_OP_FUNC(BarkMap)(unsigned* binIdxV, unsigned* cntV, unsigned bandCnt, unsigned binCnt, double srate );
  451. // Calc a set of triangle fitler masks into each row of maskMtx.
  452. // maskMtx[ bandCnt, binCnt ] - result matrix
  453. // binHz - freq resolution of the output filters.
  454. // stSpread - Semi-tone spread above and below each center frequency (stSpread*2) is the total bandwidth.
  455. // (Only used if lowHzV or uprHzV are NULL)
  456. // lowHz[ bandCnt ] - set of upper frequency limits for each band.
  457. // ctrHz[ bandCnt ] set to the center value in Hz for each band
  458. // uprHz[ bandCnt ] - set of lower frequency limits for each band.
  459. // Note if lowHz[] and uprHz[] are set to NULL then stSpread is used to set the bandwidth of each band.
  460. VECT_OP_TYPE* VECT_OP_FUNC(TriangleMask)(VECT_OP_TYPE* maskMtx, unsigned bandCnt, unsigned binCnt, const VECT_OP_TYPE* ctrHzV, VECT_OP_TYPE binHz, VECT_OP_TYPE stSpread, const VECT_OP_TYPE* lowHzV, const VECT_OP_TYPE* uprHzV );
  461. // Calculate a set of Bark band triangle filters into maskMtx.
  462. // Each row of maskMtx contains the filter for one band.
  463. // maskMtx[ bandCnt, binCnt ]
  464. // bandCnt - the number of triangle bankds. If bandCnt is > 24 it will be reduced to 24.
  465. // binCnt - the number of bins in the filters.
  466. // binHz - the width of each bin in Hz.
  467. VECT_OP_TYPE* VECT_OP_FUNC(BarkMask)(VECT_OP_TYPE* maskMtx, unsigned bandCnt, unsigned binCnt, double binHz );
  468. // Terhardt 1979 (Calculating virtual pitch, Hearing Research #1, pp 155-182)
  469. // See enum { kNoTtmFlags=0, kModifiedTtmFl=0x01 } in cmVectOps.h
  470. VECT_OP_TYPE* VECT_OP_FUNC(TerhardtThresholdMask)(VECT_OP_TYPE* maskV, unsigned binCnt, double srate, unsigned flags);
  471. //Schroeder et al., 1979, JASA, Optimizing digital speech coders by exploiting masking properties of the human ear
  472. VECT_OP_TYPE* VECT_OP_FUNC(ShroederSpreadingFunc)(VECT_OP_TYPE* m, unsigned bandCnt, double srate);
  473. //======================================================================================================================
  474. //)
  475. //( { label:"Machine learning" desc:"K-means clustering and Viterbi algorithms." kw:[vop] }
  476. // Assign each data point to one of k clusters using an expectation-maximization algorithm.
  477. // k gives the number of clusters to identify
  478. // Each column of sp[ srn, scn ] contains a multidimensional data point.
  479. // srn therefore defines the dimensionality of the data.
  480. // Each column of centroidV[ srn, k ] is set to the centroid of each of k clusters.
  481. // classIdxV[ scn ] assigns the index (0 to k-1) of a cluster to each soure data point
  482. // The function returns the number of iterations required for the EM process to converge.
  483. // selIdxV[ scn ] is optional and contains a list of id's assoc'd with each column of sM.
  484. // selKey is a integer value.
  485. // If selIdxV is non-NULL then only columns of sM[] where selIdxV[] == selKey will be clustered.
  486. // All columns of sM[] where the associated column in selIdxV[] do not match will be ignored.
  487. // Set 'initFromCentroidFl' to true if the initial centroids should be taken from centroidM[].
  488. // otherwise the initial centroids are selected from 'k' random data points in sp[].
  489. // The distance function distFunc(cV,dV,dN) is called to determine the distance from a
  490. // centroid the centroid 'cV[dN]' to a data point 'dV[dN]'. 'dN' is the dimensionality of the
  491. // feature vector and is therefore equal to 'srn'.
  492. unsigned VECT_OP_FUNC(Kmeans)(
  493. unsigned* classIdxV,
  494. VECT_OP_TYPE* centroidM,
  495. unsigned k,
  496. const VECT_OP_TYPE* sp,
  497. unsigned srn,
  498. unsigned scn,
  499. const unsigned* selIdxV,
  500. unsigned selKey,
  501. bool initFromCentroidFl,
  502. VECT_OP_TYPE (*distFunc)( void* userPtr, const VECT_OP_TYPE* cV, const VECT_OP_TYPE* dV, unsigned dN ),
  503. void* userDistPtr );
  504. // 'srcFunc() should return NULL if the data point located at 'frmIdx' should not be included in the clustering.
  505. // Clustering is considered to be complete after 'maxIterCnt' iterations or when
  506. // 'deltaStopCnt' or fewer data points change class on a single iteration
  507. unsigned VECT_OP_FUNC(Kmeans2)(
  508. unsigned* classIdxV, // classIdxV[scn] - data point class assignments
  509. VECT_OP_TYPE* centroidM, // centroidM[srn,K] - cluster centroids
  510. unsigned K, // count of clusters
  511. const VECT_OP_TYPE* (*srcFunc)(void* userPtr, unsigned frmIdx ),
  512. unsigned srn, // dimensionality of each data point
  513. unsigned scn, // count of data points
  514. void* userSrcPtr, // callback data for srcFunc
  515. VECT_OP_TYPE (*distFunc)( void* userPtr, const VECT_OP_TYPE* cV, const VECT_OP_TYPE* dV, unsigned dN ),
  516. void* userDistPtr, // arg. to distFunc()
  517. int iterCnt, // max. number of iterations (-1 to ignore)
  518. int deltaStopCnt); // if less than deltaStopCnt data points change classes on a given iteration then convergence occurs.
  519. // Determine the most likely state sequece stateV[timeN] given a
  520. // transition matrix a[stateN,stateN],
  521. // observation probability matrix b[stateN,timeN] and
  522. // initial state probability vector phi[stateN].
  523. // a[i,j] is the probability of transitioning from state i to state j.
  524. // b[i,t] is the probability of state i emitting the obj t.
  525. void VECT_OP_FUNC(DiscreteViterbi)(unsigned* stateV, unsigned timeN, unsigned stateN, const VECT_OP_TYPE* phi, const VECT_OP_TYPE* a, const VECT_OP_TYPE* b );
  526. //======================================================================================================================
  527. //)
  528. //( { label:"Graphics" desc:"Graphics related algorithms." kw:[vop] }
  529. // Generate the set of coordinates which describe a circle with a center at x,y.
  530. // dbp[dn,2] must contain 2*dn elements. The first column holds the x coord and and the second holds the y coord.
  531. VECT_OP_TYPE* VECT_OP_FUNC(CircleCoords)( VECT_OP_TYPE* dbp, unsigned dn, VECT_OP_TYPE x, VECT_OP_TYPE y, VECT_OP_TYPE varX, VECT_OP_TYPE varY );
  532. // Clip the line defined by x0,y0 to x1,y1 into the rect defined by xMin,yMin xMax,yMax.
  533. bool VECT_OP_FUNC(ClipLine)( VECT_OP_TYPE* x0, VECT_OP_TYPE* y0, VECT_OP_TYPE* x1, VECT_OP_TYPE* y1, VECT_OP_TYPE xMin, VECT_OP_TYPE yMin, VECT_OP_TYPE xMax, VECT_OP_TYPE yMax );
  534. // Return true if the line defined by x0,y0 to x1,y1 intersects with
  535. // the rectangle formed by xMin,yMin - xMax,yMax
  536. bool VECT_OP_FUNC(IsLineInRect)( VECT_OP_TYPE x0, VECT_OP_TYPE y0, VECT_OP_TYPE x1, VECT_OP_TYPE y1, VECT_OP_TYPE xMin, VECT_OP_TYPE yMin, VECT_OP_TYPE xMax, VECT_OP_TYPE yMax );
  537. // Return the perpendicular distance from the line formed by x0,y0 and x1,y1
  538. // and the point px,py
  539. VECT_OP_TYPE VECT_OP_FUNC(PtToLineDistance)( VECT_OP_TYPE x0, VECT_OP_TYPE y0, VECT_OP_TYPE x1, VECT_OP_TYPE y1, VECT_OP_TYPE px, VECT_OP_TYPE py);
  540. //======================================================================================================================
  541. //)
  542. //( { label:"Miscellaneous DSP" desc:"Common DSP algorithms." kw:[vop] }
  543. // Compute the complex transient detection function from successive spectral frames.
  544. // The spectral magntidue mag0V precedes mag1V and the phase (radians) spectrum phs0V precedes the phs1V which precedes phs2V.
  545. // binCnt gives the length of each of the spectral vectors.
  546. VECT_OP_TYPE VECT_OP_FUNC(ComplexDetect)(const VECT_OP_TYPE* mag0V, const VECT_OP_TYPE* mag1V, const VECT_OP_TYPE* phs0V, const VECT_OP_TYPE* phs1V, const VECT_OP_TYPE* phs2V, unsigned binCnt );
  547. // Compute a set of DCT-II coefficients. Result dp[ coeffCnt, filtCnt ]
  548. VECT_OP_TYPE* VECT_OP_FUNC(DctMatrix)( VECT_OP_TYPE* dp, unsigned coeffCnt, unsigned filtCnt );
  549. // Set the indexes of local peaks greater than threshold in dbp[].
  550. // Returns the number of peaks in dbp[]
  551. // The maximum number of peaks from n source values is max(0,floor((n-1)/2)).
  552. // Note that peaks will never be found at index 0 or index sn-1.
  553. unsigned VECT_OP_FUNC(PeakIndexes)( unsigned* dbp, unsigned dn, const VECT_OP_TYPE* sbp, unsigned sn, VECT_OP_TYPE threshold );
  554. // Return the index of the bin containing v otherwise return kInvalidIdx if v is below sbp[0] or above sbp[ n-1 ]
  555. // The bin limits are contained in sbp[].
  556. // The value in spb[] are therefore expected to be in increasing order.
  557. // The value returned will be in the range 0:sn-1.
  558. unsigned VECT_OP_FUNC(BinIndex)( const VECT_OP_TYPE* sbp, unsigned sn, VECT_OP_TYPE v );
  559. // Given the points x0[xy0N],y0[xy0N] fill y1[i] with the interpolated value of y0[] at
  560. // x1[i]. Note that x0[] and x1[] must be increasing monotonic.
  561. // This function is similar to the octave interp1() function.
  562. void VECT_OP_FUNC(Interp1)(VECT_OP_TYPE* y1, const VECT_OP_TYPE* x1, unsigned xy1N, const VECT_OP_TYPE* x0, const VECT_OP_TYPE* y0, unsigned xy0N );
  563. //======================================================================================================================
  564. //)