/usr/include/dune/istl/supermatrix.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 | // -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
// vi: set et ts=4 sw=2 sts=2:
#ifndef DUNE_ISTL_SUPERMATRIX_HH
#define DUNE_ISTL_SUPERMATRIX_HH
#if HAVE_SUPERLU
#ifdef SUPERLU_POST_2005_VERSION
#ifndef SUPERLU_NTYPE
#define SUPERLU_NTYPE 1
#endif
#if SUPERLU_NTYPE==0
#include "slu_sdefs.h"
#endif
#if SUPERLU_NTYPE==1
#include "slu_ddefs.h"
#endif
#if SUPERLU_NTYPE==2
#include "slu_cdefs.h"
#endif
#if SUPERLU_NTYPE>=3
#include "slu_zdefs.h"
#endif
#else
#if SUPERLU_NTYPE==0
#include "ssp_defs.h"
#endif
#if SUPERLU_NTYPE==1
#include "dsp_defs.h"
#endif
#if SUPERLU_NTYPE==2
#include "csp_defs.h"
#endif
#if SUPERLU_NTYPE>=3
#include "zsp_defs.h"
#endif
#endif
#include "bcrsmatrix.hh"
#include "bvector.hh"
#include <dune/common/fmatrix.hh>
#include <dune/common/fvector.hh>
#include <dune/common/typetraits.hh>
#include <limits>
#include"colcompmatrix.hh"
namespace Dune
{
template<class T>
struct SuperMatrixCreateSparseChooser
{};
template<class T>
struct SuperMatrixPrinter
{};
#if SUPERLU_NTYPE==0
template<>
struct SuperMatrixCreateSparseChooser<float>
{
static void create(SuperMatrix *mat, int n, int m, int offset,
float *values, int *rowindex, int* colindex,
Stype_t stype, Dtype_t dtype, Mtype_t mtype)
{
sCreate_CompCol_Matrix(mat, n, m, offset, values, rowindex, colindex,
stype, dtype, mtype);
}
};
template<>
struct SuperMatrixPrinter<float>
{
static void print(char* name, SuperMatrix* mat)
{
sPrint_CompCol_Matrix(name, mat);
}
};
#endif
#if SUPERLU_NTYPE==1
template<>
struct SuperMatrixCreateSparseChooser<double>
{
static void create(SuperMatrix *mat, int n, int m, int offset,
double *values, int *rowindex, int* colindex,
Stype_t stype, Dtype_t dtype, Mtype_t mtype)
{
dCreate_CompCol_Matrix(mat, n, m, offset, values, rowindex, colindex,
stype, dtype, mtype);
}
};
template<>
struct SuperMatrixPrinter<double>
{
static void print(char* name, SuperMatrix* mat)
{
dPrint_CompCol_Matrix(name, mat);
}
};
#endif
#if SUPERLU_NTYPE==2
template<>
struct SuperMatrixCreateSparseChooser<std::complex<float> >
{
static void create(SuperMatrix *mat, int n, int m, int offset,
std::complex<float> *values, int *rowindex, int* colindex,
Stype_t stype, Dtype_t dtype, Mtype_t mtype)
{
cCreate_CompCol_Matrix(mat, n, m, offset, reinterpret_cast< ::complex*>(values),
rowindex, colindex, stype, dtype, mtype);
}
};
template<>
struct SuperMatrixPrinter<std::complex<float> >
{
static void print(char* name, SuperMatrix* mat)
{
cPrint_CompCol_Matrix(name, mat);
}
};
#endif
#if SUPERLU_NTYPE>=3
template<>
struct SuperMatrixCreateSparseChooser<std::complex<double> >
{
static void create(SuperMatrix *mat, int n, int m, int offset,
std::complex<double> *values, int *rowindex, int* colindex,
Stype_t stype, Dtype_t dtype, Mtype_t mtype)
{
zCreate_CompCol_Matrix(mat, n, m, offset, reinterpret_cast<doublecomplex*>(values),
rowindex, colindex, stype, dtype, mtype);
}
};
template<>
struct SuperMatrixPrinter<std::complex<double> >
{
static void print(char* name, SuperMatrix* mat)
{
zPrint_CompCol_Matrix(name, mat);
}
};
#endif
template<class T>
struct BaseGetSuperLUType
{
static const Dtype_t type;
};
template<class T>
struct GetSuperLUType
{};
template<class T>
const Dtype_t BaseGetSuperLUType<T>::type =
Dune::is_same<T,float>::value ? SLU_S :
( Dune::is_same<T,std::complex<double> >::value ? SLU_Z :
( Dune::is_same<T,std::complex<float> >::value ? SLU_C : SLU_D ));
template<>
struct GetSuperLUType<double>
: public BaseGetSuperLUType<double>
{
typedef double float_type;
};
template<>
struct GetSuperLUType<float>
: public BaseGetSuperLUType<float>
{
typedef float float_type;
};
template<>
struct GetSuperLUType<std::complex<double> >
: public BaseGetSuperLUType<std::complex<double> >
{
typedef double float_type;
};
template<>
struct GetSuperLUType<std::complex<float> >
: public BaseGetSuperLUType<std::complex<float> >
{
typedef float float_type;
};
/**
* @brief Utility class for converting an ISTL Matrix
* into a SuperLU Matrix.
*/
template<class M>
struct SuperLUMatrix
{};
template<class M>
struct SuperMatrixInitializer
{};
template<class T>
class SuperLU;
/**
* @brief Converter for BCRSMatrix to SuperLU Matrix.
*/
template<class B, class TA, int n, int m>
class SuperLUMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >
: public ColCompMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >
{
template<class M, class X, class TM, class TD, class T1>
friend class SeqOverlappingSchwarz;
friend struct SuperMatrixInitializer<BCRSMatrix<FieldMatrix<B,n,m>,TA> >;
public:
/** @brief The type of the matrix to convert. */
typedef BCRSMatrix<FieldMatrix<B,n,m>,TA> Matrix;
friend struct SeqOverlappingSchwarzAssembler<SuperLU<Matrix> >;
typedef typename Matrix::size_type size_type;
/**
* @brief Constructor that initializes the data.
* @param mat The matrix to convert.
*/
explicit SuperLUMatrix(const Matrix& mat) : ColCompMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >(mat)
{}
SuperLUMatrix() : ColCompMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >()
{}
/** @brief Destructor */
virtual ~SuperLUMatrix()
{
if (this->N_+this->M_*this->Nnz_ != 0)
free();
}
/** @brief Cast to a SuperLU Matrix */
operator SuperMatrix&()
{
return A;
}
/** @brief Cast to a SuperLU Matrix */
operator const SuperMatrix&() const
{
return A;
}
SuperLUMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >& operator=(const BCRSMatrix<FieldMatrix<B,n,m>,TA>& mat)
{
this->ColCompMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >::operator=(mat);
SuperMatrixCreateSparseChooser<B>
::create(&A, this->N_, this->M_, this->colstart[this->N_],
this->values,this->rowindex, this->colstart, SLU_NC,
static_cast<Dtype_t>(GetSuperLUType<B>::type), SLU_GE);
return *this;
}
SuperLUMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >& operator=(const SuperLUMatrix <BCRSMatrix<FieldMatrix<B,n,m>,TA> >& mat)
{
this->ColCompMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >::operator=(mat);
SuperMatrixCreateSparseChooser<B>
::create(&A, this->N_, this->M_, this->colstart[this->N_],
this->values,this->rowindex, this->colstart, SLU_NC,
static_cast<Dtype_t>(GetSuperLUType<B>::type), SLU_GE);
return *this;
}
/**
* @brief Initialize data from a given set of matrix rows and columns
* @tparam The type of the row index set.
* @param mat the matrix with the values
* @param mrs The set of row (and column) indices to represent
*/
virtual void setMatrix(const Matrix& mat, const std::set<std::size_t>& mrs)
{
if(this->N_+this->M_+this->Nnz_!=0)
free();
this->N_=mrs.size()*n;
this->M_=mrs.size()*m;
SuperMatrixInitializer<Matrix> initializer(*this);
copyToColCompMatrix(initializer, MatrixRowSubset<Matrix,std::set<std::size_t> >(mat,mrs));
}
/** @brief Initialize data from given matrix. */
virtual void setMatrix(const Matrix& mat)
{
this->N_=n*mat.N();
this->M_=m*mat.M();
SuperMatrixInitializer<Matrix> initializer(*this);
copyToColCompMatrix(initializer, MatrixRowSet<Matrix>(mat));
}
/** @brief free allocated space. */
virtual void free()
{
ColCompMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >::free();
SUPERLU_FREE(A.Store);
}
private:
SuperMatrix A;
};
template<class T, class A, int n, int m>
class SuperMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A> >
: public ColCompMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A> >
{
template<class I, class S, class D>
friend class OverlappingSchwarzInitializer;
public:
typedef BCRSMatrix<FieldMatrix<T,n,m>,A> Matrix;
typedef Dune::SuperLUMatrix<Matrix> SuperLUMatrix;
SuperMatrixInitializer(SuperLUMatrix& lum) : ColCompMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A> >(lum)
,slumat(&lum)
{}
SuperMatrixInitializer() : ColCompMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A> >()
{}
virtual void createMatrix() const
{
ColCompMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A> >::createMatrix();
SuperMatrixCreateSparseChooser<T>
::create(&slumat->A, slumat->N_, slumat->M_, slumat->colstart[this->cols],
slumat->values,slumat->rowindex, slumat->colstart, SLU_NC,
static_cast<Dtype_t>(GetSuperLUType<T>::type), SLU_GE);
}
private:
SuperLUMatrix* slumat;
};
}
#endif
#endif
|