24 #ifndef BOOST_UBLAS_NDEBUG
25 #define BOOST_UBLAS_NDEBUG
29 #include <boost/numeric/ublas/matrix_sparse.hpp>
30 #include <boost/numeric/ublas/io.hpp>
31 #include <boost/numeric/ublas/operation_sparse.hpp>
33 #define VIENNACL_WITH_UBLAS 1
62 using namespace boost::numeric;
64 #define BENCHMARK_RUNS 1
67 inline void printOps(
double num_ops,
double exec_time)
69 std::cout <<
"GFLOPs: " << num_ops / (1000000 * exec_time * 1000) << std::endl;
73 template<
typename ScalarType>
76 ublas::vector<ScalarType> v2_cpu(v2.
size());
79 for (
unsigned int i=0;i<v1.size(); ++i)
81 if (
std::max( fabs(v2_cpu[i]), fabs(v1[i]) ) > 0 )
82 v2_cpu[i] = fabs(v2_cpu[i] - v1[i]) /
std::max( fabs(v2_cpu[i]), fabs(v1[i]) );
90 template<
typename ScalarType>
93 ublas::vector<ScalarType> v2_cpu(v2.
size());
100 template<
typename MatrixType,
typename VectorType,
typename SolverTag,
typename PrecondTag>
101 void run_solver(MatrixType
const & matrix, VectorType
const & rhs, VectorType
const & ref_result, SolverTag
const & solver, PrecondTag
const & precond,
long ops)
104 VectorType result(rhs);
105 VectorType residual(rhs);
114 double exec_time = timer.
get();
115 std::cout <<
"Exec. time: " << exec_time << std::endl;
116 std::cout <<
"Est. ";
printOps(static_cast<double>(ops), exec_time / BENCHMARK_RUNS);
119 std::cout <<
"Estimated rel. residual: " << solver.error() << std::endl;
120 std::cout <<
"Iterations: " << solver.iters() << std::endl;
121 result -= ref_result;
126 template<
typename ScalarType>
132 ublas::vector<ScalarType> ublas_vec1;
133 ublas::vector<ScalarType> ublas_vec2;
134 ublas::vector<ScalarType> ublas_result;
135 unsigned int solver_iters = 100;
136 unsigned int solver_krylov_dim = 20;
137 double solver_tolerance = 1e-6;
139 ublas::compressed_matrix<ScalarType> ublas_matrix;
142 std::cout <<
"Error reading Matrix file" << std::endl;
145 std::cout <<
"done reading matrix" << std::endl;
147 ublas_result = ublas::scalar_vector<ScalarType>(ublas_matrix.size1(),
ScalarType(1.0));
148 ublas_vec1 =
ublas::prod(ublas_matrix, ublas_result);
149 ublas_vec2 = ublas_vec1;
173 std::cout <<
"------- Jacobi preconditioner ----------" << std::endl;
178 std::cout <<
"------- Row-Scaling preconditioner ----------" << std::endl;
186 std::cout <<
"------- ICHOL0 on CPU (ublas) ----------" << std::endl;
190 exec_time = timer.
get();
191 std::cout <<
"Setup time: " << exec_time << std::endl;
195 ublas_ichol0.
apply(ublas_vec1);
196 exec_time = timer.
get();
197 std::cout <<
"ublas time: " << exec_time << std::endl;
199 std::cout <<
"------- ICHOL0 with ViennaCL ----------" << std::endl;
203 exec_time = timer.
get();
204 std::cout <<
"Setup time: " << exec_time << std::endl;
209 vcl_ichol0.
apply(vcl_vec1);
211 exec_time = timer.
get();
212 std::cout <<
"ViennaCL time: " << exec_time << std::endl;
214 std::cout <<
"------- Chow-Patel parallel ICC with ViennaCL ----------" << std::endl;
219 std::cout <<
"Setup time: " << timer.
get() << std::endl;
223 vcl_chow_patel_icc.apply(vcl_vec1);
225 std::cout <<
"ViennaCL Chow-Patel-ICC substitution time: " << timer.
get() << std::endl;
231 std::cout <<
"------- ILU0 on with ublas ----------" << std::endl;
235 exec_time = timer.
get();
236 std::cout <<
"Setup time (no level scheduling): " << exec_time << std::endl;
239 ublas_ilu0.
apply(ublas_vec1);
240 exec_time = timer.
get();
241 std::cout <<
"ublas ILU0 substitution time (no level scheduling): " << exec_time << std::endl;
243 std::cout <<
"------- ILU0 with ViennaCL ----------" << std::endl;
247 exec_time = timer.
get();
248 std::cout <<
"Setup time (no level scheduling): " << exec_time << std::endl;
253 vcl_ilu0.
apply(vcl_vec1);
255 exec_time = timer.
get();
256 std::cout <<
"ViennaCL ILU0 substitution time (no level scheduling): " << exec_time << std::endl;
261 exec_time = timer.
get();
262 std::cout <<
"Setup time (with level scheduling): " << exec_time << std::endl;
267 vcl_ilu0_level_scheduling.
apply(vcl_vec1);
269 exec_time = timer.
get();
270 std::cout <<
"ViennaCL ILU0 substitution time (with level scheduling): " << exec_time << std::endl;
276 std::cout <<
"------- Block-ILU0 with ublas ----------" << std::endl;
278 ublas_vec1 = ublas_vec2;
284 exec_time = timer.
get();
285 std::cout <<
"Setup time: " << exec_time << std::endl;
289 ublas_block_ilu0.
apply(ublas_vec1);
290 exec_time = timer.
get();
291 std::cout <<
"ublas time: " << exec_time << std::endl;
293 std::cout <<
"------- Block-ILU0 with ViennaCL ----------" << std::endl;
298 exec_time = timer.
get();
299 std::cout <<
"Setup time: " << exec_time << std::endl;
305 vcl_block_ilu0.
apply(vcl_vec1);
307 exec_time = timer.
get();
308 std::cout <<
"ViennaCL time: " << exec_time << std::endl;
312 std::cout <<
"------- ILUT with ublas ----------" << std::endl;
314 ublas_vec1 = ublas_vec2;
319 exec_time = timer.
get();
320 std::cout <<
"Setup time (no level scheduling): " << exec_time << std::endl;
323 ublas_ilut.
apply(ublas_vec1);
324 exec_time = timer.
get();
325 std::cout <<
"ublas ILUT substitution time (no level scheduling): " << exec_time << std::endl;
328 std::cout <<
"------- ILUT with ViennaCL ----------" << std::endl;
332 exec_time = timer.
get();
333 std::cout <<
"Setup time (no level scheduling): " << exec_time << std::endl;
338 vcl_ilut.
apply(vcl_vec1);
340 exec_time = timer.
get();
341 std::cout <<
"ViennaCL ILUT substitution time (no level scheduling): " << exec_time << std::endl;
346 exec_time = timer.
get();
347 std::cout <<
"Setup time (with level scheduling): " << exec_time << std::endl;
352 vcl_ilut_level_scheduling.
apply(vcl_vec1);
354 exec_time = timer.
get();
355 std::cout <<
"ViennaCL ILUT substitution time (with level scheduling): " << exec_time << std::endl;
360 std::cout <<
"------- Block-ILUT with ublas ----------" << std::endl;
362 ublas_vec1 = ublas_vec2;
368 exec_time = timer.
get();
369 std::cout <<
"Setup time: " << exec_time << std::endl;
374 ublas_block_ilut.
apply(ublas_vec1);
375 exec_time = timer.
get();
376 std::cout <<
"ublas time: " << exec_time << std::endl;
378 std::cout <<
"------- Block-ILUT with ViennaCL ----------" << std::endl;
383 exec_time = timer.
get();
384 std::cout <<
"Setup time: " << exec_time << std::endl;
390 vcl_block_ilut.
apply(vcl_vec1);
392 exec_time = timer.
get();
393 std::cout <<
"ViennaCL time: " << exec_time << std::endl;
395 std::cout <<
"------- Chow-Patel parallel ILU with ViennaCL ----------" << std::endl;
400 std::cout <<
"Setup time: " << timer.
get() << std::endl;
404 vcl_chow_patel_ilu.apply(vcl_vec1);
406 std::cout <<
"ViennaCL Chow-Patel-ILU substitution time: " << timer.
get() << std::endl;
411 long cg_ops =
static_cast<long>(solver_iters * (ublas_matrix.nnz() + 6 * ublas_vec2.size()));
415 std::cout <<
"------- CG solver (no preconditioner) using ublas ----------" << std::endl;
418 std::cout <<
"------- CG solver (no preconditioner) via ViennaCL, compressed_matrix ----------" << std::endl;
421 bool is_double = (
sizeof(
ScalarType) ==
sizeof(
double));
424 std::cout <<
"------- CG solver, mixed precision (no preconditioner) via ViennaCL, compressed_matrix ----------" << std::endl;
430 std::cout <<
"------- CG solver (no preconditioner) via ViennaCL, coordinate_matrix ----------" << std::endl;
433 std::cout <<
"------- CG solver (no preconditioner) via ViennaCL, ell_matrix ----------" << std::endl;
436 std::cout <<
"------- CG solver (no preconditioner) via ViennaCL, sliced_ell_matrix ----------" << std::endl;
439 std::cout <<
"------- CG solver (no preconditioner) via ViennaCL, hyb_matrix ----------" << std::endl;
442 std::cout <<
"------- CG solver (ICHOL0 preconditioner) using ublas ----------" << std::endl;
443 run_solver(ublas_matrix, ublas_vec2, ublas_result, cg_solver, ublas_ichol0, cg_ops);
445 std::cout <<
"------- CG solver (ICHOL0 preconditioner) via ViennaCL, compressed_matrix ----------" << std::endl;
446 run_solver(vcl_compressed_matrix, vcl_vec2, vcl_result, cg_solver, vcl_ichol0, cg_ops);
448 std::cout <<
"------- CG solver (Chow-Patel ICHOL0 preconditioner) via ViennaCL, compressed_matrix ----------" << std::endl;
449 run_solver(vcl_compressed_matrix, vcl_vec2, vcl_result, cg_solver, vcl_chow_patel_icc, cg_ops);
451 std::cout <<
"------- CG solver (ILU0 preconditioner) using ublas ----------" << std::endl;
452 run_solver(ublas_matrix, ublas_vec2, ublas_result, cg_solver, ublas_ilu0, cg_ops);
454 std::cout <<
"------- CG solver (ILU0 preconditioner) via ViennaCL, compressed_matrix ----------" << std::endl;
455 run_solver(vcl_compressed_matrix, vcl_vec2, vcl_result, cg_solver, vcl_ilu0, cg_ops);
458 std::cout <<
"------- CG solver (Block-ILU0 preconditioner) using ublas ----------" << std::endl;
459 run_solver(ublas_matrix, ublas_vec2, ublas_result, cg_solver, ublas_block_ilu0, cg_ops);
461 std::cout <<
"------- CG solver (Block-ILU0 preconditioner) via ViennaCL, compressed_matrix ----------" << std::endl;
462 run_solver(vcl_compressed_matrix, vcl_vec2, vcl_result, cg_solver, vcl_block_ilu0, cg_ops);
464 std::cout <<
"------- CG solver (ILUT preconditioner) using ublas ----------" << std::endl;
465 run_solver(ublas_matrix, ublas_vec2, ublas_result, cg_solver, ublas_ilut, cg_ops);
467 std::cout <<
"------- CG solver (ILUT preconditioner) via ViennaCL, compressed_matrix ----------" << std::endl;
468 run_solver(vcl_compressed_matrix, vcl_vec2, vcl_result, cg_solver, vcl_ilut, cg_ops);
470 std::cout <<
"------- CG solver (ILUT preconditioner) via ViennaCL, coordinate_matrix ----------" << std::endl;
471 run_solver(vcl_coordinate_matrix, vcl_vec2, vcl_result, cg_solver, vcl_ilut, cg_ops);
473 std::cout <<
"------- CG solver (Block-ILUT preconditioner) using ublas ----------" << std::endl;
474 run_solver(ublas_matrix, ublas_vec2, ublas_result, cg_solver, ublas_block_ilut, cg_ops);
476 std::cout <<
"------- CG solver (Block-ILUT preconditioner) via ViennaCL, compressed_matrix ----------" << std::endl;
477 run_solver(vcl_compressed_matrix, vcl_vec2, vcl_result, cg_solver, vcl_block_ilut, cg_ops);
479 std::cout <<
"------- CG solver (Jacobi preconditioner) using ublas ----------" << std::endl;
480 run_solver(ublas_matrix, ublas_vec2, ublas_result, cg_solver, ublas_jacobi, cg_ops);
482 std::cout <<
"------- CG solver (Jacobi preconditioner) via ViennaCL, compressed_matrix ----------" << std::endl;
483 run_solver(vcl_compressed_matrix, vcl_vec2, vcl_result, cg_solver, vcl_jacobi_csr, cg_ops);
485 std::cout <<
"------- CG solver (Jacobi preconditioner) via ViennaCL, coordinate_matrix ----------" << std::endl;
486 run_solver(vcl_coordinate_matrix, vcl_vec2, vcl_result, cg_solver, vcl_jacobi_coo, cg_ops);
489 std::cout <<
"------- CG solver (row scaling preconditioner) using ublas ----------" << std::endl;
490 run_solver(ublas_matrix, ublas_vec2, ublas_result, cg_solver, ublas_row_scaling, cg_ops);
492 std::cout <<
"------- CG solver (row scaling preconditioner) via ViennaCL, compressed_matrix ----------" << std::endl;
493 run_solver(vcl_compressed_matrix, vcl_vec2, vcl_result, cg_solver, vcl_row_scaling_csr, cg_ops);
495 std::cout <<
"------- CG solver (row scaling preconditioner) via ViennaCL, coordinate_matrix ----------" << std::endl;
496 run_solver(vcl_coordinate_matrix, vcl_vec2, vcl_result, cg_solver, vcl_row_scaling_coo, cg_ops);
502 long bicgstab_ops =
static_cast<long>(solver_iters * (2 * ublas_matrix.nnz() + 13 * ublas_vec2.size()));
506 std::cout <<
"------- BiCGStab solver (no preconditioner) using ublas ----------" << std::endl;
509 std::cout <<
"------- BiCGStab solver (no preconditioner) via ViennaCL, compressed_matrix ----------" << std::endl;
512 std::cout <<
"------- BiCGStab solver (ILU0 preconditioner) via ViennaCL, compressed_matrix ----------" << std::endl;
513 run_solver(vcl_compressed_matrix, vcl_vec2, vcl_result, bicgstab_solver, vcl_ilu0, bicgstab_ops);
515 std::cout <<
"------- BiCGStab solver (Chow-Patel-ILU preconditioner) via ViennaCL, compressed_matrix ----------" << std::endl;
516 run_solver(vcl_compressed_matrix, vcl_vec2, vcl_result, bicgstab_solver, vcl_chow_patel_ilu, bicgstab_ops);
518 std::cout <<
"------- BiCGStab solver (no preconditioner) via ViennaCL, coordinate_matrix ----------" << std::endl;
521 std::cout <<
"------- BiCGStab solver (no preconditioner) via ViennaCL, ell_matrix ----------" << std::endl;
524 std::cout <<
"------- BiCGStab solver (no preconditioner) via ViennaCL, sliced_ell_matrix ----------" << std::endl;
527 std::cout <<
"------- BiCGStab solver (no preconditioner) via ViennaCL, hyb_matrix ----------" << std::endl;
530 std::cout <<
"------- BiCGStab solver (ILUT preconditioner) using ublas ----------" << std::endl;
531 run_solver(ublas_matrix, ublas_vec2, ublas_result, bicgstab_solver, ublas_ilut, bicgstab_ops);
533 std::cout <<
"------- BiCGStab solver (ILUT preconditioner) via ViennaCL, compressed_matrix ----------" << std::endl;
534 run_solver(vcl_compressed_matrix, vcl_vec2, vcl_result, bicgstab_solver, vcl_ilut, bicgstab_ops);
536 std::cout <<
"------- BiCGStab solver (Block-ILUT preconditioner) using ublas ----------" << std::endl;
537 run_solver(ublas_matrix, ublas_vec2, ublas_result, bicgstab_solver, ublas_block_ilut, bicgstab_ops);
539 #ifdef VIENNACL_WITH_OPENCL
540 std::cout <<
"------- BiCGStab solver (Block-ILUT preconditioner) via ViennaCL, compressed_matrix ----------" << std::endl;
541 run_solver(vcl_compressed_matrix, vcl_vec2, vcl_result, bicgstab_solver, vcl_block_ilut, bicgstab_ops);
547 std::cout <<
"------- BiCGStab solver (Jacobi preconditioner) using ublas ----------" << std::endl;
548 run_solver(ublas_matrix, ublas_vec2, ublas_result, bicgstab_solver, ublas_jacobi, bicgstab_ops);
550 std::cout <<
"------- BiCGStab solver (Jacobi preconditioner) via ViennaCL, compressed_matrix ----------" << std::endl;
551 run_solver(vcl_compressed_matrix, vcl_vec2, vcl_result, bicgstab_solver, vcl_jacobi_csr, bicgstab_ops);
553 std::cout <<
"------- BiCGStab solver (Jacobi preconditioner) via ViennaCL, coordinate_matrix ----------" << std::endl;
554 run_solver(vcl_coordinate_matrix, vcl_vec2, vcl_result, bicgstab_solver, vcl_jacobi_coo, bicgstab_ops);
557 std::cout <<
"------- BiCGStab solver (row scaling preconditioner) using ublas ----------" << std::endl;
558 run_solver(ublas_matrix, ublas_vec2, ublas_result, bicgstab_solver, ublas_row_scaling, bicgstab_ops);
560 std::cout <<
"------- BiCGStab solver (row scaling preconditioner) via ViennaCL, compressed_matrix ----------" << std::endl;
561 run_solver(vcl_compressed_matrix, vcl_vec2, vcl_result, bicgstab_solver, vcl_row_scaling_csr, bicgstab_ops);
563 std::cout <<
"------- BiCGStab solver (row scaling preconditioner) via ViennaCL, coordinate_matrix ----------" << std::endl;
564 run_solver(vcl_coordinate_matrix, vcl_vec2, vcl_result, bicgstab_solver, vcl_row_scaling_coo, bicgstab_ops);
571 long gmres_ops =
static_cast<long>(solver_iters * (ublas_matrix.nnz() + (solver_iters * 2 + 7) * ublas_vec2.size()));
575 std::cout <<
"------- GMRES solver (no preconditioner) using ublas ----------" << std::endl;
578 std::cout <<
"------- GMRES solver (no preconditioner) via ViennaCL, compressed_matrix ----------" << std::endl;
581 std::cout <<
"------- GMRES solver (no preconditioner) on GPU, coordinate_matrix ----------" << std::endl;
584 std::cout <<
"------- GMRES solver (no preconditioner) on GPU, ell_matrix ----------" << std::endl;
587 std::cout <<
"------- GMRES solver (no preconditioner) on GPU, sliced_ell_matrix ----------" << std::endl;
590 std::cout <<
"------- GMRES solver (no preconditioner) on GPU, hyb_matrix ----------" << std::endl;
593 std::cout <<
"------- GMRES solver (ILUT preconditioner) using ublas ----------" << std::endl;
594 run_solver(ublas_matrix, ublas_vec2, ublas_result, gmres_solver, ublas_ilut, gmres_ops);
596 std::cout <<
"------- GMRES solver (ILUT preconditioner) via ViennaCL, compressed_matrix ----------" << std::endl;
597 run_solver(vcl_compressed_matrix, vcl_vec2, vcl_result, gmres_solver, vcl_ilut, gmres_ops);
599 std::cout <<
"------- GMRES solver (ILUT preconditioner) via ViennaCL, coordinate_matrix ----------" << std::endl;
600 run_solver(vcl_coordinate_matrix, vcl_vec2, vcl_result, gmres_solver, vcl_ilut, gmres_ops);
603 std::cout <<
"------- GMRES solver (Jacobi preconditioner) using ublas ----------" << std::endl;
604 run_solver(ublas_matrix, ublas_vec2, ublas_result, gmres_solver, ublas_jacobi, gmres_ops);
606 std::cout <<
"------- GMRES solver (Jacobi preconditioner) via ViennaCL, compressed_matrix ----------" << std::endl;
607 run_solver(vcl_compressed_matrix, vcl_vec2, vcl_result, gmres_solver, vcl_jacobi_csr, gmres_ops);
609 std::cout <<
"------- GMRES solver (Jacobi preconditioner) via ViennaCL, coordinate_matrix ----------" << std::endl;
610 run_solver(vcl_coordinate_matrix, vcl_vec2, vcl_result, gmres_solver, vcl_jacobi_coo, gmres_ops);
613 std::cout <<
"------- GMRES solver (row scaling preconditioner) using ublas ----------" << std::endl;
614 run_solver(ublas_matrix, ublas_vec2, ublas_result, gmres_solver, ublas_row_scaling, gmres_ops);
616 std::cout <<
"------- GMRES solver (row scaling preconditioner) via ViennaCL, compressed_matrix ----------" << std::endl;
617 run_solver(vcl_compressed_matrix, vcl_vec2, vcl_result, gmres_solver, vcl_row_scaling_csr, gmres_ops);
619 std::cout <<
"------- GMRES solver (row scaling preconditioner) via ViennaCL, coordinate_matrix ----------" << std::endl;
620 run_solver(vcl_coordinate_matrix, vcl_vec2, vcl_result, gmres_solver, vcl_row_scaling_coo, gmres_ops);
627 std::cout << std::endl;
628 std::cout <<
"----------------------------------------------" << std::endl;
629 std::cout <<
" Device Info" << std::endl;
630 std::cout <<
"----------------------------------------------" << std::endl;
632 #ifdef VIENNACL_WITH_OPENCL
634 std::vector<viennacl::ocl::device>
const & devices = pf.
devices();
640 if (devices.size() > 1)
651 std::cout <<
"---------------------------------------------------------------------------" << std::endl;
652 std::cout <<
"---------------------------------------------------------------------------" << std::endl;
653 std::cout <<
" Benchmark for Execution Times of Iterative Solvers provided with ViennaCL " << std::endl;
654 std::cout <<
"---------------------------------------------------------------------------" << std::endl;
655 std::cout <<
" Note that the purpose of this benchmark is not to run solvers until" << std::endl;
656 std::cout <<
" convergence. Instead, only the execution times of a few iterations are" << std::endl;
657 std::cout <<
" recorded. Residual errors are only printed for information." << std::endl << std::endl;
660 std::cout << std::endl;
661 std::cout <<
"----------------------------------------------" << std::endl;
662 std::cout <<
"----------------------------------------------" << std::endl;
663 std::cout <<
"## Benchmark :: Solver" << std::endl;
664 std::cout <<
"----------------------------------------------" << std::endl;
665 std::cout << std::endl;
666 std::cout <<
" -------------------------------" << std::endl;
667 std::cout <<
" # benchmarking single-precision" << std::endl;
668 std::cout <<
" -------------------------------" << std::endl;
669 run_benchmark<float>(ctx);
670 #ifdef VIENNACL_WITH_OPENCL
674 std::cout << std::endl;
675 std::cout <<
" -------------------------------" << std::endl;
676 std::cout <<
" # benchmarking double-precision" << std::endl;
677 std::cout <<
" -------------------------------" << std::endl;
678 run_benchmark<double>(ctx);
Sparse matrix class using a hybrid format composed of the ELL and CSR format for storing the nonzeros...
A tag for incomplete LU and incomplete Cholesky factorization with static pattern (Parallel-ILU0...
T norm_2(std::vector< T, A > const &v1)
A reader and writer for the matrix market format is implemented here.
std::vector< platform > get_platforms()
ILU0 preconditioner class, can be supplied to solve()-routines.
Incomplete Cholesky preconditioner class with static pattern (ICHOL0), can be supplied to solve()-rou...
int run_benchmark(viennacl::context ctx)
A tag for incomplete Cholesky factorization with static pattern (ILU0)
Jacobi preconditioner class, can be supplied to solve()-routines. Generic version for non-ViennaCL ma...
A tag for incomplete LU factorization with static pattern (ILU0)
The stabilized bi-conjugate gradient method is implemented here.
ScalarType diff_2(ublas::vector< ScalarType > &v1, viennacl::vector< ScalarType > &v2)
Implementations of incomplete Cholesky factorization preconditioners with static nonzero pattern...
ScalarType diff_inf(ublas::vector< ScalarType > &v1, viennacl::vector< ScalarType > &v2)
void finish()
Synchronizes the execution. finish() will only return after all compute kernels (CUDA, OpenCL) have completed.
bool use_level_scheduling() const
A tag for a jacobi preconditioner.
Parallel Chow-Patel ILU preconditioner class, can be supplied to solve()-routines.
A block ILU preconditioner class, can be supplied to solve()-routines.
T max(const T &lhs, const T &rhs)
Maximum.
Implementation of a simple Jacobi preconditioner.
Jacobi-type preconditioner class, can be supplied to solve()-routines. This is a diagonal preconditio...
viennacl::ocl::device const & current_device()
Convenience function for returning the active device in the current context.
A tag for the solver GMRES. Used for supplying solver parameters and for dispatching the solve() func...
Implementation of the coordinate_matrix class.
VectorT solve(MatrixT const &matrix, VectorT const &rhs, bicgstab_tag const &tag, PreconditionerT const &precond)
Implementation of a OpenCL-like context, which serves as a unification of {OpenMP, CUDA, OpenCL} at the user API.
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...
Represents a generic 'context' similar to an OpenCL context, but is backend-agnostic and thus also su...
void apply(VectorT &vec) const
viennacl::vector< float > v1
Implementation of the hyb_matrix class.
VectorT prod(std::vector< std::vector< T, A1 >, A2 > const &matrix, VectorT const &vector)
Implementations of the generalized minimum residual method are in this file.
Sparse matrix class using the ELLPACK format for storing the nonzeros.
iterator begin()
Returns an iterator pointing to the beginning of the vector (STL like)
A tag class representing the use of no preconditioner.
Implementations of incomplete factorization preconditioners. Convenience header file.
A tag for incomplete LU factorization with threshold (ILUT)
Sparse matrix class using the sliced ELLPACK with parameters C, .
Implementation of the compressed_matrix class.
Implementation of the sliced_ell_matrix class.
A tag for a row scaling preconditioner which merely normalizes the equation system such that each row...
void run_solver(MatrixType const &matrix, VectorType const &rhs, VectorType const &ref_result, SolverTag const &solver, PrecondTag const &precond, long ops)
bool double_support() const
ViennaCL convenience function: Returns true if the device supports double precision.
ILUT preconditioner class, can be supplied to solve()-routines.
The conjugate gradient method is implemented here.
void printOps(double num_ops, double exec_time)
void apply(VectorT &vec) const
Implementation of the ell_matrix class.
void apply(VectorT &vec) const
A row normalization preconditioner is implemented here.
A simple, yet (mostly) sufficiently accurate timer for benchmarking and profiling.
void apply(VectorT &vec) const
viennacl::vector< int > v2
bool use_level_scheduling() const
void prod(std::vector< std::map< IndexT, NumericT > > const &stl_A, std::vector< std::map< IndexT, NumericT > > const &stl_B, std::vector< std::map< IndexT, NumericT > > &stl_C)
The vector type with operator-overloads and proxy classes is defined here. Linear algebra operations ...
T norm_inf(std::vector< T, A > const &v1)
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) ...
size_type size() const
Returns the length of the vector (cf. std::vector)
A tag for the conjugate gradient Used for supplying solver parameters and for dispatching the solve()...
A sparse square matrix in compressed sparse rows format.
A tag for the stabilized Bi-conjugate gradient solver. Used for supplying solver parameters and for d...
long read_matrix_market_file(MatrixT &mat, const char *file, long index_base=1)
Reads a sparse matrix from a file (MatrixMarket format)
iterator end()
Returns an iterator pointing to the end of the vector (STL like)
viennacl::ocl::context & get_context(long i)
Convenience function for returning the current context.
A tag for the conjugate gradient Used for supplying solver parameters and for dispatching the solve()...
Implementation of the ViennaCL scalar class.
void setup_context(long i, std::vector< cl_device_id > const &devices)
Convenience function for setting devices for a context.
The conjugate gradient method using mixed precision is implemented here. Experimental.
Parallel Chow-Patel ILU preconditioner class, can be supplied to solve()-routines.
A sparse square matrix, where entries are stored as triplets (i,j, val), where i and j are the row an...