Eigen_Tools.cc 22.9 KB
Newer Older
1
#include <Base/Security/Mandatory.hh>
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
/*===========================================================================*\
 *                                                                           *
 *                               CoMISo                                      *
 *      Copyright (C) 2008-2009 by Computer Graphics Group, RWTH Aachen      *
 *                           www.rwth-graphics.de                            *
 *                                                                           *
 *---------------------------------------------------------------------------* 
 *  This file is part of CoMISo.                                             *
 *                                                                           *
 *  CoMISo is free software: you can redistribute it and/or modify           *
 *  it under the terms of the GNU General Public License as published by     *
 *  the Free Software Foundation, either version 3 of the License, or        *
 *  (at your option) any later version.                                      *
 *                                                                           *
 *  CoMISo is distributed in the hope that it will be useful,                *
 *  but WITHOUT ANY WARRANTY; without even the implied warranty of           *
 *  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the            *
 *  GNU General Public License for more details.                             *
 *                                                                           *
 *  You should have received a copy of the GNU General Public License        *
 *  along with CoMISo.  If not, see <http://www.gnu.org/licenses/>.          *
 *                                                                           *
\*===========================================================================*/ 



//=============================================================================
//
//  CLASS Eigen_Tools - IMPLEMENTATION
//
//=============================================================================

//== COMPILE-TIME PACKAGE REQUIREMENTS ========================================
#include <CoMISo/Config/config.hh>
#if COMISO_EIGEN3_AVAILABLE


#define COMISO_Eigen_TOOLS_C

//== INCLUDES =================================================================

#include "Eigen_Tools.hh"
#include <queue>
#include <CoMISo/Utils/StopWatch.hh>
#include <CoMISo/Utils/VSToolsT.hh>
47
DISABLE_BANNED_COMPLIANCE
48
#include <gmm/gmm.h>
49
ENABLE_BANNED_COMPLIANCE
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782

//== NAMESPACES ===============================================================

