OPAL (Object Oriented Parallel Accelerator Library) 2022.1
OPAL
Amesos2BottomSolver.hpp
Go to the documentation of this file.
1//
2// Class Amesos2BottomSolver
3// Interface to Amesos2 solvers of the Trilinos package.
4//
5// Copyright (c) 2017 - 2020, Matthias Frey, Paul Scherrer Institut, Villigen PSI, Switzerland
6// All rights reserved
7//
8// Implemented as part of the PhD thesis
9// "Precise Simulations of Multibunches in High Intensity Cyclotrons"
10//
11// This file is part of OPAL.
12//
13// OPAL is free software: you can redistribute it and/or modify
14// it under the terms of the GNU General Public License as published by
15// the Free Software Foundation, either version 3 of the License, or
16// (at your option) any later version.
17//
18// You should have received a copy of the GNU General Public License
19// along with OPAL. If not, see <https://www.gnu.org/licenses/>.
20//
21template <class Level>
23 : BottomSolver<Teuchos::RCP<amr::matrix_t>,
24 Teuchos::RCP<amr::multivector_t>,
25 Level>()
26 , solvertype_m(solvertype)
27 , solver_mp(Teuchos::null)
28{ }
29
30
31template <class Level>
32void Amesos2BottomSolver<Level>::solve(const Teuchos::RCP<mv_t>& x,
33 const Teuchos::RCP<mv_t>& b)
34{
35 /*
36 * solve linear system Ax = b
37 */
38 solver_mp->solve(x.get(), b.get());
39}
40
41
42template <class Level>
43void Amesos2BottomSolver<Level>::setOperator(const Teuchos::RCP<matrix_t>& A,
44 Level* /*level_p*/)
45{
46 try {
47 solver_mp = Amesos2::create<matrix_t, mv_t>(solvertype_m, A);
48 } catch(const std::invalid_argument& ex) {
49 *gmsg << ex.what() << endl;
50 }
51
52 solver_mp->symbolicFactorization();
53 solver_mp->numericFactorization();
54
55 this->isInitialized_m = true;
56}
57
58
59template <class Level>
61 return 1; // direct solvers do only one step
62}
Inform & endl(Inform &inf)
Definition: Inform.cpp:42
Some AMR types used a lot.
Definition: AmrDefs.h:33
Tpetra::MultiVector< scalar_t, local_ordinal_t, global_ordinal_t, node_t > multivector_t
void setOperator(const Teuchos::RCP< matrix_t > &A, Level *level_p=nullptr)
void solve(const Teuchos::RCP< mv_t > &x, const Teuchos::RCP< mv_t > &b)
Amesos2BottomSolver(std::string solvertype="klu2")
Inform * gmsg
Definition: Main.cpp:61