TriangularSolverMatrix.h 13.6 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_TRIANGULAR_SOLVER_MATRIX_H
#define EIGEN_TRIANGULAR_SOLVER_MATRIX_H

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

LM's avatar
LM committed
15 16 17 18 19 20
namespace internal {

// if the rhs is row major, let's transpose the product
template <typename Scalar, typename Index, int Side, int Mode, bool Conjugate, int TriStorageOrder>
struct triangular_solve_matrix<Scalar,Index,Side,Mode,Conjugate,TriStorageOrder,RowMajor>
{
Don Gagne's avatar
Don Gagne committed
21
  static void run(
LM's avatar
LM committed
22 23
    Index size, Index cols,
    const Scalar*  tri, Index triStride,
Don Gagne's avatar
Don Gagne committed
24 25
    Scalar* _other, Index otherStride,
    level3_blocking<Scalar,Scalar>& blocking)
LM's avatar
LM committed
26 27 28 29 30 31
  {
    triangular_solve_matrix<
      Scalar, Index, Side==OnTheLeft?OnTheRight:OnTheLeft,
      (Mode&UnitDiag) | ((Mode&Upper) ? Lower : Upper),
      NumTraits<Scalar>::IsComplex && Conjugate,
      TriStorageOrder==RowMajor ? ColMajor : RowMajor, ColMajor>
Don Gagne's avatar
Don Gagne committed
32
      ::run(size, cols, tri, triStride, _other, otherStride, blocking);
LM's avatar
LM committed
33 34 35 36 37 38 39 40 41 42 43
  }
};

/* Optimized triangular solver with multiple right hand side and the triangular matrix on the left
 */
template <typename Scalar, typename Index, int Mode, bool Conjugate, int TriStorageOrder>
struct triangular_solve_matrix<Scalar,Index,OnTheLeft,Mode,Conjugate,TriStorageOrder,ColMajor>
{
  static EIGEN_DONT_INLINE void run(
    Index size, Index otherSize,
    const Scalar* _tri, Index triStride,
Don Gagne's avatar
Don Gagne committed
44 45 46 47 48 49 50 51 52
    Scalar* _other, Index otherStride,
    level3_blocking<Scalar,Scalar>& blocking);
};
template <typename Scalar, typename Index, int Mode, bool Conjugate, int TriStorageOrder>
EIGEN_DONT_INLINE void triangular_solve_matrix<Scalar,Index,OnTheLeft,Mode,Conjugate,TriStorageOrder,ColMajor>::run(
    Index size, Index otherSize,
    const Scalar* _tri, Index triStride,
    Scalar* _other, Index otherStride,
    level3_blocking<Scalar,Scalar>& blocking)
LM's avatar
LM committed
53 54 55 56 57 58 59 60 61 62 63
  {
    Index cols = otherSize;
    const_blas_data_mapper<Scalar, Index, TriStorageOrder> tri(_tri,triStride);
    blas_data_mapper<Scalar, Index, ColMajor> other(_other,otherStride);

    typedef gebp_traits<Scalar,Scalar> Traits;
    enum {
      SmallPanelWidth   = EIGEN_PLAIN_ENUM_MAX(Traits::mr,Traits::nr),
      IsLower = (Mode&Lower) == Lower
    };

Don Gagne's avatar
Don Gagne committed
64 65
    Index kc = blocking.kc();                   // cache block size along the K direction
    Index mc = (std::min)(size,blocking.mc());  // cache block size along the M direction
LM's avatar
LM committed
66

Don Gagne's avatar
Don Gagne committed
67 68
    std::size_t sizeA = kc*mc;
    std::size_t sizeB = kc*cols;
LM's avatar
LM committed
69
    std::size_t sizeW = kc*Traits::WorkSpaceFactor;
Don Gagne's avatar
Don Gagne committed
70 71 72 73

    ei_declare_aligned_stack_constructed_variable(Scalar, blockA, sizeA, blocking.blockA());
    ei_declare_aligned_stack_constructed_variable(Scalar, blockB, sizeB, blocking.blockB());
    ei_declare_aligned_stack_constructed_variable(Scalar, blockW, sizeW, blocking.blockW());
LM's avatar
LM committed
74 75 76 77 78 79

    conj_if<Conjugate> conj;
    gebp_kernel<Scalar, Scalar, Index, Traits::mr, Traits::nr, Conjugate, false> gebp_kernel;
    gemm_pack_lhs<Scalar, Index, Traits::mr, Traits::LhsProgress, TriStorageOrder> pack_lhs;
    gemm_pack_rhs<Scalar, Index, Traits::nr, ColMajor, false, true> pack_rhs;

Don Gagne's avatar
Don Gagne committed
80 81 82 83 84 85 86
    // the goal here is to subdivise the Rhs panels such that we keep some cache
    // coherence when accessing the rhs elements
    std::ptrdiff_t l1, l2;
    manage_caching_sizes(GetAction, &l1, &l2);
    Index subcols = cols>0 ? l2/(4 * sizeof(Scalar) * otherStride) : 0;
    subcols = std::max<Index>((subcols/Traits::nr)*Traits::nr, Traits::nr);

LM's avatar
LM committed
87 88 89 90
    for(Index k2=IsLower ? 0 : size;
        IsLower ? k2<size : k2>0;
        IsLower ? k2+=kc : k2-=kc)
    {
Don Gagne's avatar
Don Gagne committed
91
      const Index actual_kc = (std::min)(IsLower ? size-k2 : k2, kc);
LM's avatar
LM committed
92 93 94 95 96 97

      // We have selected and packed a big horizontal panel R1 of rhs. Let B be the packed copy of this panel,
      // and R2 the remaining part of rhs. The corresponding vertical panel of lhs is split into
      // A11 (the triangular part) and A21 the remaining rectangular part.
      // Then the high level algorithm is:
      //  - B = R1                    => general block copy (done during the next step)
Don Gagne's avatar
Don Gagne committed
98
      //  - R1 = A11^-1 B             => tricky part
LM's avatar
LM committed
99
      //  - update B from the new R1  => actually this has to be performed continuously during the above step
Don Gagne's avatar
Don Gagne committed
100
      //  - R2 -= A21 * B             => GEPP
LM's avatar
LM committed
101

Don Gagne's avatar
Don Gagne committed
102 103 104 105 106
      // The tricky part: compute R1 = A11^-1 B while updating B from R1
      // The idea is to split A11 into multiple small vertical panels.
      // Each panel can be split into a small triangular part T1k which is processed without optimization,
      // and the remaining small part T2k which is processed using gebp with appropriate block strides
      for(Index j2=0; j2<cols; j2+=subcols)
LM's avatar
LM committed
107
      {
Don Gagne's avatar
Don Gagne committed
108 109
        Index actual_cols = (std::min)(cols-j2,subcols);
        // for each small vertical panels [T1k^T, T2k^T]^T of lhs
LM's avatar
LM committed
110 111 112 113 114 115 116 117 118 119 120 121
        for (Index k1=0; k1<actual_kc; k1+=SmallPanelWidth)
        {
          Index actualPanelWidth = std::min<Index>(actual_kc-k1, SmallPanelWidth);
          // tr solve
          for (Index k=0; k<actualPanelWidth; ++k)
          {
            // TODO write a small kernel handling this (can be shared with trsv)
            Index i  = IsLower ? k2+k1+k : k2-k1-k-1;
            Index s  = IsLower ? k2+k1 : i+1;
            Index rs = actualPanelWidth - k - 1; // remaining size

            Scalar a = (Mode & UnitDiag) ? Scalar(1) : Scalar(1)/conj(tri(i,i));
Don Gagne's avatar
Don Gagne committed
122
            for (Index j=j2; j<j2+actual_cols; ++j)
LM's avatar
LM committed
123 124 125
            {
              if (TriStorageOrder==RowMajor)
              {
Don Gagne's avatar
Don Gagne committed
126
                Scalar b(0);
LM's avatar
LM committed
127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150
                const Scalar* l = &tri(i,s);
                Scalar* r = &other(s,j);
                for (Index i3=0; i3<k; ++i3)
                  b += conj(l[i3]) * r[i3];

                other(i,j) = (other(i,j) - b)*a;
              }
              else
              {
                Index s = IsLower ? i+1 : i-rs;
                Scalar b = (other(i,j) *= a);
                Scalar* r = &other(s,j);
                const Scalar* l = &tri(s,i);
                for (Index i3=0;i3<rs;++i3)
                  r[i3] -= b * conj(l[i3]);
              }
            }
          }

          Index lengthTarget = actual_kc-k1-actualPanelWidth;
          Index startBlock   = IsLower ? k2+k1 : k2-k1-actualPanelWidth;
          Index blockBOffset = IsLower ? k1 : lengthTarget;

          // update the respective rows of B from other
Don Gagne's avatar
Don Gagne committed
151
          pack_rhs(blockB+actual_kc*j2, &other(startBlock,j2), otherStride, actualPanelWidth, actual_cols, actual_kc, blockBOffset);
LM's avatar
LM committed
152 153 154 155 156 157 158 159

          // GEBP
          if (lengthTarget>0)
          {
            Index startTarget  = IsLower ? k2+k1+actualPanelWidth : k2-actual_kc;

            pack_lhs(blockA, &tri(startTarget,startBlock), triStride, actualPanelWidth, lengthTarget);

Don Gagne's avatar
Don Gagne committed
160 161
            gebp_kernel(&other(startTarget,j2), otherStride, blockA, blockB+actual_kc*j2, lengthTarget, actualPanelWidth, actual_cols, Scalar(-1),
                        actualPanelWidth, actual_kc, 0, blockBOffset, blockW);
LM's avatar
LM committed
162 163 164
          }
        }
      }
Don Gagne's avatar
Don Gagne committed
165 166
      
      // R2 -= A21 * B => GEPP
LM's avatar
LM committed
167 168 169 170 171
      {
        Index start = IsLower ? k2+kc : 0;
        Index end   = IsLower ? size : k2-kc;
        for(Index i2=start; i2<end; i2+=mc)
        {
Don Gagne's avatar
Don Gagne committed
172
          const Index actual_mc = (std::min)(mc,end-i2);
LM's avatar
LM committed
173 174 175 176
          if (actual_mc>0)
          {
            pack_lhs(blockA, &tri(i2, IsLower ? k2 : k2-kc), triStride, actual_kc, actual_mc);

Don Gagne's avatar
Don Gagne committed
177
            gebp_kernel(_other+i2, otherStride, blockA, blockB, actual_mc, actual_kc, cols, Scalar(-1), -1, -1, 0, 0, blockW);
LM's avatar
LM committed
178 179 180 181 182 183 184 185 186 187 188 189 190 191
          }
        }
      }
    }
  }

/* Optimized triangular solver with multiple left hand sides and the trinagular matrix on the right
 */
template <typename Scalar, typename Index, int Mode, bool Conjugate, int TriStorageOrder>
struct triangular_solve_matrix<Scalar,Index,OnTheRight,Mode,Conjugate,TriStorageOrder,ColMajor>
{
  static EIGEN_DONT_INLINE void run(
    Index size, Index otherSize,
    const Scalar* _tri, Index triStride,
Don Gagne's avatar
Don Gagne committed
192 193 194 195 196 197 198 199 200
    Scalar* _other, Index otherStride,
    level3_blocking<Scalar,Scalar>& blocking);
};
template <typename Scalar, typename Index, int Mode, bool Conjugate, int TriStorageOrder>
EIGEN_DONT_INLINE void triangular_solve_matrix<Scalar,Index,OnTheRight,Mode,Conjugate,TriStorageOrder,ColMajor>::run(
    Index size, Index otherSize,
    const Scalar* _tri, Index triStride,
    Scalar* _other, Index otherStride,
    level3_blocking<Scalar,Scalar>& blocking)
LM's avatar
LM committed
201 202 203 204 205 206 207 208 209 210 211 212
  {
    Index rows = otherSize;
    const_blas_data_mapper<Scalar, Index, TriStorageOrder> rhs(_tri,triStride);
    blas_data_mapper<Scalar, Index, ColMajor> lhs(_other,otherStride);

    typedef gebp_traits<Scalar,Scalar> Traits;
    enum {
      RhsStorageOrder   = TriStorageOrder,
      SmallPanelWidth   = EIGEN_PLAIN_ENUM_MAX(Traits::mr,Traits::nr),
      IsLower = (Mode&Lower) == Lower
    };

Don Gagne's avatar
Don Gagne committed
213 214
    Index kc = blocking.kc();                   // cache block size along the K direction
    Index mc = (std::min)(rows,blocking.mc());  // cache block size along the M direction
LM's avatar
LM committed
215

Don Gagne's avatar
Don Gagne committed
216 217
    std::size_t sizeA = kc*mc;
    std::size_t sizeB = kc*size;
LM's avatar
LM committed
218
    std::size_t sizeW = kc*Traits::WorkSpaceFactor;
Don Gagne's avatar
Don Gagne committed
219 220 221 222

    ei_declare_aligned_stack_constructed_variable(Scalar, blockA, sizeA, blocking.blockA());
    ei_declare_aligned_stack_constructed_variable(Scalar, blockB, sizeB, blocking.blockB());
    ei_declare_aligned_stack_constructed_variable(Scalar, blockW, sizeW, blocking.blockW());
LM's avatar
LM committed
223 224 225 226 227 228 229 230 231 232 233

    conj_if<Conjugate> conj;
    gebp_kernel<Scalar,Scalar, Index, Traits::mr, Traits::nr, false, Conjugate> gebp_kernel;
    gemm_pack_rhs<Scalar, Index, Traits::nr,RhsStorageOrder> pack_rhs;
    gemm_pack_rhs<Scalar, Index, Traits::nr,RhsStorageOrder,false,true> pack_rhs_panel;
    gemm_pack_lhs<Scalar, Index, Traits::mr, Traits::LhsProgress, ColMajor, false, true> pack_lhs_panel;

    for(Index k2=IsLower ? size : 0;
        IsLower ? k2>0 : k2<size;
        IsLower ? k2-=kc : k2+=kc)
    {
Don Gagne's avatar
Don Gagne committed
234
      const Index actual_kc = (std::min)(IsLower ? k2 : size-k2, kc);
LM's avatar
LM committed
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
      Index actual_k2 = IsLower ? k2-actual_kc : k2 ;

      Index startPanel = IsLower ? 0 : k2+actual_kc;
      Index rs = IsLower ? actual_k2 : size - actual_k2 - actual_kc;
      Scalar* geb = blockB+actual_kc*actual_kc;

      if (rs>0) pack_rhs(geb, &rhs(actual_k2,startPanel), triStride, actual_kc, rs);

      // triangular packing (we only pack the panels off the diagonal,
      // neglecting the blocks overlapping the diagonal
      {
        for (Index j2=0; j2<actual_kc; j2+=SmallPanelWidth)
        {
          Index actualPanelWidth = std::min<Index>(actual_kc-j2, SmallPanelWidth);
          Index actual_j2 = actual_k2 + j2;
          Index panelOffset = IsLower ? j2+actualPanelWidth : 0;
          Index panelLength = IsLower ? actual_kc-j2-actualPanelWidth : j2;

          if (panelLength>0)
          pack_rhs_panel(blockB+j2*actual_kc,
                         &rhs(actual_k2+panelOffset, actual_j2), triStride,
                         panelLength, actualPanelWidth,
                         actual_kc, panelOffset);
        }
      }

      for(Index i2=0; i2<rows; i2+=mc)
      {
Don Gagne's avatar
Don Gagne committed
263
        const Index actual_mc = (std::min)(mc,rows-i2);
LM's avatar
LM committed
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

        // triangular solver kernel
        {
          // for each small block of the diagonal (=> vertical panels of rhs)
          for (Index j2 = IsLower
                      ? (actual_kc - ((actual_kc%SmallPanelWidth) ? Index(actual_kc%SmallPanelWidth)
                                                                  : Index(SmallPanelWidth)))
                      : 0;
               IsLower ? j2>=0 : j2<actual_kc;
               IsLower ? j2-=SmallPanelWidth : j2+=SmallPanelWidth)
          {
            Index actualPanelWidth = std::min<Index>(actual_kc-j2, SmallPanelWidth);
            Index absolute_j2 = actual_k2 + j2;
            Index panelOffset = IsLower ? j2+actualPanelWidth : 0;
            Index panelLength = IsLower ? actual_kc - j2 - actualPanelWidth : j2;

            // GEBP
            if(panelLength>0)
            {
              gebp_kernel(&lhs(i2,absolute_j2), otherStride,
                          blockA, blockB+j2*actual_kc,
                          actual_mc, panelLength, actualPanelWidth,
                          Scalar(-1),
                          actual_kc, actual_kc, // strides
                          panelOffset, panelOffset, // offsets
Don Gagne's avatar
Don Gagne committed
289
                          blockW);  // workspace
LM's avatar
LM committed
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
            }

            // unblocked triangular solve
            for (Index k=0; k<actualPanelWidth; ++k)
            {
              Index j = IsLower ? absolute_j2+actualPanelWidth-k-1 : absolute_j2+k;

              Scalar* r = &lhs(i2,j);
              for (Index k3=0; k3<k; ++k3)
              {
                Scalar b = conj(rhs(IsLower ? j+1+k3 : absolute_j2+k3,j));
                Scalar* a = &lhs(i2,IsLower ? j+1+k3 : absolute_j2+k3);
                for (Index i=0; i<actual_mc; ++i)
                  r[i] -= a[i] * b;
              }
              Scalar b = (Mode & UnitDiag) ? Scalar(1) : Scalar(1)/conj(rhs(j,j));
              for (Index i=0; i<actual_mc; ++i)
                r[i] *= b;
            }

            // pack the just computed part of lhs to A
            pack_lhs_panel(blockA, _other+absolute_j2*otherStride+i2, otherStride,
                           actualPanelWidth, actual_mc,
                           actual_kc, j2);
          }
        }

        if (rs>0)
          gebp_kernel(_other+i2+startPanel*otherStride, otherStride, blockA, geb,
                      actual_mc, actual_kc, rs, Scalar(-1),
Don Gagne's avatar
Don Gagne committed
320
                      -1, -1, 0, 0, blockW);
LM's avatar
LM committed
321 322 323 324 325 326
      }
    }
  }

} // end namespace internal

Don Gagne's avatar
Don Gagne committed
327 328
} // end namespace Eigen

LM's avatar
LM committed
329
#endif // EIGEN_TRIANGULAR_SOLVER_MATRIX_H