1 #ifndef DUNE_FEM_UMFPACKSOLVER_HH 2 #define DUNE_FEM_UMFPACKSOLVER_HH 19 #include <dune/istl/umfpack.hh> 44 template<
class DF,
class Op,
bool symmetric=false>
45 class UMFPACKOp:
public Operator<DF, DF>
48 typedef DF DiscreteFunctionType;
49 typedef Op OperatorType;
52 typedef ColCompMatrix<typename OperatorType::MatrixType::MatrixBaseType> CCSMatrixType;
53 typedef typename DiscreteFunctionType::DofType DofType;
54 typedef typename DiscreteFunctionType::DiscreteFunctionSpaceType DiscreteFunctionSpaceType;
63 UMFPACKOp(
const OperatorType& op,
const double& redEps,
const double& absLimit,
const int& maxIter,
const bool& verbose,
65 op_(op), verbose_(verbose), ccsmat_(), isloaded_(false)
67 Caller::defaults(UMF_Control);
68 UMF_Control[UMFPACK_PRL] = 4;
77 UMFPACKOp(
const OperatorType& op,
const double& redEps,
const double& absLimit,
79 op_(op), verbose_(parameter.getValue<bool>(
"fem.solver.verbose",false)), ccsmat_(), isloaded_(false)
81 Caller::defaults(UMF_Control);
82 UMF_Control[UMFPACK_PRL] = 4;
86 op_(op), verbose_(parameter.getValue<bool>(
"fem.solver.verbose",false)), ccsmat_(), isloaded_(false)
88 Caller::defaults(UMF_Control);
89 UMF_Control[UMFPACK_PRL] = 4;
100 void operator()(
const DiscreteFunctionType& arg, DiscreteFunctionType& dest)
const 108 template<
typename... A>
109 void prepare(A... )
const 113 ccsmat_ = op_.systemMatrix().matrix();
120 inline void finalize()
const 125 Caller::free_symbolic(&UMF_Symbolic);
126 Caller::free_numeric(&UMF_Numeric);
137 void apply(
const DofType*& arg, DofType*& dest)
const 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);
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;
160 void apply(
const AdaptiveDiscreteFunction<DiscreteFunctionSpaceType>& arg,
161 AdaptiveDiscreteFunction<DiscreteFunctionSpaceType>& dest)
const 163 const DofType* argPtr(arg.leakPointer());
164 DofType* destPtr(dest.leakPointer());
165 apply(argPtr,destPtr);
174 void apply(
const ISTLBlockVectorDiscreteFunction<DiscreteFunctionSpaceType>& arg,
175 ISTLBlockVectorDiscreteFunction<DiscreteFunctionSpaceType>& dest)
const 178 std::vector<DofType> vecArg(arg.size());
179 std::copy(arg.dbegin(),arg.dend(),vecArg.begin());
180 std::vector<DofType> vecDest(dest.size());
182 const DofType* argDataPtr(vecArg.data());
183 DofType* destDataPtr(vecDest.data());
184 apply(argDataPtr,destDataPtr);
186 std::copy(vecDest.begin(),vecDest.end(),dest.dbegin());
189 inline void printTexInfo(std::ostream& out)
const 191 out<<
"Solver: UMFPACK direct solver";
195 inline double averageCommTime()
const 200 inline int iterations()
const 208 inline CCSMatrixType& getCCSMatrix()
214 inline void printDecompositionInfo()
const 216 Caller::report_info(UMF_Control,UMF_Decomposition_Info);
220 const OperatorType& op_;
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];
229 typedef typename Dune::UMFPackMethodChooser<DofType> Caller;
232 void decompose()
const 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);
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;
258 #endif // #if HAVE_DUNE_ISTL 260 #endif // #if HAVE_UMFPACK 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