This file is indexed.

/usr/include/dune/istl/spqr.hh is in libdune-istl-dev 2.5.1-1.

This file is owned by root:root, with mode 0o644.

The actual contents of the file can be viewed below.

  1
  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
 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
 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
// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
// vi: set et ts=4 sw=2 sts=2:
#ifndef DUNE_ISTL_SPQR_HH
#define DUNE_ISTL_SPQR_HH

#if HAVE_SUITESPARSE_SPQR || defined DOXYGEN

#include <complex>
#include <type_traits>

#include <SuiteSparseQR.hpp>

#include <dune/common/exceptions.hh>
#include <dune/common/unused.hh>

#include <dune/istl/colcompmatrix.hh>
#include <dune/istl/solvers.hh>
#include <dune/istl/solvertype.hh>

namespace Dune {
  /**
   * @addtogroup ISTL
   *
   * @{
   */
  /**
   * @file
   * @author Marco Agnese, Andrea Sacconi
   * @brief Class for using SPQR with ISTL matrices.
   */

  // forward declarations
  template<class M, class T, class TM, class TD, class TA>
  class SeqOverlappingSchwarz;

  template<class T, bool tag>
  struct SeqOverlappingSchwarzAssemblerHelper;

  /** @brief Use the %SPQR package to directly solve linear systems -- empty default class
   * @tparam Matrix the matrix type defining the system
   * Details on SPQR can be found on
   * http://www.cise.ufl.edu/research/sparse/spqr/
   */
  template<class Matrix>
  class SPQR
  {};

  /** @brief The %SPQR direct sparse solver for matrices of type BCRSMatrix
   *
   * Specialization for the Dune::BCRSMatrix. %SPQR will always go double
   * precision and supports complex numbers
   * too (use std::complex<double> for that).
   *
   * \tparam T Number type.  Only double and std::complex<double> is supported
   * \tparam A STL-compatible allocator type
   * \tparam n Number of rows in a matrix block
   * \tparam m Number of columns in a matrix block
   *
   * \note This will only work if dune-istl has been configured to use SPQR
   */
  template<typename T, typename A, int n, int m>
  class SPQR<BCRSMatrix<FieldMatrix<T,n,m>,A > >
    : public InverseOperator<BlockVector<FieldVector<T,m>, typename A::template rebind<FieldVector<T,m> >::other>,
                             BlockVector<FieldVector<T,n>, typename A::template rebind<FieldVector<T,n> >::other> >
  {
    public:
    /** @brief The matrix type. */
    typedef Dune::BCRSMatrix<FieldMatrix<T,n,m>,A> Matrix;
    typedef Dune::BCRSMatrix<FieldMatrix<T,n,m>,A> matrix_type;
    /** @brief The corresponding SuperLU Matrix type.*/
    typedef Dune::ColCompMatrix<Matrix> SPQRMatrix;
    /** @brief Type of an associated initializer class. */
    typedef ColCompMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A> > MatrixInitializer;
    /** @brief The type of the domain of the solver. */
    typedef Dune::BlockVector<FieldVector<T,m>, typename A::template rebind<FieldVector<T,m> >::other> domain_type;
    /** @brief The type of the range of the solver. */
    typedef Dune::BlockVector<FieldVector<T,n>, typename A::template rebind<FieldVector<T,n> >::other> range_type;

    /** @brief Construct a solver object from a BCRSMatrix
     *
     * This computes the matrix decomposition, and may take a long time
     * (and use a lot of memory).
     *
     *  @param matrix the matrix to solve for
     *  @param verbose, 0 or 1, set the verbosity level, defaults to 0
     */
    SPQR(const Matrix& matrix, int verbose=0) : matrixIsLoaded_(false), verbose_(verbose)
    {
      //check whether T is a supported type
      static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
                    "Unsupported Type in SPQR (only double and std::complex<double> supported)");
      cc_ = new cholmod_common();
      cholmod_l_start(cc_);
      setMatrix(matrix);
    }

