This file is indexed.

/usr/include/dune/istl/umfpack.hh is in libdune-istl-dev 2.4.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
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
// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
// vi: set et ts=4 sw=2 sts=2:
#ifndef DUNE_ISTL_UMFPACK_HH
#define DUNE_ISTL_UMFPACK_HH

// some trickery to make the preprocessor define setup
// of the backported UMFPACK test work with old code
#if ENABLE_UMFPACK && not defined(ENABLE_SUITESPARSE)
#define ENABLE_SUITESPARSE ENABLE_UMFPACK
#endif

#if HAVE_UMFPACK || defined DOXYGEN

#include<complex>
#include<type_traits>

#include<umfpack.h>

#include<dune/common/exceptions.hh>
#include<dune/common/fmatrix.hh>
#include<dune/common/fvector.hh>
#include<dune/common/unused.hh>
#include<dune/istl/bcrsmatrix.hh>
#include<dune/istl/solvers.hh>
#include<dune/istl/solvertype.hh>

#include"colcompmatrix.hh"


namespace Dune {
  /**
   * @addtogroup ISTL
   *
   * @{
   */
  /**
   * @file
   * @author Dominic Kempf
   * @brief Classes for using UMFPack 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 %UMFPack package to directly solve linear systems -- empty default class
   * @tparam Matrix the matrix type defining the system
   * Details on UMFPack can be found on
   * http://www.cise.ufl.edu/research/sparse/umfpack/
   */
  template<class Matrix>
  class UMFPack
  {};

  // wrapper class for C-Function Calls in the backend. Choose the right function namespace
  // depending on the template parameter used.
  template<typename T>
  struct UMFPackMethodChooser
  {};

  template<>
  struct UMFPackMethodChooser<double>
  {
    template<typename... A>
    static void defaults(A... args)
    {
      umfpack_di_defaults(args...);
    }
    template<typename... A>
    static void free_numeric(A... args)
    {
      umfpack_di_free_numeric(args...);
    }
    template<typename... A>
    static void free_symbolic(A... args)
    {
      umfpack_di_free_symbolic(args...);
    }
    template<typename... A>
    static int load_numeric(A... args)
    {
      return umfpack_di_load_numeric(args...);
    }
    template<typename... A>
    static void numeric(A... args)
    {
      umfpack_di_numeric(args...);
    }
    template<typename... A>
    static void report_info(A... args)
    {
      umfpack_di_report_info(args...);
    }
    template<typename... A>
    static void report_status(A... args)
    {
      umfpack_di_report_status(args...);
    }
    template<typename... A>
    static int save_numeric(A... args)
    {
      return umfpack_di_save_numeric(args...);
    }
    template<typename... A>
    static void solve(A... args)
    {
      umfpack_di_solve(args...);
    }
    template<typename... A>
    static void symbolic(A... args)
    {
      umfpack_di_symbolic(args...);
    }
  };

  template<>
  struct UMFPackMethodChooser<std::complex<double> >
  {
    template<typename... A>
    static void defaults(A... args)
    {
      umfpack_zi_defaults(args...);
    }
    template<typename... A>
    static void free_numeric(A... args)
    {
      umfpack_zi_free_numeric(args...);
    }
    template<typename... A>
    static void free_symbolic(A... args)
    {
      umfpack_zi_free_symbolic(args...);
    }
    template<typename... A>
    static int load_numeric(A... args)
    {
      return umfpack_zi_load_numeric(args...);
    }
    template<typename... A>
    static void numeric(const int* cs, const int* ri, const double* val, A... args)
    {
      umfpack_zi_numeric(cs,ri,val,NULL,args...);
    }
    template<typename... A>
    static void report_info(A... args)
    {
      umfpack_zi_report_info(args...);
    }
    template<typename... A>
    static void report_status(A... args)
    {
      umfpack_zi_report_status(args...);
    }
    template<typename... A>
    static int save_numeric(A... args)
    {
      return umfpack_zi_save_numeric(args...);
    }
    template<typename... A>
    static void solve(int m, const int* cs, const int* ri, std::complex<double>* val, double* x, const double* b,A... args)
    {
      const double* cval = reinterpret_cast<const double*>(val);
      umfpack_zi_solve(m,cs,ri,cval,NULL,x,NULL,b,NULL,args...);
    }
    template<typename... A>
    static void symbolic(int m, int n, const int* cs, const int* ri, const double* val, A... args)
    {
      umfpack_zi_symbolic(m,n,cs,ri,val,NULL,args...);
    }
  };