namespace COMISO_EIGEN
{

//== IMPLEMENTATION ==========================================================


//-----------------------------------------------------------------------------

template<class MatrixT, class REALT, class INTT>
void get_ccs_symmetric_data( const MatrixT&      _mat,
                             const char          _uplo,
                             std::vector<REALT>& _values,
                             std::vector<INTT>&  _rowind,
                             std::vector<INTT>&  _colptr )
{
  // Assumes col major
  
   int m = _mat.innerSize();
   int n = _mat.outerSize();

   _values.resize( 0 );
   _rowind.resize( 0 );
   _colptr.resize( 0 );

   INTT iv( 0 );

   typedef typename MatrixT::InnerIterator It;

   switch ( _uplo )
   {
      case 'l':
      case 'L':
         // for all columns
         for ( int i=0; i<n; ++i )
         {
            _colptr.push_back( iv );

            // row it
            It it(_mat, i);

            for( ; it; ++it)
            {
              if( it.index() >= i )
              {
                _values.push_back( it.value());
                _rowind.push_back( it.index());
                ++iv;
              }
            }
         }
         _colptr.push_back( iv );
         break;

      case 'u':
      case 'U':
         // for all columns
         for ( int i=0; i<n; ++i )
         {
            _colptr.push_back( iv );

            // row it
            It it(_mat, i);

            for( ; it; ++it)
            {
              if( it.index() <= i )
              {
                _values.push_back( it.value());
                _rowind.push_back( it.index());
                ++iv;
              }
            }
         }
         _colptr.push_back( iv );
         break;

      case 'c':
      case 'C':
         // for all columns
         for ( int i=0; i<n; ++i )
         {
            _colptr.push_back( iv );

            // row it
            It it(_mat, i);

            for( ; it; ++it)
            {
              _values.push_back( it.value());
              _rowind.push_back( it.index());
              ++iv;
            }
         }
         _colptr.push_back( iv );
         break;

      default:
         std::cerr << "ERROR: parameter uplo must bei either 'U' or 'L' or 'C'!!!\n";
         break;
   }
}


//-----------------------------------------------------------------------------


// inspect the matrix: dimension, symmetry, zero_rows, zero_cols, nnz, max, min, max_abs, min_abs, NAN, INF
template<class MatrixT>
void inspect_matrix( const MatrixT& _A)
{


  std::cerr << "################### INSPECT MATRIX ##################\n";
  std::cerr << "#outer size  : " << _A.outerSize() << std::endl;
  std::cerr << "#inner size  : " << _A.innerSize() << std::endl;
  std::cerr << "#rows        : " << _A.rows() << std::endl;
  std::cerr << "#cols        : " << _A.cols() << std::endl;
  std::cerr << "#nonzeros    : " << _A.nonZeros() << std::endl;
  std::cerr << "#nonzeros/row: " << (double(_A.nonZeros())/double(_A.rows())) << std::endl;
  std::cerr << "symmetric    : " << is_symmetric( _A) << std::endl;

  MatrixT trans( _A.transpose());

  int zero_rows = 0;
  int zero_cols = 0;

  for(int i=0; i<_A.outerSize(); ++i)
  {
    typename MatrixT::InnerIterator it(_A, i); 
    if( !it) ++zero_rows;
  }

  for(int i=0; i<trans.outerSize(); ++i)
  {
    typename MatrixT::InnerIterator it(trans, i); 
    if( !it) ++zero_cols;
  }

  std::cerr << "zero rows    : " << zero_rows << std::endl;
  std::cerr << "zero cols    : " << zero_cols << std::endl;

  typedef typename MatrixT::Scalar Scalar;
  Scalar vmin     = std::numeric_limits<Scalar>::max();
  Scalar vmax     = std::numeric_limits<Scalar>::min();
  Scalar vmin_abs = std::numeric_limits<Scalar>::max();
  Scalar vmax_abs = 0;

  int n_nan = 0;
  int n_inf = 0;
  
  // inspect elements
  for(int i=0; i<_A.outerSize(); ++i)
  {
    typename MatrixT::InnerIterator it( _A, i);
    
    for(; it ; ++it)
    {
      if( it.value() < vmin ) vmin = it.value();
      if( it.value() > vmax ) vmax = it.value();

      if( fabs(it.value()) < vmin_abs) vmin_abs = fabs(it.value());
      if( fabs(it.value()) > vmax_abs) vmax_abs = fabs(it.value());

      if( std::isnan(it.value())) ++n_nan;
      if( std::isinf(it.value())) ++n_inf;
    }
  }
  
  std::cerr << "min  val     : " << vmin << std::endl;
  std::cerr << "max  val     : " << vmax << std::endl;
  std::cerr << "min |val|    : " << vmin_abs << std::endl;
  std::cerr << "max |val|    : " << vmax_abs << std::endl;
  std::cerr << "#nan         : " << n_nan << std::endl;
  std::cerr << "#inf         : " << n_inf << std::endl;
  
  std::cerr << "min eval     : " << "..." << std::endl;
  std::cerr << "max eval     : " << "..." << std::endl;
  std::cerr << "min|eval|    : " << "..." << std::endl;
  std::cerr << "max|eval|    : " << "..." << std::endl;
}

//-----------------------------------------------------------------------------


// symmetric ?
template<class MatrixT>
bool is_symmetric( const MatrixT& _A)
{
  typedef typename MatrixT::InnerIterator It;
  typedef typename MatrixT::Scalar Scalar;

  int nouter( _A.outerSize());
  int ninner( _A.innerSize());
  
  if( nouter != ninner )
    return false;

  bool symmetric(true);

  for( int c = 0; c < nouter; ++c)
  {
    for( It it(_A,c); it; ++it)
    {
      int r(it.index());

      Scalar val(it.value());

      // find diagonal partner element
      bool found(false);
      for( It dit(_A,r); dit; ++dit)
      {
        if( dit.index() < c )
        {}
        else if( dit.index() == c)
        {
          if( dit.value() == val) 
            found = true;
          break;
        }
        else 
        {
          break;
        }
      }
      if( !found) 
      {
        symmetric = false;
        break;
      }
    }
  }
  return symmetric;
}

//-----------------------------------------------------------------------------


template< class Eigen_MatrixT, class IntT >
void permute( const Eigen_MatrixT& _QR, const std::vector< IntT>& _Pvec, Eigen_MatrixT& _A)
{
#ifdef COMISO_EIGEN3_AVAILABLE
  typedef typename Eigen_MatrixT::Scalar Scalar;

  int m = _QR.innerSize();
  int n = _QR.outerSize();
  
  if( _Pvec.size() == 0)
  {
    _A = _QR;
    return;
  }
    
  if( _Pvec.size() != (size_t)_QR.cols() && _Pvec.size() != 0)
  {
    std::cerr << __FUNCTION__ << " wrong size of permutation vector, should have #cols length (or zero)" << std::endl;
  }

  // build sparse permutation matrix
  typedef Eigen::Triplet< Scalar > Triplet;
  std::vector< Triplet > triplets;
  triplets.reserve(_QR.nonZeros());
  _A = Eigen_MatrixT( m, n);

  typedef typename Eigen_MatrixT::InnerIterator It;

  
  for( int c = 0; c < n; ++c) // cols
  {
    for( It it(_QR,c); it; ++it) // rows
    {
      int r(it.index());

      Scalar val(it.value());

      int newcol( _Pvec[c]);

      triplets.push_back( Triplet( r, newcol, val));
    }
  }
  _A.setFromTriplets( triplets.begin(), triplets.end());
#endif
}

//-----------------------------------------------------------------------------

#if COMISO_SUITESPARSE_AVAILABLE

/// Eigen to Cholmod_sparse interface
template<class MatrixT>
void cholmod_to_eigen( const cholmod_sparse& _AC, MatrixT& _A)
{
#ifdef COMISO_EIGEN3_AVAILABLE
  // initialize dimensions
  typedef typename MatrixT::Scalar Scalar;
  typedef Eigen::Triplet< Scalar > Triplet;
  size_t nzmax( _AC.nzmax);
  std::cerr << __FUNCTION__ << " row " << _AC.nrow << " col " << _AC.ncol << " stype " << _AC.stype << std::endl;
  _A = MatrixT(_AC.nrow, _AC.ncol);
  std::vector< Triplet > triplets;
  triplets.reserve(nzmax);

  if(!_AC.packed)
  {
    std::cerr << "Warning: " << __FUNCTION__ << " does not support unpacked matrices yet!!!" << std::endl;
    return;
  }

  // Pointer to data
  double* X((double*)_AC.x);

  // complete matrix stored
  if(_AC.stype == 0)
  {
    // which integer type?
    if(_AC.itype == CHOLMOD_LONG)
    {
      UF_long* P((UF_long*)_AC.p);
      UF_long* I((UF_long*)_AC.i);

      for(UF_long i=0; i<(UF_long)_AC.ncol; ++i)
        for(UF_long j= P[i]; j< P[i+1]; ++j)
          //_A( I[j], i) += X[j]; // += really needed?
          triplets.push_back( Triplet( I[j], i, X[j]));
    }
    else
    {
      int* P((int*)_AC.p);
      int* I((int*)_AC.i);

      for(int i=0; i<(int)_AC.ncol; ++i)
        for(int j= P[i]; j< P[i+1]; ++j)
          triplets.push_back( Triplet( I[j], i, X[j]));
      //_A( I[j], i) += X[j];
    }

  }
  else // only upper or lower diagonal stored
  {
    // which integer type?
    if(_AC.itype == CHOLMOD_LONG)
    {
      UF_long* P((UF_long*)_AC.p);
      UF_long* I((UF_long*)_AC.i);

      for(UF_long i=0; i<(UF_long)_AC.ncol; ++i)
        for(UF_long j=P[i]; j<P[i+1]; ++j)
        {
          //_A(I[j], i) += X[j];
          triplets.push_back( Triplet( I[j], i, X[j]));

          // add up symmetric part
          if( I[j] != i)
            triplets.push_back( Triplet( i, I[j], X[j]));
          //_A(i,I[j]) += X[j];
        }
    }
    else
    {
      int* P((int*)_AC.p);
      int* I((int*)_AC.i);

      for(int i=0; i<(int)_AC.ncol; ++i)
        for(int j=P[i]; j<P[i+1]; ++j)
        {
          //_A(I[j], i) += X[j];
          triplets.push_back( Triplet( I[j], i, X[j]));

          // add up symmetric part
          if( I[j] != i)
            //  _A(i,I[j]) += X[j];
            triplets.push_back( Triplet( i, I[j], X[j]));
        }
    }
  }
  _A.setFromTriplets( triplets.begin(), triplets.end());
#endif
}


/// GMM to Cholmod_sparse interface
template<class MatrixT>
void eigen_to_cholmod( const MatrixT& _A, cholmod_sparse* &_AC, cholmod_common* _common, int _sparsity_type, bool _long_int)
{
  /* _sparsity_type
          * 0:  matrix is "unsymmetric": use both upper and lower triangular parts
          *     (the matrix may actually be symmetric in pattern and value, but
          *     both parts are explicitly stored and used).  May be square or
          *     rectangular.
          * >0: matrix is square and symmetric, use upper triangular part.
          *     Entries in the lower triangular part are ignored.
          * <0: matrix is square and symmetric, use lower triangular part.
          *     Entries in the upper triangular part are ignored. */

  int m = _A.innerSize();
  int n = _A.outerSize();

  // get upper or lower
  char uplo = 'c';
  if(_sparsity_type < 0) uplo = 'l';
  if(_sparsity_type > 0) uplo = 'u';


  if( _long_int) // long int version
  {
    std::vector<double> values;
    std::vector<UF_long> rowind;
    std::vector<UF_long> colptr;

    // get data of gmm matrix
    COMISO_EIGEN::get_ccs_symmetric_data( _A, uplo, values, rowind, colptr);

    // allocate cholmod matrix
    _AC = cholmod_l_allocate_sparse(m,n,values.size(),true,true,_sparsity_type,CHOLMOD_REAL, _common);

    // copy data to cholmod matrix
    for(UF_long i=0; i<(UF_long)values.size(); ++i)
    {
      ((double*) (_AC->x))[i] = values[i];
      ((UF_long*)(_AC->i))[i] = rowind[i];
    }

    for(UF_long i=0; i<(UF_long)colptr.size(); ++i)
      ((UF_long*)(_AC->p))[i] = colptr[i];
  }
  else // int version
  {
     std::vector<double> values;
     std::vector<int> rowind;
     std::vector<int> colptr;

     // get data of gmm matrix
     COMISO_EIGEN::get_ccs_symmetric_data( _A, uplo, values, rowind, colptr);

     // allocate cholmod matrix
     _AC = cholmod_allocate_sparse(m,n,values.size(),true,true,_sparsity_type,CHOLMOD_REAL, _common);

     // copy data to cholmod matrix
     for(unsigned int i=0; i<values.size(); ++i)
     {
       ((double*)(_AC->x))[i] = values[i];
       ((int*)   (_AC->i))[i] = rowind[i];
     }
     for(unsigned int i=0; i<colptr.size(); ++i)
       ((int*)(_AC->p))[i] = colptr[i];
  }

}
#endif

/*
/// Eigen to Cholmod_dense interface
template<class MatrixT>
void cholmod_to_eigen_dense( const cholmod_dense& _AC, MatrixT& _A)
{
  // initialize dimensions
  typedef typename MatrixT::Scalar Scalar;
  typedef Eigen::Triplet< Scalar > Triplet;
  size_t nzmax( _AC.nzmax);
  _A = MatrixT(_AC.nrow, _AC.ncol);
  std::vector< Triplet > triplets;
  triplets.reserve(nzmax);

  if(!_AC.packed)
  {
    std::cerr << "Warning: " << __FUNCTION__ << " does not support unpacked matrices yet!!!" << std::endl;
    return;
  }

  // Pointer to data
  double* X((double*)_AC.x);

  // complete matrix stored
  if(_AC.stype == 0)
  {
    // which integer type?
    if(_AC.itype == CHOLMOD_LONG)
    {
      UF_long* P((UF_long*)_AC.p);
      UF_long* I((UF_long*)_AC.i);

      for(UF_long i=0; i<(UF_long)_AC.ncol; ++i)
        for(UF_long j= P[i]; j< P[i+1]; ++j)
          //_A( I[j], i) += X[j]; // += really needed?
          triplets.push_back( Triplet( I[j], i, X[j]));
    }
    else
    {
      int* P((int*)_AC.p);
      int* I((int*)_AC.i);

      for(int i=0; i<(int)_AC.ncol; ++i)
        for(int j= P[i]; j< P[i+1]; ++j)
          triplets.push_back( Triplet( I[j], i, X[j]));
      //_A( I[j], i) += X[j];
    }

  }
  else // only upper or lower diagonal stored
  {
    // which integer type?
    if(_AC.itype == CHOLMOD_LONG)
    {
      UF_long* P((UF_long*)_AC.p);
      UF_long* I((UF_long*)_AC.i);

      for(UF_long i=0; i<(UF_long)_AC.ncol; ++i)
        for(UF_long j=P[i]; j<P[i+1]; ++j)
        {
          //_A(I[j], i) += X[j];
          triplets.push_back( Triplet( I[j], i, X[j]));

          // add up symmetric part
          if( I[j] != i)
            triplets.push_back( Triplet( i, I[j], X[j]));
          //_A(i,I[j]) += X[j];
        }
    }
    else
    {
      int* P((int*)_AC.p);
      int* I((int*)_AC.i);

      for(int i=0; i<(int)_AC.ncol; ++i)
        for(int j=P[i]; j<P[i+1]; ++j)
        {
          //_A(I[j], i) += X[j];
          triplets.push_back( Triplet( I[j], i, X[j]));

          // add up symmetric part
          if( I[j] != i)
            //  _A(i,I[j]) += X[j];
            triplets.push_back( Triplet( i, I[j], X[j]));
        }
    }
  }
  _A.setFromTriplets( triplets.begin(), triplets.end());
}


/// GMM to Cholmod_sparse interface
template<class MatrixT>
void eigen_to_cholmod_dense( const MatrixT& _A, cholmod_dense* &_AC, cholmod_common* _common, bool _long_int)
{
  int m = _A.innerSize();
  int n = _A.outerSize();

  // allocate cholmod matrix
  _AC = cholmod_l_allocate_sparse(m,n,values.size(),true,true,_sparsity_type,CHOLMOD_REAL, _common);
  _AC = cholmod_l_allocate_dense (m,n,n, xtype, cc) 

  if( _long_int) // long int version
  {
    std::vector<double> values;
    std::vector<UF_long> rowind;
    std::vector<UF_long> colptr;

    // get data of gmm matrix
    COMISO_EIGEN::get_ccs_symmetric_data( _A, uplo, values, rowind, colptr);

    // allocate cholmod matrix
    _AC = cholmod_l_allocate_sparse(m,n,values.size(),true,true,_sparsity_type,CHOLMOD_REAL, _common);

    // copy data to cholmod matrix
    for(UF_long i=0; i<(UF_long)values.size(); ++i)
    {
      ((double*) (_AC->x))[i] = values[i];
      ((UF_long*)(_AC->i))[i] = rowind[i];
    }

    for(UF_long i=0; i<(UF_long)colptr.size(); ++i)
      ((UF_long*)(_AC->p))[i] = colptr[i];
  }
  else // int version
  {
     std::vector<double> values;
     std::vector<int> rowind;
     std::vector<int> colptr;

     // get data of gmm matrix
     COMISO_EIGEN::get_ccs_symmetric_data( _A, uplo, values, rowind, colptr);

     // allocate cholmod matrix
     _AC = cholmod_allocate_sparse(m,n,values.size(),true,true,_sparsity_type,CHOLMOD_REAL, _common);

     // copy data to cholmod matrix
     for(unsigned int i=0; i<values.size(); ++i)
     {
       ((double*)(_AC->x))[i] = values[i];
       ((int*)   (_AC->i))[i] = rowind[i];
     }
     for(unsigned int i=0; i<colptr.size(); ++i)
       ((int*)(_AC->p))[i] = colptr[i];
  }

}*/

