dune-fem  2.4.1-rc
umfpacksolver.hh
Go to the documentation of this file.
1 #ifndef DUNE_FEM_UMFPACKSOLVER_HH
2 #define DUNE_FEM_UMFPACKSOLVER_HH
3 
4 #include <limits>
5 #include <iostream>
6 #include <type_traits>
7 #include <vector>
8 #include <algorithm>
9 
12 
14 #include <dune/fem/io/parameter.hh>
16 
17 
18 #if HAVE_DUNE_ISTL
19 #include <dune/istl/umfpack.hh>
20 
21 #ifdef ENABLE_UMFPACK
22 
23 namespace Dune
24 {
25 namespace Fem
26 {
27 
44 template<class DF, class Op, bool symmetric=false>
45 class UMFPACKOp:public Operator<DF, DF>
46 {
47  public:
48  typedef DF DiscreteFunctionType;
49  typedef Op OperatorType;
50 
51  // \brief The column-compressed matrix type.
52  typedef ColCompMatrix<typename OperatorType::MatrixType::MatrixBaseType> CCSMatrixType;
53  typedef typename DiscreteFunctionType::DofType DofType;
54  typedef typename DiscreteFunctionType::DiscreteFunctionSpaceType DiscreteFunctionSpaceType;
55 
63  UMFPACKOp(const OperatorType& op, const double& redEps, const double& absLimit, const int& maxIter, const bool& verbose,
64  const ParameterReader &parameter = Parameter::container() ) :
65  op_(op), verbose_(verbose), ccsmat_(), isloaded_(false)
66  {
67  Caller::defaults(UMF_Control);
68  UMF_Control[UMFPACK_PRL] = 4;
69  }
70 
77  UMFPACKOp(const OperatorType& op, const double& redEps, const double& absLimit,
78  const int& maxIter, const ParameterReader& parameter = Parameter::container() ) :
79  op_(op), verbose_(parameter.getValue<bool>("fem.solver.verbose",false)), ccsmat_(), isloaded_(false)
80  {
81  Caller::defaults(UMF_Control);
82  UMF_Control[UMFPACK_PRL] = 4;
83  }
84 
85  UMFPACKOp(const OperatorType& op, const ParameterReader& parameter = Parameter::container() ) :
86  op_(op), verbose_(parameter.getValue<bool>("fem.solver.verbose",false)), ccsmat_(), isloaded_(false)
87  {
88  Caller::defaults(UMF_Control);
89  UMF_Control[UMFPACK_PRL] = 4;
90  }
91 
92  // \brief Destructor.
93  ~UMFPACKOp()
94  {}
95 
100  void operator()(const DiscreteFunctionType& arg, DiscreteFunctionType& dest) const
101  {
102  prepare();
103  apply(arg,dest);
104  finalize();
105  }
106 
107  // \brief Decompose matrix.
108  template<typename... A>
109  void prepare(A... ) const
110  {
111  if(!isloaded_)
112  {
113  ccsmat_ = op_.systemMatrix().matrix();
114  decompose();
115  isloaded_ = true;
116  }
117  }
118 
119  // \brief Free allocated memory.
120  inline void finalize() const
121  {
122  if(isloaded_)
123  {
124  ccsmat_.free();
125  Caller::free_symbolic(&UMF_Symbolic);
126  Caller::free_numeric(&UMF_Numeric);
127  isloaded_ = false;
128  }
129  }
130 
137  void apply(const DofType*& arg, DofType*& dest) const
138  {
139  double UMF_Apply_Info[UMFPACK_INFO];
140  Caller::solve(UMFPACK_A, ccsmat_.getColStart(), ccsmat_.getRowIndex(), ccsmat_.getValues(),
141  dest, const_cast<DofType*>(arg), UMF_Numeric, UMF_Control, UMF_Apply_Info);
142  if(verbose_)
143  {
144  Caller::report_status(UMF_Control, UMF_Apply_Info[UMFPACK_STATUS]);
145  std::cout <<"[UMFPack Solve]" << std::endl;
146  std::cout << "Wallclock Time: " << UMF_Apply_Info[UMFPACK_SOLVE_WALLTIME]
147  << " (CPU Time: " << UMF_Apply_Info[UMFPACK_SOLVE_TIME] << ")" << std::endl;
148  std::cout << "Flops Taken: " << UMF_Apply_Info[UMFPACK_SOLVE_FLOPS] << std::endl;
149  std::cout << "Iterative Refinement steps taken: " << UMF_Apply_Info[UMFPACK_IR_TAKEN] << std::endl;
150  std::cout << "Error Estimate: " << UMF_Apply_Info[UMFPACK_OMEGA1] << " resp. " << UMF_Apply_Info[UMFPACK_OMEGA2] << std::endl;
151  }
152  }
153 
160  void apply(const AdaptiveDiscreteFunction<DiscreteFunctionSpaceType>& arg,
161  AdaptiveDiscreteFunction<DiscreteFunctionSpaceType>& dest) const
162  {
163  const DofType* argPtr(arg.leakPointer());
164  DofType* destPtr(dest.leakPointer());
165  apply(argPtr,destPtr);
166  }
167 
174  void apply(const ISTLBlockVectorDiscreteFunction<DiscreteFunctionSpaceType>& arg,
175  ISTLBlockVectorDiscreteFunction<DiscreteFunctionSpaceType>& dest) const
176  {
177  // copy DOF arg into a consecutive vector
178  std::vector<DofType> vecArg(arg.size());
179  std::copy(arg.dbegin(),arg.dend(),vecArg.begin());
180  std::vector<DofType> vecDest(dest.size());
181  // apply operator
182  const DofType* argDataPtr(vecArg.data());
183  DofType* destDataPtr(vecDest.data());
184  apply(argDataPtr,destDataPtr);
185  // copy back solution into dest
186  std::copy(vecDest.begin(),vecDest.end(),dest.dbegin());
187  }
188 
189  inline void printTexInfo(std::ostream& out) const
190  {
191  out<<"Solver: UMFPACK direct solver";
192  out<<"\\\\ \n";
193  }
194 
195  inline double averageCommTime() const
196  {
197  return 0.0;
198  }
199 
200  inline int iterations() const
201  {
202  return 0;
203  }
204 
208  inline CCSMatrixType& getCCSMatrix()
209  {
210  return ccsmat_;
211  }
212 
213  // /brief Print some statistics about the UMFPACK decomposition.
214  inline void printDecompositionInfo() const
215  {
216  Caller::report_info(UMF_Control,UMF_Decomposition_Info);
217  }
218 
219  private:
220  const OperatorType& op_;
221  const bool verbose_;
222  mutable CCSMatrixType ccsmat_;
223  mutable bool isloaded_;
224  mutable void *UMF_Symbolic;
225  mutable void *UMF_Numeric;
226  mutable double UMF_Control[UMFPACK_CONTROL];
227  mutable double UMF_Decomposition_Info[UMFPACK_INFO];
228 
229  typedef typename Dune::UMFPackMethodChooser<DofType> Caller;
230 
231  // /brief Computes the UMFPACK decomposition.
232  void decompose() const
233  {
234  const std::size_t dimMat(ccsmat_.N());
235  Caller::symbolic(static_cast<int>(dimMat), static_cast<int>(dimMat), ccsmat_.getColStart(), ccsmat_.getRowIndex(),
236  reinterpret_cast<double*>(ccsmat_.getValues()), &UMF_Symbolic, UMF_Control, UMF_Decomposition_Info);
237  Caller::numeric(ccsmat_.getColStart(), ccsmat_.getRowIndex(), reinterpret_cast<double*>(ccsmat_.getValues()),
238  UMF_Symbolic, &UMF_Numeric, UMF_Control, UMF_Decomposition_Info);
239  if(verbose_)
240  {
241  Caller::report_status(UMF_Control,UMF_Decomposition_Info[UMFPACK_STATUS]);
242  std::cout << "[UMFPack Decomposition]" << std::endl;
243  std::cout << "Wallclock Time taken: " << UMF_Decomposition_Info[UMFPACK_NUMERIC_WALLTIME]
244  << " (CPU Time: " << UMF_Decomposition_Info[UMFPACK_NUMERIC_TIME] << ")" << std::endl;
245  std::cout << "Flops taken: " << UMF_Decomposition_Info[UMFPACK_FLOPS] << std::endl;
246  std::cout << "Peak Memory Usage: " << UMF_Decomposition_Info[UMFPACK_PEAK_MEMORY]*UMF_Decomposition_Info[UMFPACK_SIZE_OF_UNIT]
247  << " bytes" << std::endl;
248  std::cout << "Condition number estimate: " << 1./UMF_Decomposition_Info[UMFPACK_RCOND] << std::endl;
249  std::cout << "Numbers of non-zeroes in decomposition: L: " << UMF_Decomposition_Info[UMFPACK_LNZ]
250  << " U: " << UMF_Decomposition_Info[UMFPACK_UNZ] << std::endl;
251  }
252  }
253 };
254 
255 }
256 }
257 
258 #endif // #if HAVE_DUNE_ISTL
259 
260 #endif // #if HAVE_UMFPACK
261 
262 #endif // #ifndef DUNE_FEM_UMFPACKSOLVER_HH
BasicParameterReader< std::function< const std::string *(const std::string &, const std::string *) > > ParameterReader
Definition: reader.hh:264
Definition: coordinate.hh:4
static ParameterContainer & container()
Definition: io/parameter.hh:190