OPAL (Object Oriented Parallel Accelerator Library) 2022.1
OPAL
MueLuBottomSolver.h
Go to the documentation of this file.
1//
2// Class MueLuBottomSolver
3// Interface to the SAAMG solver of MueLu.
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//
21#ifndef MUELU_BOTTOM_SOLVER_H
22#define MUELU_BOTTOM_SOLVER_H
23
25
26#include "Amr/AmrDefs.h"
27
28#include "Ippl.h"
29
30#include <MueLu.hpp>
31#include <MueLu_Level.hpp>
32#include <MueLu_Utilities.hpp>
33
34#include <MueLu_ParameterListInterpreter.hpp>
35
36template <class Level>
37class MueLuBottomSolver : public BottomSolver<Teuchos::RCP<amr::matrix_t>,
38 Teuchos::RCP<amr::multivector_t>,
39 Level >
40{
41public:
50
53
54// typedef amr::AmrGeometry_t AmrGeometry_t;
55
56 typedef MueLu::Hierarchy<scalar_t, lo_t, go_t, node_t> hierarchy_t;
57 typedef MueLu::Level level_t;
58 typedef Xpetra::Matrix<scalar_t, lo_t, go_t, node_t> xmatrix_t;
59 typedef Xpetra::MultiVector<scalar_t, lo_t, go_t, node_t> xmv_t;
60 typedef MueLu::Utilities<scalar_t, lo_t, go_t, node_t> util_t;
61
62 typedef MueLu::ParameterListInterpreter<scalar_t, lo_t, go_t, node_t> pListInterpreter_t;
63 typedef MueLu::HierarchyManager<scalar_t, lo_t, go_t, node_t> manager_t;
64
65public:
66 MueLuBottomSolver(const bool& rebalance,
67 const std::string& reuse);
68
69 void solve(const Teuchos::RCP<mv_t>& x,
70 const Teuchos::RCP<mv_t>& b);
71
72 void setOperator(const Teuchos::RCP<matrix_t>& A,
73 Level* level_p = nullptr);
74
75 std::size_t getNumIters();
76
77 /*
78 * MueLu reuse option.
79 * Either NONE, RP, RAP, SYMBOLIC or FULL
80 */
81 static std::string convertToMueLuReuseOption(const std::string& reuse);
82
83private:
84 void initMueLuList_m(const std::string& reuse);
85
86private:
87 Teuchos::RCP<hierarchy_t> hierarchy_mp;
88 Teuchos::RCP<manager_t> factory_mp;
89 Teuchos::RCP<xmatrix_t> A_mp;
90
92
93 Teuchos::ParameterList mueluList_m;
94
96
98};
99
101
102#endif
long global_ordinal_t
Tpetra::Operator< scalar_t, local_ordinal_t, global_ordinal_t, node_t > operator_t
KokkosClassic::DefaultNode::DefaultNodeType node_t
amrex::Box AmrBox_t
Definition: AmrDefs.h:50
amrex::IntVect AmrIntVect_t
Definition: AmrDefs.h:48
int local_ordinal_t
Tpetra::Vector< scalar_t, local_ordinal_t, global_ordinal_t, node_t > vector_t
double scalar_t
Tpetra::MultiVector< scalar_t, local_ordinal_t, global_ordinal_t, node_t > multivector_t
Tpetra::CrsMatrix< scalar_t, local_ordinal_t, global_ordinal_t, node_t > matrix_t
IpplTimings::TimerRef setupTimer_m
MueLu::Hierarchy< scalar_t, lo_t, go_t, node_t > hierarchy_t
void setOperator(const Teuchos::RCP< matrix_t > &A, Level *level_p=nullptr)
amr::matrix_t matrix_t
MueLu::ParameterListInterpreter< scalar_t, lo_t, go_t, node_t > pListInterpreter_t
MueLu::Utilities< scalar_t, lo_t, go_t, node_t > util_t
std::size_t getNumIters()
amr::global_ordinal_t go_t
amr::AmrBox_t AmrBox_t
amr::multivector_t mv_t
amr::operator_t op_t
MueLuBottomSolver(const bool &rebalance, const std::string &reuse)
void solve(const Teuchos::RCP< mv_t > &x, const Teuchos::RCP< mv_t > &b)
Teuchos::RCP< xmatrix_t > A_mp
MueLu requires Xpetra.
Xpetra::MultiVector< scalar_t, lo_t, go_t, node_t > xmv_t
amr::AmrIntVect_t AmrIntVect_t
Teuchos::RCP< manager_t > factory_mp
sets up hierarchy
amr::vector_t vector_t
static std::string convertToMueLuReuseOption(const std::string &reuse)
MueLu::Level level_t
void initMueLuList_m(const std::string &reuse)
MueLu::HierarchyManager< scalar_t, lo_t, go_t, node_t > manager_t
Teuchos::ParameterList mueluList_m
bool rebalance_m
use subcommunicators (less communication)
amr::local_ordinal_t lo_t
Teuchos::RCP< hierarchy_t > hierarchy_mp
manages the multigrid hierarchy
Xpetra::Matrix< scalar_t, lo_t, go_t, node_t > xmatrix_t
amr::scalar_t scalar_t
lo_t nSweeps_m
the number of multigrid iterations
Timing::TimerRef TimerRef
Definition: IpplTimings.h:176