fml  0.1-0
Fused Matrix Library
cpumat.hh
1 // This file is part of fml which is released under the Boost Software
2 // License, Version 1.0. See accompanying file LICENSE or copy at
3 // https://www.boost.org/LICENSE_1_0.txt
4 
5 #ifndef FML_CPU_CPUMAT_H
6 #define FML_CPU_CPUMAT_H
7 #pragma once
8 
9 
10 #include <cmath>
11 #include <cstdint>
12 #include <cstdlib>
13 #include <cstring>
14 #include <random>
15 #include <stdexcept>
16 
17 #include "../_internals/arraytools/src/arraytools.hpp"
18 
19 #include "../_internals/rand.hh"
20 #include "../_internals/omp.hh"
21 #include "../_internals/print.hh"
22 #include "../_internals/types.hh"
23 #include "../_internals/unimat.hh"
24 
25 #include "cpuvec.hh"
26 
27 
28 namespace fml
29 {
35  template <typename REAL>
36  class cpumat : public fml::unimat<REAL>
37  {
38  public:
39  cpumat();
40  cpumat(len_t nrows, len_t ncols);
41  cpumat(REAL *data, len_t nrows, len_t ncols, bool free_on_destruct=false);
42  cpumat(cpumat &&x);
43  cpumat(const cpumat &x);
44  ~cpumat();
45 
46  void resize(len_t nrows, len_t ncols);
47  void inherit(REAL *data, len_t nrows, len_t ncols, bool free_on_destruct=false);
48  void inherit(cpumat<REAL> &data);
49  cpumat<REAL> dupe() const;
50 
51  void print(uint8_t ndigits=4, bool add_final_blank=true) const;
52  void info() const;
53 
54  void fill_zero();
55  void fill_val(const REAL v);
56  void fill_linspace();
57  void fill_linspace(const REAL start, const REAL stop);
58  void fill_eye();
59  void fill_diag(const cpuvec<REAL> &v);
60  void fill_runif(const uint32_t seed, const REAL min=0, const REAL max=1);
61  void fill_runif(const REAL min=0, const REAL max=1);
62  void fill_rnorm(const uint32_t seed, const REAL mean=0, const REAL sd=1);
63  void fill_rnorm(const REAL mean=0, const REAL sd=1);
64 
65  void diag(cpuvec<REAL> &v);
66  void antidiag(cpuvec<REAL> &v);
67  void scale(const REAL s);
68  void rev_rows();
69  void rev_cols();
70 
71  bool any_inf() const;
72  bool any_nan() const;
73 
74  REAL get(const len_t i) const;
75  REAL get(const len_t i, const len_t j) const;
76  void set(const len_t i, const REAL v);
77  void set(const len_t i, const len_t j, const REAL v);
78  void get_row(const len_t i, cpuvec<REAL> &v) const;
79  void set_row(const len_t i, const cpuvec<REAL> &v);
80  void get_col(const len_t j, cpuvec<REAL> &v) const;
81  void set_col(const len_t i, const cpuvec<REAL> &v);
82 
83  bool operator==(const cpumat<REAL> &x) const;
84  bool operator!=(const cpumat<REAL> &x) const;
87 
88  protected:
89  bool free_on_destruct() const {return this->free_data;};
90  void dont_free_on_destruct() {this->free_data=false;};
91 
92  private:
93  void free();
94  void check_params(len_t nrows, len_t ncols);
95  };
96 }
97 
98 
99 
100 // -----------------------------------------------------------------------------
101 // public
102 // -----------------------------------------------------------------------------
103 
104 // constructors/destructor
105 
113 template <typename REAL>
115 {
116  this->m = 0;
117  this->n = 0;
118  this->data = NULL;
119 
120  this->free_data = true;
121 }
122 
123 
124 
137 template <typename REAL>
138 fml::cpumat<REAL>::cpumat(len_t nrows, len_t ncols)
139 {
140  check_params(nrows, ncols);
141 
142  this->m = nrows;
143  this->n = ncols;
144  this->free_data = true;
145 
146  if (this->m == 0 || this->n == 0)
147  return;
148 
149  size_t len = (size_t) nrows * ncols * sizeof(REAL);
150  this->data = (REAL*) std::malloc(len);
151  if (this->data == NULL)
152  throw std::bad_alloc();
153 }
154 
155 
156 
170 template <typename REAL>
171 fml::cpumat<REAL>::cpumat(REAL *data_, len_t nrows, len_t ncols, bool free_on_destruct)
172 {
173  check_params(nrows, ncols);
174 
175  this->m = nrows;
176  this->n = ncols;
177  this->data = data_;
178 
179  this->free_data = free_on_destruct;
180 }
181 
182 
183 
184 template <typename REAL>
186 {
187  this->m = x.nrows();
188  this->n = x.ncols();
189  this->data = x.data_ptr();
190 
191  this->free_data = x.free_on_destruct();
192  x.dont_free_on_destruct();
193 }
194 
195 
196 
197 template <typename REAL>
198 fml::cpumat<REAL>::cpumat(const cpumat<REAL> &x)
199 {
200  this->m = x.nrows();
201  this->n = x.ncols();
202  this->data.resize(this->m, this->n);
203 
204  size_t len = (size_t) this->m * this->n * sizeof(REAL);
205  std::memcpy(this->data, x.data_ptr(), len);
206 
207  this->free_data = true;
208 }
209 
210 
211 
212 template <typename REAL>
214 {
215  this->free();
216 }
217 
218 
219 
220 // memory management
221 
232 template <typename REAL>
233 void fml::cpumat<REAL>::resize(len_t nrows, len_t ncols)
234 {
235  check_params(nrows, ncols);
236 
237  const size_t len = (size_t) nrows * ncols * sizeof(REAL);
238  const size_t oldlen = (size_t) this->m * this->n * sizeof(REAL);
239 
240  if ((nrows == 0 || ncols == 0) || (len == oldlen))
241  {
242  this->m = nrows;
243  this->n = ncols;
244  return;
245  }
246 
247  void *realloc_ptr;
248  if (oldlen == 0)
249  realloc_ptr = malloc(len);
250  else
251  realloc_ptr = realloc(this->data, len);
252 
253  if (realloc_ptr == NULL)
254  throw std::bad_alloc();
255 
256  this->data = (REAL*) realloc_ptr;
257 
258  this->m = nrows;
259  this->n = ncols;
260 }
261 
262 
263 
276 template <typename REAL>
277 void fml::cpumat<REAL>::inherit(REAL *data, len_t nrows, len_t ncols, bool free_on_destruct)
278 {
279  check_params(nrows, ncols);
280 
281  this->free();
282 
283  this->m = nrows;
284  this->n = ncols;
285  this->data = data;
286 
287  this->free_data = free_on_destruct;
288 }
289 
290 
291 
293 template <typename REAL>
295 {
296  this->free();
297 
298  this->m = data_.nrows();
299  this->n = data_.ncols();
300  this->data = data_.data_ptr();
301 
302  this->free_data = false;
303 }
304 
305 
306 
308 template <typename REAL>
310 {
311  fml::cpumat<REAL> cpy(this->m, this->n);
312 
313  const size_t len = (size_t) this->m * this->n * sizeof(REAL);
314  memcpy(cpy.data_ptr(), this->data, len);
315 
316  return cpy;
317 }
318 
319 
320 
321 // printers
322 
329 template <typename REAL>
330 void fml::cpumat<REAL>::print(uint8_t ndigits, bool add_final_blank) const
331 {
332  for (len_t i=0; i<this->m; i++)
333  {
334  for (len_t j=0; j<this->n; j++)
335  this->printval(this->data[i + this->m*j], ndigits);
336 
337  fml::print::putchar('\n');
338  }
339 
340  if (add_final_blank)
341  fml::print::putchar('\n');
342 }
343 
344 
345 
347 template <typename REAL>
349 {
350  fml::print::printf("# cpumat");
351  fml::print::printf(" %dx%d", this->m, this->n);
352  fml::print::printf(" type=%s", typeid(REAL).name());
353  fml::print::printf("\n");
354 }
355 
356 
357 
358 // fillers
359 
361 template <typename REAL>
363 {
364  const size_t len = (size_t) this->m * this->n * sizeof(REAL);
365  memset(this->data, 0, len);
366 }
367 
368 
369 
375 template <typename REAL>
376 void fml::cpumat<REAL>::fill_val(const REAL v)
377 {
378  #pragma omp parallel for if((this->m)*(this->n) > fml::omp::OMP_MIN_SIZE)
379  for (len_t j=0; j<this->n; j++)
380  {
381  #pragma omp simd
382  for (len_t i=0; i<this->m; i++)
383  this->data[i + this->m*j] = v;
384  }
385 }
386 
387 
388 
395 template <typename T>
397 {
398  T start = 1;
399  T stop = (T) (this->m * this->n);
400  this->fill_linspace(start, stop);
401 }
402 
403 template <typename REAL>
404 void fml::cpumat<REAL>::fill_linspace(const REAL start, const REAL stop)
405 {
406  if (start == stop)
407  this->fill_val(start);
408  else
409  {
410  const REAL v = (stop-start)/((REAL) this->m*this->n - 1);
411 
412  #pragma omp parallel for if((this->m)*(this->n) > fml::omp::OMP_MIN_SIZE)
413  for (len_t j=0; j<this->n; j++)
414  {
415  #pragma omp simd
416  for (len_t i=0; i<this->m; i++)
417  {
418  const len_t ind = i + this->m*j;
419  this->data[ind] = v*((REAL) ind) + start;
420  }
421  }
422  }
423 }
424 
425 template <>
426 inline void fml::cpumat<int>::fill_linspace(const int start, const int stop)
427 {
428  if (start == stop)
429  this->fill_val(start);
430  else
431  {
432  const float v = (stop-start)/((float) this->m*this->n - 1);
433 
434  #pragma omp parallel for if((this->m)*(this->n) > fml::omp::OMP_MIN_SIZE)
435  for (len_t j=0; j<this->n; j++)
436  {
437  #pragma omp simd
438  for (len_t i=0; i<this->m; i++)
439  {
440  const len_t ind = i + this->m*j;
441  this->data[ind] = (int) roundf(v*((float) ind) + start);
442  }
443  }
444  }
445 }
446 
447 
448 
450 template <typename REAL>
452 {
453  cpuvec<REAL> v(1);
454  v.set(0, (REAL)1);
455  this->fill_diag(v);
456 }
457 
458 
459 
469 template <typename REAL>
471 {
472  this->fill_zero();
473 
474  REAL *v_d = v.data_ptr();
475  len_t min = std::min(this->m, this->n);
476 
477  #pragma omp for simd
478  for (len_t i=0; i<min; i++)
479  this->data[i + this->m*i] = v_d[i % v.size()];
480 }
481 
482 
483 
490 template <typename REAL>
491 void fml::cpumat<REAL>::fill_runif(const uint32_t seed, const REAL min, const REAL max)
492 {
493  std::mt19937 mt(seed);
494  static std::uniform_real_distribution<REAL> dist(min, max);
495 
496  for (len_t j=0; j<this->n; j++)
497  {
498  for (len_t i=0; i<this->m; i++)
499  this->data[i + this->m*j] = dist(mt);
500  }
501 }
502 
504 template <typename REAL>
505 void fml::cpumat<REAL>::fill_runif(const REAL min, const REAL max)
506 {
507  uint32_t seed = fml::rand::get_seed();
508  this->fill_runif(seed, min, max);
509 }
510 
511 
512 
519 template <typename REAL>
520 void fml::cpumat<REAL>::fill_rnorm(const uint32_t seed, const REAL mean, const REAL sd)
521 {
522  std::mt19937 mt(seed);
523  static std::normal_distribution<REAL> dist(mean, sd);
524 
525  for (len_t j=0; j<this->n; j++)
526  {
527  for (len_t i=0; i<this->m; i++)
528  this->data[i + this->m*j] = dist(mt);
529  }
530 }
531 
533 template <typename REAL>
534 void fml::cpumat<REAL>::fill_rnorm(const REAL mean, const REAL sd)
535 {
536  uint32_t seed = fml::rand::get_seed();
537  this->fill_rnorm(seed, mean, sd);
538 }
539 
540 
541 
555 template <typename REAL>
557 {
558  const len_t minmn = std::min(this->m, this->n);
559  v.resize(minmn);
560  REAL *v_ptr = v.data_ptr();
561 
562  #pragma omp parallel for simd if(minmn > fml::omp::OMP_MIN_SIZE)
563  for (len_t i=0; i<minmn; i++)
564  v_ptr[i] = this->data[i + this->m*i];
565 }
566 
567 
568 
583 template <typename REAL>
585 {
586  const len_t minmn = std::min(this->m, this->n);
587  v.resize(minmn);
588  REAL *v_ptr = v.data_ptr();
589 
590  #pragma omp parallel for simd if(minmn > fml::omp::OMP_MIN_SIZE)
591  for (len_t i=0; i<minmn; i++)
592  v_ptr[i] = this->data[(this->m-1-i) + this->m*i];
593 }
594 
595 
596 
602 template <typename REAL>
603 void fml::cpumat<REAL>::scale(const REAL s)
604 {
605  #pragma omp parallel for if((this->m)*(this->n) > fml::omp::OMP_MIN_SIZE)
606  for (len_t j=0; j<this->n; j++)
607  {
608  #pragma omp simd
609  for (len_t i=0; i<this->m; i++)
610  this->data[i + this->m*j] *= s;
611  }
612 }
613 
614 
615 
617 template <typename REAL>
619 {
620  for (len_t j=0; j<this->n; j++)
621  {
622  len_t last = this->m - 1;
623  for (len_t i=0; i<this->m/2; i++)
624  {
625  const REAL tmp = this->data[i + this->m*j];
626  this->data[i + this->m*j] = this->data[last + this->m*j];
627  this->data[last + this->m*j] = tmp;
628  }
629 
630  last--;
631  }
632 }
633 
634 
635 
637 template <typename REAL>
639 {
640  len_t last = this->n - 1;
641 
642  for (len_t j=0; j<this->n/2; j++)
643  {
644  for (len_t i=0; i<this->m; i++)
645  {
646  const REAL tmp = this->data[i + this->m*j];
647  this->data[i + this->m*j] = this->data[i + this->m*last];
648  this->data[i + this->m*last] = tmp;
649  }
650 
651  last--;
652  }
653 }
654 
655 
656 
658 template <typename REAL>
660 {
661  for (len_t j=0; j<this->n; j++)
662  {
663  for (len_t i=0; i<this->m; i++)
664  {
665  if (isinf(this->data[i + this->m*j]))
666  return true;
667  }
668  }
669 
670  return false;
671 }
672 
673 
674 
676 template <typename REAL>
678 {
679  for (len_t j=0; j<this->n; j++)
680  {
681  for (len_t i=0; i<this->m; i++)
682  {
683  if (isnan(this->data[i + this->m*j]))
684  return true;
685  }
686  }
687 
688  return false;
689 }
690 
691 
692 
693 // operators
694 
704 template <typename REAL>
705 REAL fml::cpumat<REAL>::get(const len_t i) const
706 {
707  this->check_index(i);
708  return this->data[i];
709 }
710 
719 template <typename REAL>
720 REAL fml::cpumat<REAL>::get(const len_t i, const len_t j) const
721 {
722  this->check_index(i, j);
723  return this->data[i + (this->m)*j];
724 }
725 
736 template <typename REAL>
737 void fml::cpumat<REAL>::set(const len_t i, const REAL v)
738 {
739  this->check_index(i);
740  this->data[i] = v;
741 }
742 
752 template <typename REAL>
753 void fml::cpumat<REAL>::set(const len_t i, const len_t j, const REAL v)
754 {
755  this->check_index(i, j);
756  this->data[i + (this->m)*j] = v;
757 }
758 
759 
760 
774 template <typename REAL>
775 void fml::cpumat<REAL>::get_row(const len_t i, cpuvec<REAL> &v) const
776 {
777  if (i < 0 || i >= this->m)
778  throw std::logic_error("invalid matrix row");
779 
780  v.resize(this->n);
781  REAL *v_d = v.data_ptr();
782 
783  #pragma omp parallel for simd if(this->n > fml::omp::OMP_MIN_SIZE)
784  for (len_t j=0; j<this->n; j++)
785  v_d[j] = this->data[i + this->m*j];
786 }
787 
788 
789 
800 template <typename REAL>
801 void fml::cpumat<REAL>::set_row(const len_t i, const cpuvec<REAL> &v)
802 {
803  if (i < 0 || i >= this->m)
804  throw std::logic_error("invalid matrix row");
805  if (v.size() != this->n)
806  throw std::runtime_error("non-conformable arguments");
807 
808  REAL *v_d = v.data_ptr();
809  #pragma omp parallel for simd if(this->n > fml::omp::OMP_MIN_SIZE)
810  for (len_t j=0; j<this->n; j++)
811  this->data[i + this->m*j] = v_d[j];
812 }
813 
814 
815 
829 template <typename REAL>
830 void fml::cpumat<REAL>::get_col(const len_t j, cpuvec<REAL> &v) const
831 {
832  if (j < 0 || j >= this->n)
833  throw std::logic_error("invalid matrix column");
834 
835  v.resize(this->m);
836  REAL *v_d = v.data_ptr();
837 
838  #pragma omp parallel for if(this->m > fml::omp::OMP_MIN_SIZE)
839  for (len_t i=0; i<this->m; i++)
840  v_d[i] = this->data[i + this->m*j];
841 }
842 
843 
844 
855 template <typename REAL>
856 void fml::cpumat<REAL>::set_col(const len_t j, const cpuvec<REAL> &v)
857 {
858  if (j < 0 || j >= this->n)
859  throw std::logic_error("invalid matrix column");
860  if (v.size() != this->m)
861  throw std::runtime_error("non-conformable arguments");
862 
863  REAL *v_d = v.data_ptr();
864  #pragma omp parallel for simd if(this->n > fml::omp::OMP_MIN_SIZE)
865  for (len_t i=0; i<this->m; i++)
866  this->data[i + this->m*j] = v_d[i];
867 }
868 
869 
870 
879 template <typename REAL>
881 {
882  if (this->m != x.nrows() || this->n != x.ncols())
883  return false;
884  else if (this->data == x.data_ptr())
885  return true;
886 
887  const REAL *x_d = x.data_ptr();
888  for (len_t j=0; j<this->n; j++)
889  {
890  for (len_t i=0; i<this->m; i++)
891  {
892  const REAL a = this->data[i + this->m*j];
893  const REAL b = x_d[i + this->m*j];
894  if (!arraytools::fltcmp::eq(a, b))
895  return false;
896  }
897  }
898 
899  return true;
900 }
901 
908 template <typename REAL>
910 {
911  return !(*this == x);
912 }
913 
914 
915 
922 template <typename REAL>
924 {
925  this->m = x.nrows();
926  this->n = x.ncols();
927  this->data = x.data_ptr();
928 
929  this->free_data = x.free_on_destruct();
930  x.dont_free_on_destruct();
931 
932  return *this;
933 }
934 
936 template <typename REAL>
938 {
939  this->m = x.nrows();
940  this->n = x.ncols();
941  this->data = x.data_ptr();
942 
943  this->free_data = false;
944 
945  return *this;
946 }
947 
948 
949 
950 // -----------------------------------------------------------------------------
951 // private
952 // -----------------------------------------------------------------------------
953 
954 template <typename REAL>
956 {
957  if (this->free_data && this->data)
958  {
959  std::free(this->data);
960  this->data = NULL;
961  }
962 }
963 
964 
965 
966 template <typename REAL>
967 void fml::cpumat<REAL>::check_params(len_t nrows, len_t ncols)
968 {
969  if (nrows < 0 || ncols < 0)
970  throw std::runtime_error("invalid dimensions");
971 }
972 
973 
974 #endif
fml::cpumat
Matrix class for data held on a single CPU.
Definition: cpumat.hh:36
fml::unimat
Base matrix class. Not meant for direct use. Instead see cpumat, gpumat, and mpimat.
Definition: unimat.hh:30
fml::cpumat::fill_eye
void fill_eye()
Set diagonal entries to 1 and non-diagonal entries to 0.
Definition: cpumat.hh:451
fml::cpumat::info
void info() const
Print some brief information about the object.
Definition: cpumat.hh:348
fml::cpumat::cpumat
cpumat()
Construct matrix object with no internal allocated storage.
Definition: cpumat.hh:114
fml::cpumat::fill_zero
void fill_zero()
Set all values to zero.
Definition: cpumat.hh:362
fml::univec::data_ptr
T * data_ptr()
Pointer to the internal array.
Definition: univec.hh:28
fml::cpumat::dupe
cpumat< REAL > dupe() const
Duplicate the object in a deep copy.
Definition: cpumat.hh:309
fml::cpumat::get_row
void get_row(const len_t i, cpuvec< REAL > &v) const
Get the specified row.
Definition: cpumat.hh:775
fml::cpumat::scale
void scale(const REAL s)
Multiply all values by the input value.
Definition: cpumat.hh:603
fml::unimat::nrows
len_t nrows() const
Number of rows.
Definition: unimat.hh:36
fml::cpuvec::resize
void resize(len_t size)
Resize the internal object storage.
Definition: cpuvec.hh:210
fml::cpumat::antidiag
void antidiag(cpuvec< REAL > &v)
Get the anti-diagonal entries, i.e. those on the bottom-left to top-right.
Definition: cpumat.hh:584
fml::cpumat::resize
void resize(len_t nrows, len_t ncols)
Resize the internal object storage.
Definition: cpumat.hh:233
fml::cpumat::get
REAL get(const len_t i) const
Get the specified value.
Definition: cpumat.hh:705
fml::cpuvec
Vector class for data held on a single CPU.
Definition: cpuvec.hh:31
fml::cpumat::fill_val
void fill_val(const REAL v)
Set all values to input value.
Definition: cpumat.hh:376
fml::unimat::ncols
len_t ncols() const
Number of columns.
Definition: unimat.hh:38
fml::cpumat::print
void print(uint8_t ndigits=4, bool add_final_blank=true) const
Print all values in the object.
Definition: cpumat.hh:330
fml::cpumat::operator=
cpumat< REAL > & operator=(cpumat< REAL > &x)
Operator that sets the LHS to a shallow copy of the input. Desctruction of the LHS object will not re...
Definition: cpumat.hh:923
fml::unimat::data_ptr
REAL * data_ptr()
Pointer to the internal array.
Definition: unimat.hh:40
fml::cpumat::fill_runif
void fill_runif(const uint32_t seed, const REAL min=0, const REAL max=1)
Set diagonal entries to 1 and non-diagonal entries to 0.
Definition: cpumat.hh:491
fml::cpumat::fill_diag
void fill_diag(const cpuvec< REAL > &v)
Set diagonal entries of the matrix to those in the vector.
Definition: cpumat.hh:470
fml::cpumat::rev_cols
void rev_cols()
Reverse the columns of the matrix.
Definition: cpumat.hh:638
fml::cpumat::operator==
bool operator==(const cpumat< REAL > &x) const
See if the two objects are the same.
Definition: cpumat.hh:880
fml
Core namespace.
Definition: dimops.hh:10
fml::univec::size
len_t size() const
Number of elements in the vector.
Definition: univec.hh:26
fml::cpuvec::set
void set(const len_t i, const T v)
Set the storage at the specified index with the provided value.
Definition: cpuvec.hh:540
fml::cpumat::set_row
void set_row(const len_t i, const cpuvec< REAL > &v)
Set the specified row.
Definition: cpumat.hh:801
fml::cpumat::fill_rnorm
void fill_rnorm(const uint32_t seed, const REAL mean=0, const REAL sd=1)
Set diagonal entries to 1 and non-diagonal entries to 0.
Definition: cpumat.hh:520
fml::cpumat::diag
void diag(cpuvec< REAL > &v)
Get the diagonal entries.
Definition: cpumat.hh:556
fml::cpumat::set
void set(const len_t i, const REAL v)
Set the storage at the specified index with the provided value.
Definition: cpumat.hh:737
fml::cpumat::fill_linspace
void fill_linspace()
Set values to linearly spaced numbers.
Definition: cpumat.hh:396
fml::cpumat::any_nan
bool any_nan() const
Are any values NaN?
Definition: cpumat.hh:677
fml::cpumat::rev_rows
void rev_rows()
Reverse the rows of the matrix.
Definition: cpumat.hh:618
fml::cpumat::any_inf
bool any_inf() const
Are any values infinite?
Definition: cpumat.hh:659
fml::cpumat::set_col
void set_col(const len_t i, const cpuvec< REAL > &v)
Set the specified column.
Definition: cpumat.hh:856
fml::cpumat::operator!=
bool operator!=(const cpumat< REAL > &x) const
See if the two objects are not the same. Uses same internal logic as the == method.
Definition: cpumat.hh:909
fml::cpumat::inherit
void inherit(REAL *data, len_t nrows, len_t ncols, bool free_on_destruct=false)
Set the internal object storage to the specified array.
Definition: cpumat.hh:277
fml::cpumat::get_col
void get_col(const len_t j, cpuvec< REAL > &v) const
Get the specified column.
Definition: cpumat.hh:830