     /** @brief Constructor for compatibility with SuperLU standard constructor
     *
     * This computes the matrix decomposition, and may take a long time
     * (and use a lot of memory).
     *
     * @param matrix the matrix to solve for
     * @param verbose, 0 or 1, set the verbosity level, defaults to 0
     */
    SPQR(const Matrix& matrix, int verbose, bool) : matrixIsLoaded_(false), verbose_(verbose)
    {
      //check whether T is a supported type
      static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
                    "Unsupported Type in SPQR (only double and std::complex<double> supported)");
      cc_ = new cholmod_common();
      cholmod_l_start(cc_);
      setMatrix(matrix);
    }

    /** @brief Default constructor. */
    SPQR() : matrixIsLoaded_(false), verbose_(0)
    {
      //check whether T is a supported type
      static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
                    "Unsupported Type in SPQR (only double and std::complex<double> supported)");
      cc_ = new cholmod_common();
      cholmod_l_start(cc_);
    }

    /** @brief Destructor. */
    virtual ~SPQR()
    {
      if ((spqrMatrix_.N() + spqrMatrix_.M() > 0) || matrixIsLoaded_)
        free();
      cholmod_l_finish(cc_);
    }

    /** \copydoc InverseOperator::apply(X&, Y&, InverseOperatorResult&) */
    virtual void apply(domain_type& x, range_type& b, InverseOperatorResult& res)
    {
      const std::size_t dimMat(spqrMatrix_.N());
      // fill B
      for(std::size_t k = 0; k != dimMat; ++k)
        (static_cast<T*>(B_->x))[k] = b[k];
      cholmod_dense* BTemp = B_;
      B_ = SuiteSparseQR_qmult<T>(0, spqrfactorization_, B_, cc_);
      cholmod_dense* X = SuiteSparseQR_solve<T>(1, spqrfactorization_, B_, cc_);
      cholmod_l_free_dense(&BTemp, cc_);
      // fill x
      for(std::size_t k = 0; k != dimMat; ++k)
        x [k] = (static_cast<T*>(X->x))[k];
      cholmod_l_free_dense(&X, cc_);
      // this is a direct solver
      res.iterations = 1;
      res.converged = true;
      if(verbose_ > 0)
      {
        std::cout<<std::endl<<"Solving with SuiteSparseQR"<<std::endl;
        std::cout<<"Flops Taken: "<<cc_->SPQR_flopcount<<std::endl;
        std::cout<<"Analysis Time: "<<cc_->SPQR_analyze_time<<" s"<<std::endl;
        std::cout<<"Factorize Time: "<<cc_->SPQR_factorize_time<<" s"<<std::endl;
        std::cout<<"Backsolve Time: "<<cc_->SPQR_solve_time<<" s"<<std::endl;
        std::cout<<"Peak Memory Usage: "<<cc_->memory_usage<<" bytes"<<std::endl;
        std::cout<<"Rank Estimate: "<<cc_->SPQR_istat[4]<<std::endl<<std::endl;
      }
    }

    /** \copydoc InverseOperator::apply(X&,Y&,double,InverseOperatorResult&) */
    virtual void apply (domain_type& x, range_type& b, double reduction, InverseOperatorResult& res)
    {
      DUNE_UNUSED_PARAMETER(reduction);
      apply(x, b, res);
    }

    void setOption(unsigned int option, double value)
    {
      DUNE_UNUSED_PARAMETER(option);
      DUNE_UNUSED_PARAMETER(value);
    }

    /** @brief Initialize data from given matrix. */
    void setMatrix(const Matrix& matrix)
    {
      if ((spqrMatrix_.N() + spqrMatrix_.M() > 0) || matrixIsLoaded_)
        free();
      spqrMatrix_ = matrix;
      decompose();
    }

    template<class S>
    void setSubMatrix(const Matrix& matrix, const S& rowIndexSet)
    {
      if ((spqrMatrix_.N() + spqrMatrix_.M() > 0) || matrixIsLoaded_)
        free();
      spqrMatrix_.setMatrix(matrix,rowIndexSet);
      decompose();
    }

    /**
     * @brief Sets the verbosity level for the solver.
     * @param v verbosity level: 0 only error messages, 1 a bit of statistics.
     */
    inline void setVerbosity(int v)
    {
      verbose_=v;
    }

    /**
     * @brief Return the matrix factorization.
     * @warning It is up to the user to keep consistency.
     */
    inline SuiteSparseQR_factorization<T>* getFactorization()
    {
      return spqrfactorization_;
    }

    /**
     * @brief Return the column coppressed matrix.
     * @warning It is up to the user to keep consistency.
     */
    inline SPQRMatrix& getInternalMatrix()
    {
      return spqrMatrix_;
    }

    /**
     * @brief Free allocated space.
     * @warning Later calling apply will result in an error.
     */
    void free()
    {
      cholmod_l_free_sparse(&A_, cc_);
      cholmod_l_free_dense(&B_, cc_);
      SuiteSparseQR_free<T>(&spqrfactorization_, cc_);
      spqrMatrix_.free();
      matrixIsLoaded_ = false;
    }

    /** @brief Get method name. */
    inline const char* name()
    {
      return "SPQR";
    }

    private:
    template<class M,class X, class TM, class TD, class T1>
    friend class SeqOverlappingSchwarz;

    friend struct SeqOverlappingSchwarzAssemblerHelper<SPQR<Matrix>,true>;

    /** @brief Computes the QR decomposition. */
    void decompose()
    {
      const std::size_t dimMat(spqrMatrix_.N());
      const std::size_t nnz(spqrMatrix_.getColStart()[dimMat]);
      // initialise the matrix A (sorted, packed, unsymmetric, real entries)
      A_ = cholmod_l_allocate_sparse(dimMat, dimMat, nnz, 1, 1, 0, 1, cc_);
      // copy all the entries of Ap, Ai, Ax
      for(std::size_t k = 0; k != (dimMat+1); ++k)
        (static_cast<long int *>(A_->p))[k] = spqrMatrix_.getColStart()[k];
      for(std::size_t k = 0; k != nnz; ++k)
      {
        (static_cast<long int*>(A_->i))[k] = spqrMatrix_.getRowIndex()[k];
        (static_cast<T*>(A_->x))[k] = spqrMatrix_.getValues()[k];
      }
      // initialise the vector B
      B_ = cholmod_l_allocate_dense(dimMat, 1, dimMat, A_->xtype, cc_);
      // compute factorization of A
      spqrfactorization_=SuiteSparseQR_factorize<T>(SPQR_ORDERING_DEFAULT,SPQR_DEFAULT_TOL,A_,cc_);
    }

    SPQRMatrix spqrMatrix_;
    bool matrixIsLoaded_;
    int verbose_;
    cholmod_common* cc_;
    cholmod_sparse* A_;
    cholmod_dense* B_;
    SuiteSparseQR_factorization<T>* spqrfactorization_;
  };

  template<typename T, typename A, int n, int m>
  struct IsDirectSolver<SPQR<BCRSMatrix<FieldMatrix<T,n,m>,A> > >
  {
    enum {value = true};
  };

  template<typename T, typename A, int n, int m>
  struct StoresColumnCompressed<SPQR<BCRSMatrix<FieldMatrix<T,n,m>,A> > >
  {
    enum {value = true};
  };

}

#endif //HAVE_SUITESPARSE_SPQR
#endif //DUNE_ISTL_SPQR_HH