umfpack_solver.hpp 2.08 KB
Newer Older
1
2
3
4
5
6
7
8
9
10
11
12
/*
 * Umfpack_solver.hpp
 *
 *  Created on: Nov 27, 2015
 *      Author: i-bird
 */

#ifndef OPENFPM_NUMERICS_SRC_SOLVERS_UMFPACK_SOLVER_HPP_
#define OPENFPM_NUMERICS_SRC_SOLVERS_UMFPACK_SOLVER_HPP_

#include "Vector/Vector.hpp"
#include "Eigen/UmfPackSupport"
13
14
#include <Eigen/SparseLU>

15
#define UMFPACK_NONE 0
16
17
18
19
20
21
22
23
24
25
26
27

template<typename T>
class umfpack_solver
{
public:

	template<typename impl> static Vector<T> solve(const SparseMatrix<T,impl> & A, const Vector<T> & b)
	{
		std::cerr << "Error Umfpack only suppor double precision" << "/n";
	}
};

28
29
30
31
#define SOLVER_NOOPTION 0
#define SOLVER_PRINT_RESIDUAL_NORM_INFINITY 1
#define SOLVER_PRINT_DETERMINANT 2

32
33
34
template<>
class umfpack_solver<double>
{
35

36
37
public:

38
39
40
41
42
43
44
45
46
47
	/*! \brief Here we invert the matrix and solve the system
	 *
	 *  \warning umfpack is not a parallel solver, this function work only with one processor
	 *
	 *  \note if you want to use umfpack in a NON parallel, but on a distributed data, use solve with triplet
	 *
	 *	\tparam impl Implementation of the SparseMatrix
	 *
	 */
	template<typename impl> static Vector<double> solve(const SparseMatrix<double,int,impl> & A, const Vector<double> & b, size_t opt = UMFPACK_NONE)
48
	{
49
50
		Vcluster & v_cl = *global_v_cluster;

51
		Vector<double> x;
52
		// only master processor solve
53

54
55
56
		if (v_cl.getProcessUnitID() == 0)
		{
			Eigen::UmfPackLU<Eigen::SparseMatrix<double,0,int> > solver;
57

58
			solver.compute(A.getMat());
59

60
61
62
63
64
65
			if(solver.info()!=Eigen::Success)
			{
				// decomposition failed
				std::cout << __FILE__ << ":" << __LINE__ << " solver failed" << "\n";
				return x;
			}
66

67
			x.getVec() = solver.solve(b.getVec());
68

69
70
71
72
			if (opt & SOLVER_PRINT_RESIDUAL_NORM_INFINITY)
			{
				Vector<double> res;
				res.getVec() = A.getMat() * x.getVec() - b.getVec();
73

74
75
				std::cout << "Infinity norm: " << res.getVec().lpNorm<Eigen::Infinity>() << "\n";
			}
76

77
78
79
80
			if (opt & SOLVER_PRINT_DETERMINANT)
			{
				std::cout << " Determinant: " << solver.determinant() << "\n";
			}
81

82
83
84
			// Vector is only on master, scatter back the information
			x.sync();
		}
85
86
87
88
89
90
91
		return x;
	}
};



#endif /* OPENFPM_NUMERICS_SRC_SOLVERS_UMFPACK_SOLVER_HPP_ */