Reference documentation for deal.II version 9.0.0
trilinos_precondition_muelu.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2008 - 2018 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 #include <deal.II/lac/trilinos_index_access.h>
17 #include <deal.II/lac/trilinos_precondition.h>
18 
19 #ifdef DEAL_II_WITH_TRILINOS
20 #if DEAL_II_TRILINOS_VERSION_GTE(11,14,0)
21 
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>
26 
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>
32 
33 # include <MueLu.hpp>
34 # include <MueLu_EpetraOperator.hpp>
35 # include <MueLu_MLParameterListInterpreter.hpp>
36 
37 DEAL_II_NAMESPACE_OPEN
38 
39 namespace TrilinosWrappers
40 {
42  AdditionalData (const bool elliptic,
43  const unsigned int n_cycles,
44  const bool w_cycle,
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)
52  :
53  elliptic (elliptic),
54  n_cycles (n_cycles),
55  w_cycle (w_cycle),
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)
63  {}
64 
65 
66 
68  {
69 #ifdef DEAL_II_WITH_64BIT_INDICES
70  AssertThrow (false, ExcMessage("PreconditionAMGMueLu does not support 64bit-indices!"));
71 #endif
72  }
73 
74 
75 
77  {
78  preconditioner.reset();
79  trilinos_matrix.reset();
80  }
81 
82 
83 
84  void
86  const AdditionalData &additional_data)
87  {
88  initialize(matrix.trilinos_matrix(), additional_data);
89  }
90 
91 
92 
93  void
94  PreconditionAMGMueLu::initialize (const Epetra_CrsMatrix &matrix,
95  const AdditionalData &additional_data)
96  {
97  // Build the AMG preconditioner.
98  Teuchos::ParameterList parameter_list;
99 
100  if (additional_data.elliptic == true)
101  ML_Epetra::SetDefaults("SA",parameter_list);
102  else
103  {
104  ML_Epetra::SetDefaults("NSSA",parameter_list);
105  parameter_list.set("aggregation: block scaling", true);
106  }
107  // MIS does not exist anymore, only choice are uncoupled and coupled. When using
108  // uncoupled, aggregates cannot span multiple processes. When using coupled
109  // aggregates can span multiple processes.
110  parameter_list.set("aggregation: type", "Uncoupled");
111 
112  parameter_list.set("smoother: type", additional_data.smoother_type);
113  parameter_list.set("coarse: type", additional_data.coarse_type);
114 
115  parameter_list.set("smoother: sweeps",
116  static_cast<int>(additional_data.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");
121  else
122  parameter_list.set("prec type", "MGV");
123 
124  parameter_list.set("smoother: Chebyshev alpha",10.);
125  parameter_list.set("smoother: ifpack overlap",
126  static_cast<int>(additional_data.smoother_overlap));
127  parameter_list.set("aggregation: threshold",
128  additional_data.aggregation_threshold);
129  parameter_list.set("coarse: max size", 2000);
130 
131  if (additional_data.output_details)
132  parameter_list.set("ML output", 10);
133  else
134  parameter_list.set("ML output", 0);
135 
136  const Epetra_Map &domain_map = matrix.OperatorDomainMap();
137 
138  const size_type constant_modes_dimension =
139  additional_data.constant_modes.size();
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);
144 
145  if (constant_modes_dimension > 0)
146  {
147  const size_type n_rows = TrilinosWrappers::n_global_rows(matrix);
148  const bool constant_modes_are_global =
149  additional_data.constant_modes[0].size() == n_rows;
150  const size_type n_relevant_rows =
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,
155  ExcDimensionMismatch(n_relevant_rows, my_size));
156  Assert (n_rows ==
157  static_cast<size_type>(TrilinosWrappers::global_length(distributed_constant_modes)),
158  ExcDimensionMismatch(n_rows,
159  TrilinosWrappers::global_length(distributed_constant_modes)));
160 
161  (void)n_relevant_rows;
162  (void)global_length;
163 
164  // Reshape null space as a contiguous vector of doubles so that
165  // Trilinos can read from it.
166  for (size_type d=0; d<constant_modes_dimension; ++d)
167  for (size_type row=0; row<my_size; ++row)
168  {
169  TrilinosWrappers::types::int_type global_row_id =
170  constant_modes_are_global ? TrilinosWrappers::global_index(domain_map,row) : row;
171  distributed_constant_modes[d][row] =
172  additional_data.constant_modes[d][global_row_id];
173  }
174 
175  parameter_list.set("null space: type", "pre-computed");
176  parameter_list.set("null space: dimension",
177  distributed_constant_modes.NumVectors());
178  if (my_size > 0)
179  parameter_list.set("null space: vectors",
180  distributed_constant_modes.Values());
181  // We need to set a valid pointer to data even if there is no data on
182  // the current processor. Therefore, pass a dummy in that case
183  else
184  parameter_list.set("null space: vectors",
185  dummy.data());
186  }
187 
188  initialize (matrix, parameter_list);
189  }
190 
191 
192 
193  void
195  Teuchos::ParameterList &muelu_parameters)
196  {
197  initialize(matrix.trilinos_matrix(), muelu_parameters);
198  }
199 
200 
201 
202  void
203  PreconditionAMGMueLu::initialize (const Epetra_CrsMatrix &matrix,
204  Teuchos::ParameterList &muelu_parameters)
205  {
206  // We cannot use MueLu::CreateEpetraOperator directly because, we cannot
207  // transfer ownership of MueLu::EpetraOperator from Teuchos::RCP to
208  // std::shared_ptr.
209 
210  // For now, just use serial node, i.e. no multithreaing or GPU.
211  typedef KokkosClassic::DefaultNode::DefaultNodeType node;
212  preconditioner.reset ();
213 
214  // Cast matrix into a MueLu::Matrix. The constness needs to be cast away.
215  // MueLu uses Teuchos::RCP which are Trilinos version of std::shared_ptr.
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));
222 
223  // Create the multigrid hierarchy using ML parameters.
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);
230 
231  // MueLu::EpetraOperator is just a wrapper around a "standard"
232  // Epetra_Operator.
233  preconditioner = std::make_shared<MueLu::EpetraOperator>(hierarchy);
234  }
235 
236 
237 
238  template <typename number>
239  void
241  initialize (const ::SparseMatrix<number> &deal_ii_sparse_matrix,
242  const AdditionalData &additional_data,
243  const double drop_tolerance,
244  const ::SparsityPattern *use_this_sparsity)
245  {
246  preconditioner.reset();
247  const size_type n_rows = deal_ii_sparse_matrix.m();
248 
249  // Init Epetra Matrix using an
250  // equidistributed map; avoid
251  // storing the nonzero
252  // elements.
253  vector_distributor = std::make_shared<Epetra_Map>
254  (static_cast<TrilinosWrappers::types::int_type>(n_rows),
255  0, communicator);
256 
257  if (trilinos_matrix.get() == nullptr)
258  trilinos_matrix = std::make_shared<SparseMatrix>();
259 
261  deal_ii_sparse_matrix, drop_tolerance, true,
262  use_this_sparsity);
263 
264  initialize (*trilinos_matrix, additional_data);
265  }
266 
267 
268 
270  {
272  trilinos_matrix.reset();
273  }
274 
275 
276 
279  {
280  unsigned int memory = sizeof(*this);
281 
282  // todo: find a way to read out ML's data
283  // sizes
284  if (trilinos_matrix.get() != nullptr)
285  memory += trilinos_matrix->memory_consumption();
286  return memory;
287  }
288 
289 
290 
291  // explicit instantiations
292  template void PreconditionAMGMueLu::initialize (const ::SparseMatrix<double> &,
293  const AdditionalData &, const double,
294  const ::SparsityPattern *);
295  template void PreconditionAMGMueLu::initialize (const ::SparseMatrix<float> &,
296  const AdditionalData &, const double,
297  const ::SparsityPattern *);
298 
299 }
300 
301 DEAL_II_NAMESPACE_CLOSE
302 
303 #endif // DEAL_II_TRILINOS_VERSION_GTE(11,14,0)
304 #endif // DEAL_II_WITH_TRILINOS
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)
std::shared_ptr< Epetra_Map > vector_distributor
#define AssertThrow(cond, exc)
Definition: exceptions.h:1221
static ::ExceptionBase & ExcMessage(std::string arg1)
TrilinosWrappers::types::int_type n_global_rows(const Epetra_CrsGraph &graph)
#define Assert(cond, exc)
Definition: exceptions.h:1142
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
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