SparseSelfAdjointView.h 18 KB
Newer Older
LM's avatar
LM committed
1 2 3 4 5
// This file is part of Eigen, a lightweight C++ template library
// for linear algebra.
//
// Copyright (C) 2009 Gael Guennebaud <gael.guennebaud@inria.fr>
//
Don Gagne's avatar
Don Gagne committed
6 7 8
// This Source Code Form is subject to the terms of the Mozilla
// Public License v. 2.0. If a copy of the MPL was not distributed
// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
LM's avatar
LM committed
9 10 11 12

#ifndef EIGEN_SPARSE_SELFADJOINTVIEW_H
#define EIGEN_SPARSE_SELFADJOINTVIEW_H

Don Gagne's avatar
Don Gagne committed
13 14 15 16
namespace Eigen { 

/** \ingroup SparseCore_Module
  * \class SparseSelfAdjointView
LM's avatar
LM committed
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 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71
  *
  * \brief Pseudo expression to manipulate a triangular sparse matrix as a selfadjoint matrix.
  *
  * \param MatrixType the type of the dense matrix storing the coefficients
  * \param UpLo can be either \c #Lower or \c #Upper
  *
  * This class is an expression of a sefladjoint matrix from a triangular part of a matrix
  * with given dense storage of the coefficients. It is the return type of MatrixBase::selfadjointView()
  * and most of the time this is the only way that it is used.
  *
  * \sa SparseMatrixBase::selfadjointView()
  */
template<typename Lhs, typename Rhs, int UpLo>
class SparseSelfAdjointTimeDenseProduct;

template<typename Lhs, typename Rhs, int UpLo>
class DenseTimeSparseSelfAdjointProduct;

namespace internal {
  
template<typename MatrixType, unsigned int UpLo>
struct traits<SparseSelfAdjointView<MatrixType,UpLo> > : traits<MatrixType> {
};

template<int SrcUpLo,int DstUpLo,typename MatrixType,int DestOrder>
void permute_symm_to_symm(const MatrixType& mat, SparseMatrix<typename MatrixType::Scalar,DestOrder,typename MatrixType::Index>& _dest, const typename MatrixType::Index* perm = 0);

template<int UpLo,typename MatrixType,int DestOrder>
void permute_symm_to_fullsymm(const MatrixType& mat, SparseMatrix<typename MatrixType::Scalar,DestOrder,typename MatrixType::Index>& _dest, const typename MatrixType::Index* perm = 0);

}

template<typename MatrixType, unsigned int UpLo> class SparseSelfAdjointView
  : public EigenBase<SparseSelfAdjointView<MatrixType,UpLo> >
{
  public:

    typedef typename MatrixType::Scalar Scalar;
    typedef typename MatrixType::Index Index;
    typedef Matrix<Index,Dynamic,1> VectorI;
    typedef typename MatrixType::Nested MatrixTypeNested;
    typedef typename internal::remove_all<MatrixTypeNested>::type _MatrixTypeNested;

    inline SparseSelfAdjointView(const MatrixType& matrix) : m_matrix(matrix)
    {
      eigen_assert(rows()==cols() && "SelfAdjointView is only for squared matrices");
    }

    inline Index rows() const { return m_matrix.rows(); }
    inline Index cols() const { return m_matrix.cols(); }

    /** \internal \returns a reference to the nested matrix */
    const _MatrixTypeNested& matrix() const { return m_matrix; }
    _MatrixTypeNested& matrix() { return m_matrix.const_cast_derived(); }

Don Gagne's avatar
Don Gagne committed
72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95
    /** \returns an expression of the matrix product between a sparse self-adjoint matrix \c *this and a sparse matrix \a rhs.
      *
      * Note that there is no algorithmic advantage of performing such a product compared to a general sparse-sparse matrix product.
      * Indeed, the SparseSelfadjointView operand is first copied into a temporary SparseMatrix before computing the product.
      */
    template<typename OtherDerived>
    SparseSparseProduct<typename OtherDerived::PlainObject, OtherDerived>
    operator*(const SparseMatrixBase<OtherDerived>& rhs) const
    {
      return SparseSparseProduct<typename OtherDerived::PlainObject, OtherDerived>(*this, rhs.derived());
    }

    /** \returns an expression of the matrix product between a sparse matrix \a lhs and a sparse self-adjoint matrix \a rhs.
      *
      * Note that there is no algorithmic advantage of performing such a product compared to a general sparse-sparse matrix product.
      * Indeed, the SparseSelfadjointView operand is first copied into a temporary SparseMatrix before computing the product.
      */
    template<typename OtherDerived> friend
    SparseSparseProduct<OtherDerived, typename OtherDerived::PlainObject >
    operator*(const SparseMatrixBase<OtherDerived>& lhs, const SparseSelfAdjointView& rhs)
    {
      return SparseSparseProduct<OtherDerived, typename OtherDerived::PlainObject>(lhs.derived(), rhs);
    }
    
LM's avatar
LM committed
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
    /** Efficient sparse self-adjoint matrix times dense vector/matrix product */
    template<typename OtherDerived>
    SparseSelfAdjointTimeDenseProduct<MatrixType,OtherDerived,UpLo>
    operator*(const MatrixBase<OtherDerived>& rhs) const
    {
      return SparseSelfAdjointTimeDenseProduct<MatrixType,OtherDerived,UpLo>(m_matrix, rhs.derived());
    }

    /** Efficient dense vector/matrix times sparse self-adjoint matrix product */
    template<typename OtherDerived> friend
    DenseTimeSparseSelfAdjointProduct<OtherDerived,MatrixType,UpLo>
    operator*(const MatrixBase<OtherDerived>& lhs, const SparseSelfAdjointView& rhs)
    {
      return DenseTimeSparseSelfAdjointProduct<OtherDerived,_MatrixTypeNested,UpLo>(lhs.derived(), rhs.m_matrix);
    }

    /** Perform a symmetric rank K update of the selfadjoint matrix \c *this:
      * \f$ this = this + \alpha ( u u^* ) \f$ where \a u is a vector or matrix.
      *
      * \returns a reference to \c *this
      *
      * To perform \f$ this = this + \alpha ( u^* u ) \f$ you can simply
      * call this function with u.adjoint().
      */
    template<typename DerivedU>
Don Gagne's avatar
Don Gagne committed
121
    SparseSelfAdjointView& rankUpdate(const SparseMatrixBase<DerivedU>& u, const Scalar& alpha = Scalar(1));
LM's avatar
LM committed
122 123
    
    /** \internal triggered by sparse_matrix = SparseSelfadjointView; */
Don Gagne's avatar
Don Gagne committed
124
    template<typename DestScalar,int StorageOrder> void evalTo(SparseMatrix<DestScalar,StorageOrder,Index>& _dest) const
LM's avatar
LM committed
125 126 127 128
    {
      internal::permute_symm_to_fullsymm<UpLo>(m_matrix, _dest);
    }
    
Don Gagne's avatar
Don Gagne committed
129
    template<typename DestScalar> void evalTo(DynamicSparseMatrix<DestScalar,ColMajor,Index>& _dest) const
LM's avatar
LM committed
130 131
    {
      // TODO directly evaluate into _dest;
Don Gagne's avatar
Don Gagne committed
132
      SparseMatrix<DestScalar,ColMajor,Index> tmp(_dest.rows(),_dest.cols());
LM's avatar
LM committed
133 134 135 136
      internal::permute_symm_to_fullsymm<UpLo>(m_matrix, tmp);
      _dest = tmp;
    }
    
Don Gagne's avatar
Don Gagne committed
137 138
    /** \returns an expression of P H P^-1 */
    SparseSymmetricPermutationProduct<_MatrixTypeNested,UpLo> twistedBy(const PermutationMatrix<Dynamic,Dynamic,Index>& perm) const
LM's avatar
LM committed
139 140 141 142 143 144 145 146 147 148
    {
      return SparseSymmetricPermutationProduct<_MatrixTypeNested,UpLo>(m_matrix, perm);
    }
    
    template<typename SrcMatrixType,int SrcUpLo>
    SparseSelfAdjointView& operator=(const SparseSymmetricPermutationProduct<SrcMatrixType,SrcUpLo>& permutedMatrix)
    {
      permutedMatrix.evalTo(*this);
      return *this;
    }
Don Gagne's avatar
Don Gagne committed
149 150 151 152 153 154 155 156 157 158 159 160 161 162


    SparseSelfAdjointView& operator=(const SparseSelfAdjointView& src)
    {
      PermutationMatrix<Dynamic> pnull;
      return *this = src.twistedBy(pnull);
    }

    template<typename SrcMatrixType,unsigned int SrcUpLo>
    SparseSelfAdjointView& operator=(const SparseSelfAdjointView<SrcMatrixType,SrcUpLo>& src)
    {
      PermutationMatrix<Dynamic> pnull;
      return *this = src.twistedBy(pnull);
    }
LM's avatar
LM committed
163 164 165 166 167 168 169
    

    // const SparseLLT<PlainObject, UpLo> llt() const;
    // const SparseLDLT<PlainObject, UpLo> ldlt() const;

  protected:

Don Gagne's avatar
Don Gagne committed
170
    typename MatrixType::Nested m_matrix;
LM's avatar
LM committed
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
    mutable VectorI m_countPerRow;
    mutable VectorI m_countPerCol;
};

/***************************************************************************
* Implementation of SparseMatrixBase methods
***************************************************************************/

template<typename Derived>
template<unsigned int UpLo>
const SparseSelfAdjointView<Derived, UpLo> SparseMatrixBase<Derived>::selfadjointView() const
{
  return derived();
}

template<typename Derived>
template<unsigned int UpLo>
SparseSelfAdjointView<Derived, UpLo> SparseMatrixBase<Derived>::selfadjointView()
{
  return derived();
}

/***************************************************************************
* Implementation of SparseSelfAdjointView methods
***************************************************************************/

template<typename MatrixType, unsigned int UpLo>
template<typename DerivedU>
SparseSelfAdjointView<MatrixType,UpLo>&
Don Gagne's avatar
Don Gagne committed
200
SparseSelfAdjointView<MatrixType,UpLo>::rankUpdate(const SparseMatrixBase<DerivedU>& u, const Scalar& alpha)
LM's avatar
LM committed
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
{
  SparseMatrix<Scalar,MatrixType::Flags&RowMajorBit?RowMajor:ColMajor> tmp = u * u.adjoint();
  if(alpha==Scalar(0))
    m_matrix.const_cast_derived() = tmp.template triangularView<UpLo>();
  else
    m_matrix.const_cast_derived() += alpha * tmp.template triangularView<UpLo>();

  return *this;
}

/***************************************************************************
* Implementation of sparse self-adjoint time dense matrix
***************************************************************************/

namespace internal {
template<typename Lhs, typename Rhs, int UpLo>
struct traits<SparseSelfAdjointTimeDenseProduct<Lhs,Rhs,UpLo> >
 : traits<ProductBase<SparseSelfAdjointTimeDenseProduct<Lhs,Rhs,UpLo>, Lhs, Rhs> >
{
  typedef Dense StorageKind;
};
}

template<typename Lhs, typename Rhs, int UpLo>
class SparseSelfAdjointTimeDenseProduct
  : public ProductBase<SparseSelfAdjointTimeDenseProduct<Lhs,Rhs,UpLo>, Lhs, Rhs>
{
  public:
    EIGEN_PRODUCT_PUBLIC_INTERFACE(SparseSelfAdjointTimeDenseProduct)

    SparseSelfAdjointTimeDenseProduct(const Lhs& lhs, const Rhs& rhs) : Base(lhs,rhs)
    {}

Don Gagne's avatar
Don Gagne committed
234
    template<typename Dest> void scaleAndAddTo(Dest& dest, const Scalar& alpha) const
LM's avatar
LM committed
235
    {
Don Gagne's avatar
Don Gagne committed
236
      EIGEN_ONLY_USED_FOR_DEBUG(alpha);
LM's avatar
LM committed
237 238 239 240 241 242 243 244 245 246 247 248 249 250 251
      // TODO use alpha
      eigen_assert(alpha==Scalar(1) && "alpha != 1 is not implemented yet, sorry");
      typedef typename internal::remove_all<Lhs>::type _Lhs;
      typedef typename _Lhs::InnerIterator LhsInnerIterator;
      enum {
        LhsIsRowMajor = (_Lhs::Flags&RowMajorBit)==RowMajorBit,
        ProcessFirstHalf =
                 ((UpLo&(Upper|Lower))==(Upper|Lower))
              || ( (UpLo&Upper) && !LhsIsRowMajor)
              || ( (UpLo&Lower) && LhsIsRowMajor),
        ProcessSecondHalf = !ProcessFirstHalf
      };
      for (Index j=0; j<m_lhs.outerSize(); ++j)
      {
        LhsInnerIterator i(m_lhs,j);
Don Gagne's avatar
Don Gagne committed
252
        if (ProcessSecondHalf)
LM's avatar
LM committed
253
        {
Don Gagne's avatar
Don Gagne committed
254 255 256 257 258 259
          while (i && i.index()<j) ++i;
          if(i && i.index()==j)
          {
            dest.row(j) += i.value() * m_rhs.row(j);
            ++i;
          }
LM's avatar
LM committed
260 261 262 263 264 265 266
        }
        for(; (ProcessFirstHalf ? i && i.index() < j : i) ; ++i)
        {
          Index a = LhsIsRowMajor ? j : i.index();
          Index b = LhsIsRowMajor ? i.index() : j;
          typename Lhs::Scalar v = i.value();
          dest.row(a) += (v) * m_rhs.row(b);
Don Gagne's avatar
Don Gagne committed
267
          dest.row(b) += numext::conj(v) * m_rhs.row(a);
LM's avatar
LM committed
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
        }
        if (ProcessFirstHalf && i && (i.index()==j))
          dest.row(j) += i.value() * m_rhs.row(j);
      }
    }

  private:
    SparseSelfAdjointTimeDenseProduct& operator=(const SparseSelfAdjointTimeDenseProduct&);
};

namespace internal {
template<typename Lhs, typename Rhs, int UpLo>
struct traits<DenseTimeSparseSelfAdjointProduct<Lhs,Rhs,UpLo> >
 : traits<ProductBase<DenseTimeSparseSelfAdjointProduct<Lhs,Rhs,UpLo>, Lhs, Rhs> >
{};
}

template<typename Lhs, typename Rhs, int UpLo>
class DenseTimeSparseSelfAdjointProduct
  : public ProductBase<DenseTimeSparseSelfAdjointProduct<Lhs,Rhs,UpLo>, Lhs, Rhs>
{
  public:
    EIGEN_PRODUCT_PUBLIC_INTERFACE(DenseTimeSparseSelfAdjointProduct)

    DenseTimeSparseSelfAdjointProduct(const Lhs& lhs, const Rhs& rhs) : Base(lhs,rhs)
    {}

Don Gagne's avatar
Don Gagne committed
295
    template<typename Dest> void scaleAndAddTo(Dest& /*dest*/, const Scalar& /*alpha*/) const
LM's avatar
LM committed
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
    {
      // TODO
    }

  private:
    DenseTimeSparseSelfAdjointProduct& operator=(const DenseTimeSparseSelfAdjointProduct&);
};

/***************************************************************************
* Implementation of symmetric copies and permutations
***************************************************************************/
namespace internal {
  
template<typename MatrixType, int UpLo>
struct traits<SparseSymmetricPermutationProduct<MatrixType,UpLo> > : traits<MatrixType> {
};

template<int UpLo,typename MatrixType,int DestOrder>
void permute_symm_to_fullsymm(const MatrixType& mat, SparseMatrix<typename MatrixType::Scalar,DestOrder,typename MatrixType::Index>& _dest, const typename MatrixType::Index* perm)
{
  typedef typename MatrixType::Index Index;
  typedef typename MatrixType::Scalar Scalar;
  typedef SparseMatrix<Scalar,DestOrder,Index> Dest;
  typedef Matrix<Index,Dynamic,1> VectorI;
  
  Dest& dest(_dest.derived());
  enum {
    StorageOrderMatch = int(Dest::IsRowMajor) == int(MatrixType::IsRowMajor)
  };
Don Gagne's avatar
Don Gagne committed
325
  
LM's avatar
LM committed
326 327 328 329 330 331 332 333 334 335 336
  Index size = mat.rows();
  VectorI count;
  count.resize(size);
  count.setZero();
  dest.resize(size,size);
  for(Index j = 0; j<size; ++j)
  {
    Index jp = perm ? perm[j] : j;
    for(typename MatrixType::InnerIterator it(mat,j); it; ++it)
    {
      Index i = it.index();
Don Gagne's avatar
Don Gagne committed
337 338
      Index r = it.row();
      Index c = it.col();
LM's avatar
LM committed
339
      Index ip = perm ? perm[i] : i;
Don Gagne's avatar
Don Gagne committed
340 341 342
      if(UpLo==(Upper|Lower))
        count[StorageOrderMatch ? jp : ip]++;
      else if(r==c)
LM's avatar
LM committed
343
        count[ip]++;
Don Gagne's avatar
Don Gagne committed
344
      else if(( UpLo==Lower && r>c) || ( UpLo==Upper && r<c))
LM's avatar
LM committed
345 346 347 348 349 350 351 352 353
      {
        count[ip]++;
        count[jp]++;
      }
    }
  }
  Index nnz = count.sum();
  
  // reserve space
Don Gagne's avatar
Don Gagne committed
354 355
  dest.resizeNonZeros(nnz);
  dest.outerIndexPtr()[0] = 0;
LM's avatar
LM committed
356
  for(Index j=0; j<size; ++j)
Don Gagne's avatar
Don Gagne committed
357
    dest.outerIndexPtr()[j+1] = dest.outerIndexPtr()[j] + count[j];
LM's avatar
LM committed
358
  for(Index j=0; j<size; ++j)
Don Gagne's avatar
Don Gagne committed
359
    count[j] = dest.outerIndexPtr()[j];
LM's avatar
LM committed
360 361 362 363 364 365 366
  
  // copy data
  for(Index j = 0; j<size; ++j)
  {
    for(typename MatrixType::InnerIterator it(mat,j); it; ++it)
    {
      Index i = it.index();
Don Gagne's avatar
Don Gagne committed
367 368 369 370
      Index r = it.row();
      Index c = it.col();
      
      Index jp = perm ? perm[j] : j;
LM's avatar
LM committed
371
      Index ip = perm ? perm[i] : i;
Don Gagne's avatar
Don Gagne committed
372 373
      
      if(UpLo==(Upper|Lower))
LM's avatar
LM committed
374
      {
Don Gagne's avatar
Don Gagne committed
375 376 377
        Index k = count[StorageOrderMatch ? jp : ip]++;
        dest.innerIndexPtr()[k] = StorageOrderMatch ? ip : jp;
        dest.valuePtr()[k] = it.value();
LM's avatar
LM committed
378
      }
Don Gagne's avatar
Don Gagne committed
379
      else if(r==c)
LM's avatar
LM committed
380
      {
Don Gagne's avatar
Don Gagne committed
381 382 383 384 385 386 387 388 389 390 391
        Index k = count[ip]++;
        dest.innerIndexPtr()[k] = ip;
        dest.valuePtr()[k] = it.value();
      }
      else if(( (UpLo&Lower)==Lower && r>c) || ( (UpLo&Upper)==Upper && r<c))
      {
        if(!StorageOrderMatch)
          std::swap(ip,jp);
        Index k = count[jp]++;
        dest.innerIndexPtr()[k] = ip;
        dest.valuePtr()[k] = it.value();
LM's avatar
LM committed
392
        k = count[ip]++;
Don Gagne's avatar
Don Gagne committed
393 394
        dest.innerIndexPtr()[k] = jp;
        dest.valuePtr()[k] = numext::conj(it.value());
LM's avatar
LM committed
395 396 397 398 399
      }
    }
  }
}

Don Gagne's avatar
Don Gagne committed
400 401
template<int _SrcUpLo,int _DstUpLo,typename MatrixType,int DstOrder>
void permute_symm_to_symm(const MatrixType& mat, SparseMatrix<typename MatrixType::Scalar,DstOrder,typename MatrixType::Index>& _dest, const typename MatrixType::Index* perm)
LM's avatar
LM committed
402 403 404
{
  typedef typename MatrixType::Index Index;
  typedef typename MatrixType::Scalar Scalar;
Don Gagne's avatar
Don Gagne committed
405
  SparseMatrix<Scalar,DstOrder,Index>& dest(_dest.derived());
LM's avatar
LM committed
406
  typedef Matrix<Index,Dynamic,1> VectorI;
Don Gagne's avatar
Don Gagne committed
407 408 409 410 411 412
  enum {
    SrcOrder = MatrixType::IsRowMajor ? RowMajor : ColMajor,
    StorageOrderMatch = int(SrcOrder) == int(DstOrder),
    DstUpLo = DstOrder==RowMajor ? (_DstUpLo==Upper ? Lower : Upper) : _DstUpLo,
    SrcUpLo = SrcOrder==RowMajor ? (_SrcUpLo==Upper ? Lower : Upper) : _SrcUpLo
  };
LM's avatar
LM committed
413 414 415 416 417 418 419 420 421 422 423
  
  Index size = mat.rows();
  VectorI count(size);
  count.setZero();
  dest.resize(size,size);
  for(Index j = 0; j<size; ++j)
  {
    Index jp = perm ? perm[j] : j;
    for(typename MatrixType::InnerIterator it(mat,j); it; ++it)
    {
      Index i = it.index();
Don Gagne's avatar
Don Gagne committed
424
      if((int(SrcUpLo)==int(Lower) && i<j) || (int(SrcUpLo)==int(Upper) && i>j))
LM's avatar
LM committed
425 426 427
        continue;
                  
      Index ip = perm ? perm[i] : i;
Don Gagne's avatar
Don Gagne committed
428
      count[int(DstUpLo)==int(Lower) ? (std::min)(ip,jp) : (std::max)(ip,jp)]++;
LM's avatar
LM committed
429 430
    }
  }
Don Gagne's avatar
Don Gagne committed
431
  dest.outerIndexPtr()[0] = 0;
LM's avatar
LM committed
432
  for(Index j=0; j<size; ++j)
Don Gagne's avatar
Don Gagne committed
433 434
    dest.outerIndexPtr()[j+1] = dest.outerIndexPtr()[j] + count[j];
  dest.resizeNonZeros(dest.outerIndexPtr()[size]);
LM's avatar
LM committed
435
  for(Index j=0; j<size; ++j)
Don Gagne's avatar
Don Gagne committed
436
    count[j] = dest.outerIndexPtr()[j];
LM's avatar
LM committed
437 438 439
  
  for(Index j = 0; j<size; ++j)
  {
Don Gagne's avatar
Don Gagne committed
440
    
LM's avatar
LM committed
441 442 443
    for(typename MatrixType::InnerIterator it(mat,j); it; ++it)
    {
      Index i = it.index();
Don Gagne's avatar
Don Gagne committed
444
      if((int(SrcUpLo)==int(Lower) && i<j) || (int(SrcUpLo)==int(Upper) && i>j))
LM's avatar
LM committed
445 446
        continue;
                  
Don Gagne's avatar
Don Gagne committed
447
      Index jp = perm ? perm[j] : j;
LM's avatar
LM committed
448 449
      Index ip = perm? perm[i] : i;
      
Don Gagne's avatar
Don Gagne committed
450 451 452 453 454 455
      Index k = count[int(DstUpLo)==int(Lower) ? (std::min)(ip,jp) : (std::max)(ip,jp)]++;
      dest.innerIndexPtr()[k] = int(DstUpLo)==int(Lower) ? (std::max)(ip,jp) : (std::min)(ip,jp);
      
      if(!StorageOrderMatch) std::swap(ip,jp);
      if( ((int(DstUpLo)==int(Lower) && ip<jp) || (int(DstUpLo)==int(Upper) && ip>jp)))
        dest.valuePtr()[k] = numext::conj(it.value());
LM's avatar
LM committed
456
      else
Don Gagne's avatar
Don Gagne committed
457
        dest.valuePtr()[k] = it.value();
LM's avatar
LM committed
458 459 460 461 462 463 464 465 466 467 468 469 470
    }
  }
}

}

template<typename MatrixType,int UpLo>
class SparseSymmetricPermutationProduct
  : public EigenBase<SparseSymmetricPermutationProduct<MatrixType,UpLo> >
{
  public:
    typedef typename MatrixType::Scalar Scalar;
    typedef typename MatrixType::Index Index;
Don Gagne's avatar
Don Gagne committed
471 472 473
  protected:
    typedef PermutationMatrix<Dynamic,Dynamic,Index> Perm;
  public:
LM's avatar
LM committed
474 475 476 477 478 479 480 481 482 483 484
    typedef Matrix<Index,Dynamic,1> VectorI;
    typedef typename MatrixType::Nested MatrixTypeNested;
    typedef typename internal::remove_all<MatrixTypeNested>::type _MatrixTypeNested;
    
    SparseSymmetricPermutationProduct(const MatrixType& mat, const Perm& perm)
      : m_matrix(mat), m_perm(perm)
    {}
    
    inline Index rows() const { return m_matrix.rows(); }
    inline Index cols() const { return m_matrix.cols(); }
    
Don Gagne's avatar
Don Gagne committed
485 486
    template<typename DestScalar, int Options, typename DstIndex>
    void evalTo(SparseMatrix<DestScalar,Options,DstIndex>& _dest) const
LM's avatar
LM committed
487
    {
Don Gagne's avatar
Don Gagne committed
488 489 490 491
//       internal::permute_symm_to_fullsymm<UpLo>(m_matrix,_dest,m_perm.indices().data());
      SparseMatrix<DestScalar,(Options&RowMajor)==RowMajor ? ColMajor : RowMajor, DstIndex> tmp;
      internal::permute_symm_to_fullsymm<UpLo>(m_matrix,tmp,m_perm.indices().data());
      _dest = tmp;
LM's avatar
LM committed
492 493 494 495 496 497 498 499
    }
    
    template<typename DestType,unsigned int DestUpLo> void evalTo(SparseSelfAdjointView<DestType,DestUpLo>& dest) const
    {
      internal::permute_symm_to_symm<UpLo,DestUpLo>(m_matrix,dest.matrix(),m_perm.indices().data());
    }
    
  protected:
Don Gagne's avatar
Don Gagne committed
500
    MatrixTypeNested m_matrix;
LM's avatar
LM committed
501 502 503 504
    const Perm& m_perm;

};

Don Gagne's avatar
Don Gagne committed
505 506
} // end namespace Eigen

LM's avatar
LM committed
507
#endif // EIGEN_SPARSE_SELFADJOINTVIEW_H