annotate armadillo-2.4.4/include/armadillo_bits/op_min_meat.hpp @ 0:8b6102e2a9b0

Armadillo Library
author maxzanoni76 <max.zanoni@eecs.qmul.ac.uk>
date Wed, 11 Apr 2012 09:27:06 +0100
parents
children
rev   line source
max@0 1 // Copyright (C) 2008-2011 NICTA (www.nicta.com.au)
max@0 2 // Copyright (C) 2008-2011 Conrad Sanderson
max@0 3 //
max@0 4 // This file is part of the Armadillo C++ library.
max@0 5 // It is provided without any warranty of fitness
max@0 6 // for any purpose. You can redistribute this file
max@0 7 // and/or modify it under the terms of the GNU
max@0 8 // Lesser General Public License (LGPL) as published
max@0 9 // by the Free Software Foundation, either version 3
max@0 10 // of the License or (at your option) any later version.
max@0 11 // (see http://www.opensource.org/licenses for more info)
max@0 12
max@0 13
max@0 14 //! \addtogroup op_min
max@0 15 //! @{
max@0 16
max@0 17
max@0 18
max@0 19 template<typename eT>
max@0 20 arma_pure
max@0 21 inline
max@0 22 eT
max@0 23 op_min::direct_min(const eT* const X, const uword n_elem)
max@0 24 {
max@0 25 arma_extra_debug_sigprint();
max@0 26
max@0 27 eT min_val = priv::most_pos<eT>();
max@0 28
max@0 29 uword i,j;
max@0 30
max@0 31 for(i=0, j=1; j<n_elem; i+=2, j+=2)
max@0 32 {
max@0 33 const eT X_i = X[i];
max@0 34 const eT X_j = X[j];
max@0 35
max@0 36 if(X_i < min_val)
max@0 37 {
max@0 38 min_val = X_i;
max@0 39 }
max@0 40
max@0 41 if(X_j < min_val)
max@0 42 {
max@0 43 min_val = X_j;
max@0 44 }
max@0 45 }
max@0 46
max@0 47
max@0 48 if(i < n_elem)
max@0 49 {
max@0 50 const eT X_i = X[i];
max@0 51
max@0 52 if(X_i < min_val)
max@0 53 {
max@0 54 min_val = X_i;
max@0 55 }
max@0 56 }
max@0 57
max@0 58 return min_val;
max@0 59 }
max@0 60
max@0 61
max@0 62
max@0 63 template<typename eT>
max@0 64 inline
max@0 65 eT
max@0 66 op_min::direct_min(const eT* const X, const uword n_elem, uword& index_of_min_val)
max@0 67 {
max@0 68 arma_extra_debug_sigprint();
max@0 69
max@0 70 eT min_val = priv::most_pos<eT>();
max@0 71
max@0 72 uword best_index = 0;
max@0 73
max@0 74 uword i,j;
max@0 75
max@0 76 for(i=0, j=1; j<n_elem; i+=2, j+=2)
max@0 77 {
max@0 78 const eT X_i = X[i];
max@0 79 const eT X_j = X[j];
max@0 80
max@0 81 if(X_i < min_val)
max@0 82 {
max@0 83 min_val = X_i;
max@0 84 best_index = i;
max@0 85 }
max@0 86
max@0 87 if(X_j < min_val)
max@0 88 {
max@0 89 min_val = X_j;
max@0 90 best_index = j;
max@0 91 }
max@0 92 }
max@0 93
max@0 94
max@0 95 if(i < n_elem)
max@0 96 {
max@0 97 const eT X_i = X[i];
max@0 98
max@0 99 if(X_i < min_val)
max@0 100 {
max@0 101 min_val = X_i;
max@0 102 best_index = i;
max@0 103 }
max@0 104 }
max@0 105
max@0 106 index_of_min_val = best_index;
max@0 107
max@0 108 return min_val;
max@0 109 }
max@0 110
max@0 111
max@0 112
max@0 113 template<typename eT>
max@0 114 inline
max@0 115 eT
max@0 116 op_min::direct_min(const Mat<eT>& X, const uword row)
max@0 117 {
max@0 118 arma_extra_debug_sigprint();
max@0 119
max@0 120 const uword X_n_cols = X.n_cols;
max@0 121
max@0 122 eT min_val = priv::most_pos<eT>();
max@0 123
max@0 124 for(uword col=0; col<X_n_cols; ++col)
max@0 125 {
max@0 126 const eT tmp_val = X.at(row,col);
max@0 127
max@0 128 if(tmp_val < min_val)
max@0 129 {
max@0 130 min_val = tmp_val;
max@0 131 }
max@0 132 }
max@0 133
max@0 134 return min_val;
max@0 135 }
max@0 136
max@0 137
max@0 138
max@0 139 template<typename eT>
max@0 140 inline
max@0 141 eT
max@0 142 op_min::direct_min(const subview<eT>& X)
max@0 143 {
max@0 144 arma_extra_debug_sigprint();
max@0 145
max@0 146 const uword X_n_elem = X.n_elem;
max@0 147
max@0 148 eT min_val = priv::most_pos<eT>();
max@0 149
max@0 150 for(uword i=0; i<X_n_elem; ++i)
max@0 151 {
max@0 152 eT tmp_val = X[i];
max@0 153
max@0 154 if(tmp_val < min_val)
max@0 155 {
max@0 156 min_val = tmp_val;
max@0 157 }
max@0 158 }
max@0 159
max@0 160 return min_val;
max@0 161 }
max@0 162
max@0 163
max@0 164
max@0 165 template<typename eT>
max@0 166 inline
max@0 167 eT
max@0 168 op_min::direct_min(const diagview<eT>& X)
max@0 169 {
max@0 170 arma_extra_debug_sigprint();
max@0 171
max@0 172 const uword X_n_elem = X.n_elem;
max@0 173
max@0 174 eT min_val = priv::most_pos<eT>();;
max@0 175
max@0 176 for(uword i=0; i<X_n_elem; ++i)
max@0 177 {
max@0 178 eT tmp_val = X[i];
max@0 179
max@0 180 if(tmp_val < min_val)
max@0 181 {
max@0 182 min_val = tmp_val;
max@0 183 }
max@0 184 }
max@0 185
max@0 186 return min_val;
max@0 187 }
max@0 188
max@0 189
max@0 190
max@0 191 //! \brief
max@0 192 //! For each row or for each column, find the minimum value.
max@0 193 //! The result is stored in a dense matrix that has either one column or one row.
max@0 194 //! The dimension, for which the minima are found, is set via the min() function.
max@0 195 template<typename T1>
max@0 196 inline void op_min::apply(Mat<typename T1::elem_type>& out, const Op<T1,op_min>& in)
max@0 197 {
max@0 198 arma_extra_debug_sigprint();
max@0 199
max@0 200 typedef typename T1::elem_type eT;
max@0 201
max@0 202 const unwrap_check<T1> tmp(in.m, out);
max@0 203 const Mat<eT>& X = tmp.M;
max@0 204
max@0 205 const uword dim = in.aux_uword_a;
max@0 206 arma_debug_check( (dim > 1), "min(): incorrect usage. dim must be 0 or 1");
max@0 207
max@0 208 const uword X_n_rows = X.n_rows;
max@0 209 const uword X_n_cols = X.n_cols;
max@0 210
max@0 211 if(dim == 0) // min in each column
max@0 212 {
max@0 213 arma_extra_debug_print("op_min::apply(), dim = 0");
max@0 214
max@0 215 arma_debug_check( (X_n_rows == 0), "min(): given object has zero rows" );
max@0 216
max@0 217 out.set_size(1, X_n_cols);
max@0 218
max@0 219 eT* out_mem = out.memptr();
max@0 220
max@0 221 for(uword col=0; col<X_n_cols; ++col)
max@0 222 {
max@0 223 out_mem[col] = op_min::direct_min( X.colptr(col), X_n_rows );
max@0 224 }
max@0 225 }
max@0 226 else
max@0 227 if(dim == 1) // min in each row
max@0 228 {
max@0 229 arma_extra_debug_print("op_min::apply(), dim = 1");
max@0 230
max@0 231 arma_debug_check( (X_n_cols == 0), "min(): given object has zero columns" );
max@0 232
max@0 233 out.set_size(X_n_rows, 1);
max@0 234
max@0 235 eT* out_mem = out.memptr();
max@0 236
max@0 237 for(uword row=0; row<X_n_rows; ++row)
max@0 238 {
max@0 239 out_mem[row] = op_min::direct_min( X, row );
max@0 240 }
max@0 241 }
max@0 242 }
max@0 243
max@0 244
max@0 245
max@0 246 template<typename T>
max@0 247 inline
max@0 248 std::complex<T>
max@0 249 op_min::direct_min(const std::complex<T>* const X, const uword n_elem)
max@0 250 {
max@0 251 arma_extra_debug_sigprint();
max@0 252
max@0 253 uword index = 0;
max@0 254 T min_val = priv::most_pos<T>();
max@0 255
max@0 256 for(uword i=0; i<n_elem; ++i)
max@0 257 {
max@0 258 const T tmp_val = std::abs(X[i]);
max@0 259
max@0 260 if(tmp_val < min_val)
max@0 261 {
max@0 262 min_val = tmp_val;
max@0 263 index = i;
max@0 264 }
max@0 265 }
max@0 266
max@0 267 return X[index];
max@0 268 }
max@0 269
max@0 270
max@0 271
max@0 272 template<typename T>
max@0 273 inline
max@0 274 std::complex<T>
max@0 275 op_min::direct_min(const std::complex<T>* const X, const uword n_elem, uword& index_of_min_val)
max@0 276 {
max@0 277 arma_extra_debug_sigprint();
max@0 278
max@0 279 uword index = 0;
max@0 280 T min_val = priv::most_pos<T>();
max@0 281
max@0 282 for(uword i=0; i<n_elem; ++i)
max@0 283 {
max@0 284 const T tmp_val = std::abs(X[i]);
max@0 285
max@0 286 if(tmp_val < min_val)
max@0 287 {
max@0 288 min_val = tmp_val;
max@0 289 index = i;
max@0 290 }
max@0 291 }
max@0 292
max@0 293 index_of_min_val = index;
max@0 294
max@0 295 return X[index];
max@0 296 }
max@0 297
max@0 298
max@0 299
max@0 300 template<typename T>
max@0 301 inline
max@0 302 std::complex<T>
max@0 303 op_min::direct_min(const Mat< std::complex<T> >& X, const uword row)
max@0 304 {
max@0 305 arma_extra_debug_sigprint();
max@0 306
max@0 307 const uword X_n_cols = X.n_cols;
max@0 308
max@0 309 uword index = 0;
max@0 310 T min_val = priv::most_pos<T>();
max@0 311
max@0 312 for(uword col=0; col<X_n_cols; ++col)
max@0 313 {
max@0 314 const T tmp_val = std::abs(X.at(row,col));
max@0 315
max@0 316 if(tmp_val < min_val)
max@0 317 {
max@0 318 min_val = tmp_val;
max@0 319 index = col;
max@0 320 }
max@0 321 }
max@0 322
max@0 323 return X.at(row,index);
max@0 324 }
max@0 325
max@0 326
max@0 327
max@0 328 template<typename T>
max@0 329 inline
max@0 330 std::complex<T>
max@0 331 op_min::direct_min(const subview< std::complex<T> >& X)
max@0 332 {
max@0 333 arma_extra_debug_sigprint();
max@0 334
max@0 335 const uword X_n_elem = X.n_elem;
max@0 336 uword index = 0;
max@0 337 T min_val = priv::most_pos<T>();
max@0 338
max@0 339 for(uword i=0; i<X_n_elem; ++i)
max@0 340 {
max@0 341 const T tmp_val = std::abs(X[i]);
max@0 342
max@0 343 if(tmp_val < min_val)
max@0 344 {
max@0 345 min_val = tmp_val;
max@0 346 index = i;
max@0 347 }
max@0 348 }
max@0 349
max@0 350 return X[index];
max@0 351 }
max@0 352
max@0 353
max@0 354
max@0 355 template<typename T>
max@0 356 inline
max@0 357 std::complex<T>
max@0 358 op_min::direct_min(const diagview< std::complex<T> >& X)
max@0 359 {
max@0 360 arma_extra_debug_sigprint();
max@0 361
max@0 362 const uword X_n_elem = X.n_elem;
max@0 363 uword index = 0;
max@0 364 T min_val = priv::most_pos<T>();
max@0 365
max@0 366 for(uword i=0; i<X_n_elem; ++i)
max@0 367 {
max@0 368 const T tmp_val = std::abs(X[i]);
max@0 369
max@0 370 if(tmp_val < min_val)
max@0 371 {
max@0 372 min_val = tmp_val;
max@0 373 index = i;
max@0 374 }
max@0 375 }
max@0 376
max@0 377 return X[index];
max@0 378 }
max@0 379
max@0 380
max@0 381
max@0 382 //! @}