  /** @brief The %UMFPack direct sparse solver for matrices of type BCRSMatrix
   *
   * Specialization for the Dune::BCRSMatrix. %UMFPack 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 UMFPack
   */
  template<typename T, typename A, int n, int m>
  class UMFPack<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> UMFPackMatrix;
    /** @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 mat_ the matrix to solve for
     *  @param verbose [0..2] set the verbosity level, defaults to 0
     */
    UMFPack(const Matrix& matrix, int verbose=0) : matrixIsLoaded_(false)
    {
      //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 UMFPack (only double and std::complex<double> supported)");
      Caller::defaults(UMF_Control);
      setVerbosity(verbose);
      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 mat_ the matrix to solve for
     * @param verbose [0..2] set the verbosity level, defaults to 0
     */
    UMFPack(const Matrix& matrix, int verbose, bool) : matrixIsLoaded_(false)
    {
      //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 UMFPack (only double and std::complex<double> supported)");
      Caller::defaults(UMF_Control);
      setVerbosity(verbose);
      setMatrix(matrix);
    }

    /** @brief default constructor
     */
    UMFPack() : matrixIsLoaded_(false), verbosity_(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 UMFPack (only double and std::complex<double> supported)");
      Caller::defaults(UMF_Control);
    }

    /** @brief Try loading a decomposition from file and do a decomposition if unsuccessful
     * @param mat_ the matrix to decompose when no decoposition file found
     * @param file the decomposition file
     * @param verbose the verbosity level
     *
     * Use saveDecomposition(char* file) for manually storing a decomposition. This constructor
     * will decompose mat_ and store the result to file if no file wasn't found in the first place.
     * Thus, if you always use this you will only compute the decomposition once (and when you manually
     * deleted the decomposition file).
     */
    UMFPack(const Matrix& mat_, const char* file, int 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 UMFPack (only double and std::complex<double> supported)");
      Caller::defaults(UMF_Control);
      setVerbosity(verbose);
      int errcode = Caller::load_numeric(&UMF_Numeric, const_cast<char*>(file));
      if ((errcode == UMFPACK_ERROR_out_of_memory) || (errcode == UMFPACK_ERROR_file_IO))
      {
        matrixIsLoaded_ = false;
        setMatrix(mat_);
        saveDecomposition(file);
      }
      else
      {
        matrixIsLoaded_ = true;
        std::cout << "UMFPack decomposition successfully loaded from " << file << std::endl;
      }
    }

    /** @brief try loading a decomposition from file
     * @param file the decomposition file
     * @param verbose the verbosity level
     * @throws Dune::Exception When not being able to load the file. Does not need knowledge of the
     * actual matrix!
     */
    UMFPack(const char* file, int 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 UMFPack (only double and std::complex<double> supported)");
      Caller::defaults(UMF_Control);
      int errcode = Caller::load_numeric(&UMF_Numeric, const_cast<char*>(file));
      if (errcode == UMFPACK_ERROR_out_of_memory)
        DUNE_THROW(Dune::Exception, "ran out of memory while loading UMFPack decomposition");
      if (errcode == UMFPACK_ERROR_file_IO)
        DUNE_THROW(Dune::Exception, "IO error while loading UMFPack decomposition");
      matrixIsLoaded_ = true;
      std::cout << "UMFPack decomposition successfully loaded from " << file << std::endl;
      setVerbosity(verbose);
    }

    virtual ~UMFPack()
    {
      if ((umfpackMatrix_.N() + umfpackMatrix_.M() > 0) || (matrixIsLoaded_))
        free();
    }

    /**
     *  \copydoc InverseOperator::apply(X&, Y&, InverseOperatorResult&)
     */
    virtual void apply(domain_type& x, range_type& b, InverseOperatorResult& res)
    {
      if (umfpackMatrix_.N() != b.dim())
        DUNE_THROW(Dune::ISTLError, "Size of right-hand-side vector b does not match the number of matrix rows!");
      if (umfpackMatrix_.M() != x.dim())
        DUNE_THROW(Dune::ISTLError, "Size of solution vector x does not match the number of matrix columns!");

      double UMF_Apply_Info[UMFPACK_INFO];
      Caller::solve(UMFPACK_A,
                    umfpackMatrix_.getColStart(),
                    umfpackMatrix_.getRowIndex(),
                    umfpackMatrix_.getValues(),
                    reinterpret_cast<double*>(&x[0]),
                    reinterpret_cast<double*>(&b[0]),
                    UMF_Numeric,
                    UMF_Control,
                    UMF_Apply_Info);

      //this is a direct solver
      res.iterations = 1;
      res.converged = true;
      res.elapsed = UMF_Apply_Info[UMFPACK_SOLVE_WALLTIME];

      printOnApply(UMF_Apply_Info);
    }

    /**
     *  \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);
    }

    /**
     * @brief additional apply method with c-arrays in analogy to superlu
     * @param x solution array
     * @param b rhs array
     */
    void apply(T* x, T* b)
    {
      double UMF_Apply_Info[UMFPACK_INFO];
      Caller::solve(UMFPACK_A,
                    umfpackMatrix_.getColStart(),
                    umfpackMatrix_.getRowIndex(),
                    umfpackMatrix_.getValues(),
                    x,
                    b,
                    UMF_Numeric,
                    UMF_Control,
                    UMF_Apply_Info);
      printOnApply(UMF_Apply_Info);
    }

    /** @brief Set UMFPack-specific options
     *
     * This method allows to set various options that control the UMFPack solver.
     * More specifically, it allows to set values in the UMF_Control array.
     * Please see the UMFPack documentation for a list of possible options and values.
     *
     * \param option Entry in the UMF_Control array, e.g., UMFPACK_IRSTEP
     * \param value Corresponding value
     *
     * \throws RangeError If nonexisting option was requested
     */
    void setOption(unsigned int option, double value)
    {
      if (option >= UMFPACK_CONTROL)
        DUNE_THROW(RangeError, "Requested non-existing UMFPack option");

      UMF_Control[option] = value;
    }

    /** @brief saves a decomposition to a file
     * @param file the filename to save to
     */
    void saveDecomposition(const char* file)
    {
      int errcode = Caller::save_numeric(UMF_Numeric, const_cast<char*>(file));
      if (errcode != UMFPACK_OK)
        DUNE_THROW(Dune::Exception,"IO ERROR while trying to save UMFPack decomposition");
    }

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

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

    /** @brief sets the verbosity level for the UMFPack solver
     * @param v verbosity level
     * The following levels are implemented:
     * 0 - only error messages
     * 1 - a bit of statistics on decomposition and solution
     * 2 - lots of statistics on decomposition and solution
     */
    void setVerbosity(int v)
    {
      verbosity_ = v;
      // set the verbosity level in UMFPack
      if (verbosity_ == 0)
        UMF_Control[UMFPACK_PRL] = 1;
      if (verbosity_ == 1)
        UMF_Control[UMFPACK_PRL] = 2;
      if (verbosity_ == 2)
        UMF_Control[UMFPACK_PRL] = 4;
    }

    /**
     * @brief free allocated space.
     * @warning later calling apply will result in an error.
     */
    void free()
    {
      if (!matrixIsLoaded_)
      {
        Caller::free_symbolic(&UMF_Symbolic);
        umfpackMatrix_.free();
      }
      Caller::free_numeric(&UMF_Numeric);
      matrixIsLoaded_ = false;
    }

    const char* name() { return "UMFPACK"; }

    private:
    typedef typename Dune::UMFPackMethodChooser<T> Caller;

    template<class M,class X, class TM, class TD, class T1>
    friend class SeqOverlappingSchwarz;
    friend struct SeqOverlappingSchwarzAssemblerHelper<UMFPack<Matrix>,true>;

    /** @brief computes the LU Decomposition */
    void decompose()
    {
      double UMF_Decomposition_Info[UMFPACK_INFO];
      Caller::symbolic(static_cast<int>(umfpackMatrix_.N()),
                       static_cast<int>(umfpackMatrix_.N()),
                       umfpackMatrix_.getColStart(),
                       umfpackMatrix_.getRowIndex(),
                       reinterpret_cast<double*>(umfpackMatrix_.getValues()),
                       &UMF_Symbolic,
                       UMF_Control,
                       UMF_Decomposition_Info);
      Caller::numeric(umfpackMatrix_.getColStart(),
                      umfpackMatrix_.getRowIndex(),
                      reinterpret_cast<double*>(umfpackMatrix_.getValues()),
                      UMF_Symbolic,
                      &UMF_Numeric,
                      UMF_Control,
                      UMF_Decomposition_Info);
      Caller::report_status(UMF_Control,UMF_Decomposition_Info[UMFPACK_STATUS]);
      if (verbosity_ == 1)
      {
        std::cout << "[UMFPack Decomposition]" << std::endl;
        std::cout << "Wallclock Time taken: " << UMF_Decomposition_Info[UMFPACK_NUMERIC_WALLTIME] << " (CPU Time: " << UMF_Decomposition_Info[UMFPACK_NUMERIC_TIME] << ")" << std::endl;
        std::cout << "Flops taken: " << UMF_Decomposition_Info[UMFPACK_FLOPS] << std::endl;
        std::cout << "Peak Memory Usage: " << UMF_Decomposition_Info[UMFPACK_PEAK_MEMORY]*UMF_Decomposition_Info[UMFPACK_SIZE_OF_UNIT] << " bytes" << std::endl;
        std::cout << "Condition number estimate: " << 1./UMF_Decomposition_Info[UMFPACK_RCOND] << std::endl;
        std::cout << "Numbers of non-zeroes in decomposition: L: " << UMF_Decomposition_Info[UMFPACK_LNZ] << " U: " << UMF_Decomposition_Info[UMFPACK_UNZ] << std::endl;
      }
      if (verbosity_ == 2)
      {
        Caller::report_info(UMF_Control,UMF_Decomposition_Info);
      }
    }

    void printOnApply(double* UMF_Info)
    {
      Caller::report_status(UMF_Control,UMF_Info[UMFPACK_STATUS]);
      if (verbosity_ > 0)
      {
        std::cout << "[UMFPack Solve]" << std::endl;
        std::cout << "Wallclock Time: " << UMF_Info[UMFPACK_SOLVE_WALLTIME] << " (CPU Time: " << UMF_Info[UMFPACK_SOLVE_TIME] << ")" << std::endl;
        std::cout << "Flops Taken: " << UMF_Info[UMFPACK_SOLVE_FLOPS] << std::endl;
        std::cout << "Iterative Refinement steps taken: " << UMF_Info[UMFPACK_IR_TAKEN] << std::endl;
        std::cout << "Error Estimate: " << UMF_Info[UMFPACK_OMEGA1] << " resp. " << UMF_Info[UMFPACK_OMEGA2] << std::endl;
      }
    }

    UMFPackMatrix& getInternalMatrix() { return umfpackMatrix_; }

    UMFPackMatrix umfpackMatrix_;
    bool matrixIsLoaded_;
    int verbosity_;
    void *UMF_Symbolic;
    void *UMF_Numeric;
    double UMF_Control[UMFPACK_CONTROL];
  };

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

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

#endif //HAVE_UMFPACK

#endif //DUNE_ISTL_UMFPACK_HH