16 #include <deal.II/lac/trilinos_index_access.h> 17 #include <deal.II/lac/trilinos_precondition.h> 19 #ifdef DEAL_II_WITH_TRILINOS 20 #if DEAL_II_TRILINOS_VERSION_GTE(11,14,0) 22 # include <deal.II/lac/vector.h> 23 # include <deal.II/lac/sparse_matrix.h> 24 # include <deal.II/lac/trilinos_index_access.h> 25 # include <deal.II/lac/trilinos_sparse_matrix.h> 27 # include <Teuchos_ParameterList.hpp> 28 # include <Teuchos_RCP.hpp> 29 # include <Epetra_MultiVector.h> 30 # include <ml_include.h> 31 # include <ml_MultiLevelPreconditioner.h> 34 # include <MueLu_EpetraOperator.hpp> 35 # include <MueLu_MLParameterListInterpreter.hpp> 37 DEAL_II_NAMESPACE_OPEN
43 const unsigned int n_cycles,
45 const double aggregation_threshold,
46 const std::vector<std::vector<bool> > &constant_modes,
47 const unsigned int smoother_sweeps,
48 const unsigned int smoother_overlap,
49 const bool output_details,
50 const char *smoother_type,
51 const char *coarse_type)
56 aggregation_threshold (aggregation_threshold),
57 constant_modes (constant_modes),
58 smoother_sweeps (smoother_sweeps),
59 smoother_overlap (smoother_overlap),
60 output_details (output_details),
61 smoother_type (smoother_type),
62 coarse_type (coarse_type)
69 #ifdef DEAL_II_WITH_64BIT_INDICES 88 initialize(matrix.trilinos_matrix(), additional_data);
98 Teuchos::ParameterList parameter_list;
100 if (additional_data.
elliptic ==
true)
101 ML_Epetra::SetDefaults(
"SA",parameter_list);
104 ML_Epetra::SetDefaults(
"NSSA",parameter_list);
105 parameter_list.set(
"aggregation: block scaling",
true);
110 parameter_list.set(
"aggregation: type",
"Uncoupled");
112 parameter_list.set(
"smoother: type", additional_data.
smoother_type);
113 parameter_list.set(
"coarse: type", additional_data.
coarse_type);
115 parameter_list.set(
"smoother: sweeps",
117 parameter_list.set(
"cycle applications",
118 static_cast<int>(additional_data.
n_cycles));
119 if (additional_data.
w_cycle ==
true)
120 parameter_list.set(
"prec type",
"MGW");
122 parameter_list.set(
"prec type",
"MGV");
124 parameter_list.set(
"smoother: Chebyshev alpha",10.);
125 parameter_list.set(
"smoother: ifpack overlap",
127 parameter_list.set(
"aggregation: threshold",
129 parameter_list.set(
"coarse: max size", 2000);
132 parameter_list.set(
"ML output", 10);
134 parameter_list.set(
"ML output", 0);
136 const Epetra_Map &domain_map = matrix.OperatorDomainMap();
138 const size_type constant_modes_dimension =
140 Epetra_MultiVector distributed_constant_modes (domain_map,
141 constant_modes_dimension > 0 ?
142 constant_modes_dimension : 1);
143 std::vector<double> dummy (constant_modes_dimension);
145 if (constant_modes_dimension > 0)
148 const bool constant_modes_are_global =
151 constant_modes_are_global ? n_rows : additional_data.
constant_modes[0].size();
152 const size_type my_size = domain_map.NumMyElements();
153 if (constant_modes_are_global ==
false)
154 Assert (n_relevant_rows == my_size,
161 (void)n_relevant_rows;
166 for (
size_type d=0; d<constant_modes_dimension; ++d)
167 for (
size_type row=0; row<my_size; ++row)
169 TrilinosWrappers::types::int_type global_row_id =
171 distributed_constant_modes[d][row] =
175 parameter_list.set(
"null space: type",
"pre-computed");
176 parameter_list.set(
"null space: dimension",
177 distributed_constant_modes.NumVectors());
179 parameter_list.set(
"null space: vectors",
180 distributed_constant_modes.Values());
184 parameter_list.set(
"null space: vectors",
195 Teuchos::ParameterList &muelu_parameters)
197 initialize(matrix.trilinos_matrix(), muelu_parameters);
204 Teuchos::ParameterList &muelu_parameters)
211 typedef KokkosClassic::DefaultNode::DefaultNodeType node;
216 Teuchos::RCP<Epetra_CrsMatrix> rcp_matrix = Teuchos::rcpFromRef(
217 *(const_cast<Epetra_CrsMatrix *>(&matrix)));
218 Teuchos::RCP<Xpetra::CrsMatrix<double,int,int,node> > muelu_crs_matrix =
219 Teuchos::rcp(
new Xpetra::EpetraCrsMatrix (rcp_matrix));
220 Teuchos::RCP<Xpetra::Matrix<double,int,int,node> > muelu_matrix =
221 Teuchos::rcp(
new Xpetra::CrsMatrixWrap<double,int,int,node> (muelu_crs_matrix));
224 Teuchos::RCP<MueLu::HierarchyManager<double,int,int,node> > hierarchy_factory;
225 hierarchy_factory = Teuchos::rcp(
226 new MueLu::MLParameterListInterpreter<double,int,int,node> (muelu_parameters));
227 Teuchos::RCP<MueLu::Hierarchy<double,int,int,node> > hierarchy = hierarchy_factory->CreateHierarchy();
228 hierarchy->GetLevel(0)->Set(
"A",muelu_matrix);
229 hierarchy_factory->SetupHierarchy(*hierarchy);
233 preconditioner = std::make_shared<MueLu::EpetraOperator>(hierarchy);
238 template <
typename number>
241 initialize (const ::SparseMatrix<number> &deal_ii_sparse_matrix,
243 const double drop_tolerance,
244 const ::SparsityPattern *use_this_sparsity)
247 const size_type n_rows = deal_ii_sparse_matrix.m();
254 (
static_cast<TrilinosWrappers::types::int_type
>(n_rows),
261 deal_ii_sparse_matrix, drop_tolerance,
true,
280 unsigned int memory =
sizeof(*this);
293 const AdditionalData &,
const double,
294 const ::SparsityPattern *);
296 const AdditionalData &,
const double,
297 const ::SparsityPattern *);
301 DEAL_II_NAMESPACE_CLOSE
303 #endif // DEAL_II_TRILINOS_VERSION_GTE(11,14,0) 304 #endif // DEAL_II_WITH_TRILINOS unsigned int smoother_sweeps
void initialize(const SparseMatrix &matrix, const AdditionalData &additional_data=AdditionalData())
AdditionalData(const bool elliptic=true, const unsigned int n_cycles=1, const bool w_cyle=false, const double aggregation_threshold=1e-4, const std::vector< std::vector< bool > > &constant_modes=std::vector< std::vector< bool > >(0), const unsigned int smoother_sweeps=2, const unsigned int smoother_overlap=0, const bool output_details=false, const char *smoother_type="Chebyshev", const char *coarse_type="Amesos-KLU")
TrilinosWrappers::types::int_type global_length(const Epetra_MultiVector &vector)
const char * smoother_type
std::shared_ptr< Epetra_Map > vector_distributor
#define AssertThrow(cond, exc)
size_type memory_consumption() const
static ::ExceptionBase & ExcMessage(std::string arg1)
TrilinosWrappers::types::int_type n_global_rows(const Epetra_CrsGraph &graph)
double aggregation_threshold
#define Assert(cond, exc)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
std::vector< std::vector< bool > > constant_modes
std::shared_ptr< SparseMatrix > trilinos_matrix
TrilinosWrappers::types::int_type global_index(const Epetra_BlockMap &map, const ::types::global_dof_index i)
std::shared_ptr< Epetra_Operator > preconditioner
unsigned int smoother_overlap
::types::global_dof_index size_type
Epetra_MpiComm communicator