	/*
// convert a gmm col-sparse matrix into an eigen sparse matrix
template<class GMM_MatrixT, class EIGEN_MatrixT>
void gmm_to_eigen( const GMM_MatrixT& _G, EIGEN_MatrixT& _E)
{
#ifdef COMISO_EIGEN3_AVAILABLE
  typedef typename EIGEN_MatrixT::Scalar Scalar;

  typedef typename gmm::linalg_traits<GMM_MatrixT>::const_sub_col_type ColT;
  typedef typename gmm::linalg_traits<ColT>::const_iterator CIter;

  // build matrix triplets
  typedef Eigen::Triplet< Scalar > Triplet;
  std::vector< Triplet > triplets;
  triplets.reserve(gmm::nnz(_G));

  for(unsigned int i=0; i<gmm::mat_ncols(_G); ++i)
  {
     ColT col = mat_const_col( _G, i );

     CIter it  = gmm::vect_const_begin( col );
     CIter ite = gmm::vect_const_end( col );
     for ( ; it!=ite; ++it )
       triplets.push_back( Triplet( it.index(), i, *it));

  }

  // generate eigen matrix
  _E = EIGEN_MatrixT( gmm::mat_nrows(_G), gmm::mat_ncols(_G));
  _E.setFromTriplets( triplets.begin(), triplets.end());
#endif
}
*/

// convert a gmm col-sparse matrix into an eigen sparse matrix
template<class GMM_VectorT, class EIGEN_MatrixT>
void gmm_to_eigen( const gmm::col_matrix<GMM_VectorT>& _G, EIGEN_MatrixT& _E)
{
#ifdef COMISO_EIGEN3_AVAILABLE
  typedef typename EIGEN_MatrixT::Scalar Scalar;

  typedef typename gmm::col_matrix<GMM_VectorT> GMM_MatrixT;
  typedef typename gmm::linalg_traits<GMM_MatrixT>::const_sub_col_type ColT;
  typedef typename gmm::linalg_traits<ColT>::const_iterator CIter;

  // build matrix triplets
  typedef Eigen::Triplet< Scalar > Triplet;
  std::vector< Triplet > triplets;
  triplets.reserve(gmm::nnz(_G));

  for(unsigned int i=0; i<gmm::mat_ncols(_G); ++i)
  {
     ColT col = mat_const_col( _G, i );

     CIter it  = gmm::vect_const_begin( col );
     CIter ite = gmm::vect_const_end( col );
     for ( ; it!=ite; ++it )
       triplets.push_back( Triplet( it.index(), i, *it));

  }

  // generate eigen matrix
  _E = EIGEN_MatrixT( gmm::mat_nrows(_G), gmm::mat_ncols(_G));
  _E.setFromTriplets( triplets.begin(), triplets.end());
#endif
}


// convert a gmm row-sparse matrix into an eigen sparse matrix
template<class GMM_VectorT, class EIGEN_MatrixT>
void gmm_to_eigen( const gmm::row_matrix<GMM_VectorT>& _G, EIGEN_MatrixT& _E)
{
#ifdef COMISO_EIGEN3_AVAILABLE
  typedef typename EIGEN_MatrixT::Scalar Scalar;

  typedef typename gmm::row_matrix<GMM_VectorT> GMM_MatrixT;
  typedef typename gmm::linalg_traits<GMM_MatrixT>::const_sub_row_type RowT;
  typedef typename gmm::linalg_traits<RowT>::const_iterator CIter;

  // build matrix triplets
  typedef Eigen::Triplet< Scalar > Triplet;
  std::vector< Triplet > triplets;
  triplets.reserve(gmm::nnz(_G));

  for(unsigned int i=0; i<gmm::mat_nrows(_G); ++i)
  {
     RowT row = mat_const_row( _G, i );

     CIter it  = gmm::vect_const_begin( row );
     CIter ite = gmm::vect_const_end( row );
     for ( ; it!=ite; ++it )
       triplets.push_back( Triplet( i, it.index(), *it));

  }

  // generate eigen matrix
  _E = EIGEN_MatrixT( gmm::mat_nrows(_G), gmm::mat_ncols(_G));
  _E.setFromTriplets( triplets.begin(), triplets.end());
#endif
}

// convert a gmm col-sparse matrix into an eigen sparse matrix
template<class GMM_RealT, class EIGEN_MatrixT>
void gmm_to_eigen( const gmm::csc_matrix<GMM_RealT,0>& _G, EIGEN_MatrixT& _E)
{
#ifdef COMISO_EIGEN3_AVAILABLE
  typedef typename EIGEN_MatrixT::Scalar Scalar;

  typedef typename gmm::csc_matrix<GMM_RealT,0> GMM_MatrixT;
  typedef typename gmm::linalg_traits<GMM_MatrixT>::const_sub_col_type ColT;
  typedef typename gmm::linalg_traits<ColT>::const_iterator CIter;

  // build matrix triplets
  typedef Eigen::Triplet< Scalar > Triplet;
  std::vector< Triplet > triplets;
  triplets.reserve(gmm::nnz(_G));

  for(unsigned int i=0; i<gmm::mat_ncols(_G); ++i)
  {
     ColT col = mat_const_col( _G, i );

     CIter it  = gmm::vect_const_begin( col );
     CIter ite = gmm::vect_const_end( col );
     for ( ; it!=ite; ++it )
       triplets.push_back( Triplet( it.index(), i, *it));

  }

  // generate eigen matrix
  _E = EIGEN_MatrixT( gmm::mat_nrows(_G), gmm::mat_ncols(_G));
  _E.setFromTriplets( triplets.begin(), triplets.end());
#endif
}

783
784
785
786
787
788
789
790
791
792
793
794
795
796
// this explicit instantiation does not match the call signature due to the 
// partial specializations above
//template void gmm_to_eigen(const gmm::csc_matrix<double>&,
//  Eigen::SparseMatrix<double>& );

// hence make a partial specialization that matches the main template signature
void gmm_to_eigen( const gmm::csc_matrix<double>& _G, 
  Eigen::SparseMatrix<double>& _E)
{
  // and redirect to the above partial specialization 
  gmm_to_eigen< double, Eigen::SparseMatrix<double> > (_G, _E);
}


797
//=============================================================================
798
} // namespace COMISO_EIGEN
799
800
801
802
803
//=============================================================================

//=============================================================================
#endif // COMISO_EIGEN3_AVAILABLE
//=============================================================================