umfpack_solver.hpp 3.66 KB
Newer Older
1
2
3
4
5
6
7
8
9
10
/*
 * 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_

11
12
13
14
15
16
17
18
19
20
#define UMFPACK_NONE 0

#define SOLVER_NOOPTION 0
#define SOLVER_PRINT_RESIDUAL_NORM_INFINITY 1
#define SOLVER_PRINT_DETERMINANT 2

#ifdef HAVE_EIGEN

/////// Compiled with EIGEN support

21
22
#include "Vector/Vector.hpp"
#include "Eigen/UmfPackSupport"
23
24
#include <Eigen/SparseLU>

25
26
27
28
29
30

template<typename T>
class umfpack_solver
{
public:

31
	template<unsigned int impl, typename id_type> static Vector<T> solve(const SparseMatrix<T,id_type,impl> & A, const Vector<T> & b)
32
	{
33
		std::cerr << "Error Umfpack only support double precision, and int ad id type" << "/n";
34
	}
35
36
37
38
39

	void best_solve()
	{
		std::cerr << "Error Umfpack only support double precision, and int ad id type" << "/n";
	}
40
41
};

42

43
44
45
template<>
class umfpack_solver<double>
{
46

47
48
public:

Pietro Incardona's avatar
Pietro Incardona committed
49
50
51
52
53
54
	/*! \brief No nothing
	 *
	 *
	 */
	void best_solve() {};

55
56
57
58
59
60
61
62
63
	/*! \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
	 *
	 */
Pietro Incardona's avatar
Pietro Incardona committed
64
	static Vector<double,EIGEN_BASE> solve(SparseMatrix<double,int,EIGEN_BASE> & A, const Vector<double,EIGEN_BASE> & b, size_t opt = UMFPACK_NONE)
65
	{
66
		Vcluster & vcl = create_vcluster();
67

68
		Vector<double> x;
69

70
		// only master processor solve
71
		Eigen::UmfPackLU<Eigen::SparseMatrix<double,0,int> > solver;
72

73
74
75
76
77
78
79
		// Collect the matrix on master
		auto mat_ei = A.getMat();

		Eigen::Matrix<double, Eigen::Dynamic, 1> x_ei;

		// Collect the vector on master
		auto b_ei = b.getVec();
80

81
82
83
84
85
86
		// Copy b into x, this also copy the information on how to scatter back the information on x
		x = b;

		if (vcl.getProcessUnitID() == 0)
		{
			solver.compute(mat_ei);
87

88
89
90
91
			if(solver.info()!=Eigen::Success)
			{
				// decomposition failed
				std::cout << __FILE__ << ":" << __LINE__ << " solver failed" << "\n";
92
93
94

				x.scatter();

95
96
				return x;
			}
97

98
			x_ei = solver.solve(b_ei);
99

100
101
			if (opt & SOLVER_PRINT_RESIDUAL_NORM_INFINITY)
			{
102
103
				Eigen::Matrix<double, Eigen::Dynamic, 1> res;
				res = mat_ei * x_ei - b_ei;
104

105
				std::cout << "Infinity norm: " << res.lpNorm<Eigen::Infinity>() << "\n";
106
			}
107

108
109
110
111
			if (opt & SOLVER_PRINT_DETERMINANT)
			{
				std::cout << " Determinant: " << solver.determinant() << "\n";
			}
112

113
			x = x_ei;
114
		}
115
116
117
118

		// Vector is only on master, scatter back the information
		x.scatter();

119
120
121
122
		return x;
	}
};

123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
#else

/////// Compiled without EIGEN support

#include "Vector/Vector.hpp"

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 << __FILE__ << ":" << __LINE__ << " Error Umfpack only support double precision" << "/n";
	}
138
139
140
141
142

	void best_solve()
	{
		std::cerr << __FILE__ << ":" << __LINE__ << " Error Umfpack only support double precision" << "/n";
	}
143
144
145
146
147
148
149
150
151
};


template<>
class umfpack_solver<double>
{

public:

152
	template<unsigned int impl, typename id_type> static Vector<double> solve(SparseMatrix<double,id_type,impl> & A, const Vector<double> & b, size_t opt = UMFPACK_NONE)
153
154
155
156
157
158
159
	{
		std::cerr << __FILE__ << ":" << __LINE__ << " Error in order to use umfpack you must compile OpenFPM with linear algebra support" << "/n";

		Vector<double> x;

		return x;
	}
160
161
162
163
164

	void best_solve()
	{
		std::cerr << __FILE__ << ":" << __LINE__ << " Error in order to use umfpack you must compile OpenFPM with linear algebra support" << "/n";
	}
165
166
167
};

#endif
168
169
170


#endif /* OPENFPM_NUMERICS_SRC_SOLVERS_UMFPACK_SOLVER_HPP_ */