ViennaCL - The Vienna Computing Library  1.6.0
Free open-source GPU-accelerated linear algebra and solver library.
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
amg.cpp
Go to the documentation of this file.
1 /* =========================================================================
2  Copyright (c) 2010-2014, Institute for Microelectronics,
3  Institute for Analysis and Scientific Computing,
4  TU Wien.
5  Portions of this software are copyright by UChicago Argonne, LLC.
6 
7  -----------------
8  ViennaCL - The Vienna Computing Library
9  -----------------
10 
11  Project Head: Karl Rupp rupp@iue.tuwien.ac.at
12 
13  (A list of authors and contributors can be found in the PDF manual)
14 
15  License: MIT (X11), see file LICENSE in the base directory
16 ============================================================================= */
17 
18 
28 #ifndef NDEBUG //without NDEBUG the performance of sparse ublas matrices is poor.
29  #define NDEBUG
30 #endif
31 
32 #include <boost/numeric/ublas/matrix_sparse.hpp>
33 #include <boost/numeric/ublas/operation_sparse.hpp>
34 
35 #define VIENNACL_WITH_UBLAS 1
36 
37 #include "viennacl/vector.hpp"
40 #include "viennacl/linalg/ilu.hpp"
41 #include "viennacl/linalg/cg.hpp"
45 
49 #include "viennacl/linalg/amg.hpp"
50 
54 #include <iostream>
55 #include <vector>
56 #include <ctime>
57 #include "vector-io.hpp"
58 
59 
65 template<typename MatrixType, typename VectorType, typename SolverTag, typename PrecondTag>
66 void run_solver(MatrixType const & matrix, VectorType const & rhs, VectorType const & ref_result, SolverTag const & solver, PrecondTag const & precond)
67 {
68  VectorType result(rhs);
69  VectorType residual(rhs);
70 
71  result = viennacl::linalg::solve(matrix, rhs, solver, precond);
72  residual -= viennacl::linalg::prod(matrix, result);
73  std::cout << " > Relative residual: " << viennacl::linalg::norm_2(residual) / viennacl::linalg::norm_2(rhs) << std::endl;
74  std::cout << " > Iterations: " << solver.iters() << std::endl;
75  result -= ref_result;
76  std::cout << " > Relative deviation from result: " << viennacl::linalg::norm_2(result) / viennacl::linalg::norm_2(ref_result) << std::endl;
77 }
78 
84 template<typename ScalarType>
85 void run_amg(viennacl::linalg::cg_tag & cg_solver,
86  boost::numeric::ublas::vector<ScalarType> & /*ublas_vec*/,
87  boost::numeric::ublas::vector<ScalarType> & /*ublas_result*/,
88  boost::numeric::ublas::compressed_matrix<ScalarType> & ublas_matrix,
90  viennacl::vector<ScalarType> & vcl_result,
91  viennacl::compressed_matrix<ScalarType> & vcl_compressed_matrix,
92  std::string info,
94 {
95 
97  boost::numeric::ublas::vector<ScalarType> avgstencil;
98  unsigned int coarselevels = amg_tag.get_coarselevels();
99 
100  std::cout << "-- CG with AMG preconditioner, " << info << " --" << std::endl;
101 
102  std::cout << " * Setup phase (ublas types)..." << std::endl;
103 
104  // Coarse level measure might have been changed during setup. Reload!
105  ublas_amg.tag().set_coarselevels(coarselevels);
106  ublas_amg.setup();
107 
108  std::cout << " * Operator complexity: " << ublas_amg.calc_complexity(avgstencil) << std::endl;
109 
110  amg_tag.set_coarselevels(coarselevels);
112  std::cout << " * Setup phase (ViennaCL types)..." << std::endl;
113  vcl_amg.tag().set_coarselevels(coarselevels);
114  vcl_amg.setup();
115 
116  std::cout << " * CG solver (ublas types)..." << std::endl;
117  //run_solver(ublas_matrix, ublas_vec, ublas_result, cg_solver, ublas_amg);
118 
119  std::cout << " * CG solver (ViennaCL types)..." << std::endl;
120  run_solver(vcl_compressed_matrix, vcl_vec, vcl_result, cg_solver, vcl_amg);
121 
122 }
123 
129 int main()
130 {
134  std::cout << std::endl;
135  std::cout << "----------------------------------------------" << std::endl;
136  std::cout << " Device Info" << std::endl;
137  std::cout << "----------------------------------------------" << std::endl;
138 
139 #ifdef VIENNACL_WITH_OPENCL
140  // Optional: Customize OpenCL backend
142  std::vector<viennacl::ocl::device> const & devices = pf.devices();
143 
144  // Optional: Set first device to first context:
145  viennacl::ocl::setup_context(0, devices[0]);
146 
147  // Optional: Set second device for second context (use the same device for the second context if only one device available):
148  if (devices.size() > 1)
149  viennacl::ocl::setup_context(1, devices[1]);
150  else
151  viennacl::ocl::setup_context(1, devices[0]);
152 
153  std::cout << viennacl::ocl::current_device().info() << std::endl;
155 #else
156  viennacl::context ctx;
157 #endif
158 
159  typedef float ScalarType; // feel free to change this to double if supported by your device
160 
161 
165  boost::numeric::ublas::vector<ScalarType> ublas_vec, ublas_result;
166  boost::numeric::ublas::compressed_matrix<ScalarType> ublas_matrix;
167 
168  // Read matrix
169  if (!viennacl::io::read_matrix_market_file(ublas_matrix, "../examples/testdata/mat65k.mtx"))
170  {
171  std::cout << "Error reading Matrix file" << std::endl;
172  return EXIT_FAILURE;
173  }
174 
175  // Set up rhs and result vector
176  if (!readVectorFromFile("../examples/testdata/rhs65025.txt", ublas_vec))
177  {
178  std::cout << "Error reading RHS file" << std::endl;
179  return 0;
180  }
181 
182  if (!readVectorFromFile("../examples/testdata/result65025.txt", ublas_result))
183  {
184  std::cout << "Error reading Result file" << std::endl;
185  return 0;
186  }
187 
188  viennacl::vector<ScalarType> vcl_vec(ublas_vec.size(), ctx);
189  viennacl::vector<ScalarType> vcl_result(ublas_vec.size(), ctx);
190  viennacl::compressed_matrix<ScalarType> vcl_compressed_matrix(ublas_vec.size(), ublas_vec.size(), ctx);
191 
192  // Copy to GPU
193  viennacl::copy(ublas_matrix, vcl_compressed_matrix);
194  viennacl::copy(ublas_vec, vcl_vec);
195  viennacl::copy(ublas_result, vcl_result);
196 
200  viennacl::linalg::cg_tag cg_solver;
202 
207  std::cout << "-- CG solver (CPU, no preconditioner) --" << std::endl;
208  run_solver(ublas_matrix, ublas_vec, ublas_result, cg_solver, viennacl::linalg::no_precond());
209 
210  std::cout << "-- CG solver (GPU, no preconditioner) --" << std::endl;
211  run_solver(vcl_compressed_matrix, vcl_vec, vcl_result, cg_solver, viennacl::linalg::no_precond());
212 
216  amg_tag = viennacl::linalg::amg_tag(VIENNACL_AMG_COARSE_RS, // coarsening strategy
217  VIENNACL_AMG_INTERPOL_DIRECT, // interpolation strategy
218  0.25, // strength of dependence threshold
219  0.2, // interpolation weight
220  0.67, // jacobi smoother weight
221  3, // presmoothing steps
222  3, // postsmoothing steps
223  0); // number of coarse levels to be used (0: automatically use as many as reasonable)
224  run_amg (cg_solver, ublas_vec, ublas_result, ublas_matrix, vcl_vec, vcl_result, vcl_compressed_matrix, "RS COARSENING, DIRECT INTERPOLATION", amg_tag);
225 
230  run_amg ( cg_solver, ublas_vec, ublas_result, ublas_matrix, vcl_vec, vcl_result, vcl_compressed_matrix, "RS COARSENING, CLASSIC INTERPOLATION", amg_tag);
231 
236  run_amg (cg_solver, ublas_vec, ublas_result, ublas_matrix, vcl_vec, vcl_result, vcl_compressed_matrix, "ONEPASS COARSENING, DIRECT INTERPOLATION", amg_tag);
237 
242  run_amg (cg_solver, ublas_vec, ublas_result, ublas_matrix, vcl_vec, vcl_result, vcl_compressed_matrix, "RS0 COARSENING, DIRECT INTERPOLATION", amg_tag);
243 
248  run_amg (cg_solver, ublas_vec, ublas_result, ublas_matrix, vcl_vec, vcl_result, vcl_compressed_matrix, "RS3 COARSENING, DIRECT INTERPOLATION", amg_tag);
249 
254  run_amg (cg_solver, ublas_vec, ublas_result, ublas_matrix, vcl_vec, vcl_result, vcl_compressed_matrix, "AG COARSENING, AG INTERPOLATION", amg_tag);
255 
259  amg_tag = viennacl::linalg::amg_tag(VIENNACL_AMG_COARSE_AG, VIENNACL_AMG_INTERPOL_SA, 0.08, 0.67, 0.67, 3, 3, 0);
260  run_amg (cg_solver, ublas_vec, ublas_result, ublas_matrix, vcl_vec, vcl_result, vcl_compressed_matrix, "AG COARSENING, SA INTERPOLATION",amg_tag);
261 
262 
266  std::cout << "!!!! TUTORIAL COMPLETED SUCCESSFULLY !!!!" << std::endl;
267 
268  return EXIT_SUCCESS;
269 }
270 
#define VIENNACL_AMG_INTERPOL_SA
Definition: amg_base.hpp:49
unsigned int get_coarselevels() const
Definition: amg_base.hpp:111
T norm_2(std::vector< T, A > const &v1)
Definition: norm_2.hpp:86
A reader and writer for the matrix market format is implemented here.
std::vector< platform > get_platforms()
Definition: platform.hpp:124
#define VIENNACL_AMG_INTERPOL_AG
Definition: amg_base.hpp:48
Generic interface for the l^2-norm. See viennacl/linalg/vector_operations.hpp for implementations...
#define VIENNACL_AMG_INTERPOL_DIRECT
Definition: amg_base.hpp:46
#define VIENNACL_AMG_COARSE_RS
Definition: amg_base.hpp:41
Wrapper class for an OpenCL platform.
Definition: platform.hpp:45
int main()
Definition: bisect.cpp:160
std::vector< device > devices(cl_device_type dtype=CL_DEVICE_TYPE_DEFAULT)
Returns the available devices of the supplied device type.
Definition: platform.hpp:91
#define VIENNACL_AMG_COARSE_RS3
Definition: amg_base.hpp:44
The stabilized bi-conjugate gradient method is implemented here.
viennacl::ocl::device const & current_device()
Convenience function for returning the active device in the current context.
Definition: backend.hpp:351
Implementation of the coordinate_matrix class.
#define VIENNACL_AMG_COARSE_RS0
Definition: amg_base.hpp:43
#define VIENNACL_AMG_COARSE_AG
Definition: amg_base.hpp:45
std::string info(vcl_size_t indent=0, char indent_char= ' ') const
Returns an info string with a few properties of the device. Use full_info() to get all details...
Definition: device.hpp:995
Represents a generic 'context' similar to an OpenCL context, but is backend-agnostic and thus also su...
Definition: context.hpp:39
VectorT prod(std::vector< std::vector< T, A1 >, A2 > const &matrix, VectorT const &vector)
Definition: prod.hpp:91
A tag for algebraic multigrid (AMG). Used to transport information from the user to the implementatio...
Definition: amg_base.hpp:61
A tag class representing the use of no preconditioner.
Definition: forwards.h:833
Implementations of incomplete factorization preconditioners. Convenience header file.
Implementation of the compressed_matrix class.
bool readVectorFromFile(const std::string &filename, VectorType &vec)
Definition: vector-io.hpp:104
void run_solver(MatrixType const &matrix, VectorType const &rhs, VectorType const &ref_result, SolverTag const &solver, PrecondTag const &precond, long ops)
Definition: solver.cpp:97
Main include file for algebraic multigrid (AMG) preconditioners. Experimental.
AMG preconditioner class, can be supplied to solve()-routines.
Definition: amg.hpp:339
The conjugate gradient method is implemented here.
void setup()
Start setup phase for this class and copy data structures.
Definition: amg.hpp:387
NumericType calc_complexity(VectorT &avgstencil)
Returns complexity measures.
Definition: amg.hpp:417
The vector type with operator-overloads and proxy classes is defined here. Linear algebra operations ...
void copy(std::vector< NumericT > &cpu_vec, circulant_matrix< NumericT, AlignmentV > &gpu_mat)
Copies a circulant matrix from the std::vector to the OpenCL device (either GPU or multi-core CPU) ...
#define VIENNACL_AMG_COARSE_ONEPASS
Definition: amg_base.hpp:42
float ScalarType
Definition: fft_1d.cpp:42
#define VIENNACL_AMG_INTERPOL_CLASSIC
Definition: amg_base.hpp:47
A tag for the conjugate gradient Used for supplying solver parameters and for dispatching the solve()...
Definition: cg.hpp:48
A sparse square matrix in compressed sparse rows format.
long read_matrix_market_file(MatrixT &mat, const char *file, long index_base=1)
Reads a sparse matrix from a file (MatrixMarket format)
detail::amg::amg_tag amg_tag
Definition: amg.hpp:61
viennacl::ocl::context & get_context(long i)
Convenience function for returning the current context.
Definition: backend.hpp:225
void setup_context(long i, std::vector< cl_device_id > const &devices)
Convenience function for setting devices for a context.
Definition: backend.hpp:231
void set_coarselevels(unsigned int coarselevels)
Definition: amg_base.hpp:110
viennacl::vector< NumericT > solve(MatrixT const &A, viennacl::vector_base< NumericT > const &rhs, bicgstab_tag const &tag, viennacl::linalg::no_precond)
Implementation of a pipelined stabilized Bi-conjugate gradient solver.
Definition: bicgstab.hpp:88