Reference documentation for deal.II version 9.0.0
trilinos_precondition_ml.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 
21 # include <deal.II/lac/vector.h>
22 # include <deal.II/lac/sparse_matrix.h>
23 # include <deal.II/lac/trilinos_sparse_matrix.h>
24 # include <deal.II/lac/trilinos_index_access.h>
25 
26 # include <Ifpack.h>
27 # include <Ifpack_Chebyshev.h>
28 # include <Teuchos_ParameterList.hpp>
29 # include <Teuchos_RCP.hpp>
30 # include <Epetra_MultiVector.h>
31 # include <ml_include.h>
32 # include <ml_MultiLevelPreconditioner.h>
33 
34 DEAL_II_NAMESPACE_OPEN
35 
36 namespace TrilinosWrappers
37 {
38  /* -------------------------- PreconditionAMG -------------------------- */
39 
41  AdditionalData (const bool elliptic,
42  const bool higher_order_elements,
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  higher_order_elements (higher_order_elements),
55  n_cycles (n_cycles),
56  w_cycle (w_cycle),
57  aggregation_threshold (aggregation_threshold),
58  constant_modes (constant_modes),
59  smoother_sweeps (smoother_sweeps),
60  smoother_overlap (smoother_overlap),
61  output_details (output_details),
62  smoother_type (smoother_type),
63  coarse_type (coarse_type)
64  {}
65 
66 
68  {
69  preconditioner.reset();
70  trilinos_matrix.reset();
71  }
72 
73 
74 
75  void
77  const AdditionalData &additional_data)
78  {
79  initialize(matrix.trilinos_matrix(), additional_data);
80  }
81 
82 
83 
84  void
85  PreconditionAMG::initialize (const Epetra_RowMatrix &matrix,
86  const AdditionalData &additional_data)
87  {
88  // Build the AMG preconditioner.
89  Teuchos::ParameterList parameter_list;
90 
91  if (additional_data.elliptic == true)
92  {
93  ML_Epetra::SetDefaults("SA",parameter_list);
94 
95  // uncoupled mode can give a lot of warnings or even fail when there
96  // are too many entries per row and aggreggation gets complicated, but
97  // MIS does not work if too few elements are located on one
98  // processor. work around these warnings by choosing the different
99  // strategies in different situations: for low order, always use the
100  // standard choice uncoupled. if higher order, right now we also just
101  // use Uncoupled, but we should be aware that maybe MIS might be
102  // needed
103  if (additional_data.higher_order_elements)
104  parameter_list.set("aggregation: type", "Uncoupled");
105  }
106  else
107  {
108  ML_Epetra::SetDefaults("NSSA",parameter_list);
109  parameter_list.set("aggregation: type", "Uncoupled");
110  parameter_list.set("aggregation: block scaling", true);
111  }
112 
113  parameter_list.set("smoother: type", additional_data.smoother_type);
114  parameter_list.set("coarse: type", additional_data.coarse_type);
115 
116  // Force re-initialization of the random seed to make ML deterministic
117  // (only supported in trilinos >12.2):
118 #if DEAL_II_TRILINOS_VERSION_GTE(12,4,0)
119  parameter_list.set("initialize random seed", true);
120 #endif
121 
122  parameter_list.set("smoother: sweeps",
123  static_cast<int>(additional_data.smoother_sweeps));
124  parameter_list.set("cycle applications",
125  static_cast<int>(additional_data.n_cycles));
126  if (additional_data.w_cycle == true)
127  parameter_list.set("prec type", "MGW");
128  else
129  parameter_list.set("prec type", "MGV");
130 
131  parameter_list.set("smoother: Chebyshev alpha",10.);
132  parameter_list.set("smoother: ifpack overlap",
133  static_cast<int>(additional_data.smoother_overlap));
134  parameter_list.set("aggregation: threshold",
135  additional_data.aggregation_threshold);
136  parameter_list.set("coarse: max size", 2000);
137 
138  if (additional_data.output_details)
139  parameter_list.set("ML output", 10);
140  else
141  parameter_list.set("ML output", 0);
142 
143  const Epetra_Map &domain_map = matrix.OperatorDomainMap();
144 
145  const size_type constant_modes_dimension =
146  additional_data.constant_modes.size();
147  Epetra_MultiVector distributed_constant_modes (domain_map,
148  constant_modes_dimension > 0 ?
149  constant_modes_dimension : 1);
150  std::vector<double> dummy (constant_modes_dimension);
151 
152  if (constant_modes_dimension > 0)
153  {
154  const size_type global_size = TrilinosWrappers::n_global_rows(matrix);
155  (void)global_length; // work around compiler warning about unused function in release mode
156  Assert (global_size ==
157  static_cast<size_type>(TrilinosWrappers::global_length(distributed_constant_modes)),
158  ExcDimensionMismatch(global_size,
159  TrilinosWrappers::global_length(distributed_constant_modes)));
160  const bool constant_modes_are_global
161  = additional_data.constant_modes[0].size() == global_size;
162  const size_type my_size = domain_map.NumMyElements();
163 
164  // Reshape null space as a contiguous vector of doubles so that
165  // Trilinos can read from it.
166  const size_type expected_mode_size =
167  constant_modes_are_global ? global_size : my_size;
168  for (size_type d=0; d<constant_modes_dimension; ++d)
169  {
170  Assert (additional_data.constant_modes[d].size() == expected_mode_size,
171  ExcDimensionMismatch(additional_data.constant_modes[d].size(), expected_mode_size));
172  for (size_type row=0; row<my_size; ++row)
173  {
174  const TrilinosWrappers::types::int_type mode_index =
175  constant_modes_are_global ? TrilinosWrappers::global_index(domain_map,row) : row;
176  distributed_constant_modes[d][row] =
177  additional_data.constant_modes[d][mode_index];
178  }
179  }
180  (void)expected_mode_size;
181 
182  parameter_list.set("null space: type", "pre-computed");
183  parameter_list.set("null space: dimension",
184  distributed_constant_modes.NumVectors());
185  if (my_size > 0)
186  parameter_list.set("null space: vectors",
187  distributed_constant_modes.Values());
188  // We need to set a valid pointer to data even if there is no data on
189  // the current processor. Therefore, pass a dummy in that case
190  else
191  parameter_list.set("null space: vectors",
192  dummy.data());
193  }
194 
195  initialize (matrix, parameter_list);
196 
197  if (additional_data.output_details)
198  {
199  ML_Epetra::MultiLevelPreconditioner *multilevel_operator =
200  dynamic_cast<ML_Epetra::MultiLevelPreconditioner *> (preconditioner.get());
201  Assert (multilevel_operator != nullptr,
202  ExcMessage ("Preconditioner setup failed."));
203  multilevel_operator->PrintUnused(0);
204  }
205  }
206 
207 
208 
209  void
211  const Teuchos::ParameterList &ml_parameters)
212  {
213  initialize(matrix.trilinos_matrix(), ml_parameters);
214  }
215 
216 
217 
218  void
219  PreconditionAMG::initialize (const Epetra_RowMatrix &matrix,
220  const Teuchos::ParameterList &ml_parameters)
221  {
222  preconditioner.reset ();
223  preconditioner = std::make_shared<ML_Epetra::MultiLevelPreconditioner>
224  (matrix, ml_parameters);
225  }
226 
227 
228 
229  template <typename number>
230  void
232  initialize (const ::SparseMatrix<number> &deal_ii_sparse_matrix,
233  const AdditionalData &additional_data,
234  const double drop_tolerance,
235  const ::SparsityPattern *use_this_sparsity)
236  {
237  preconditioner.reset();
238  const size_type n_rows = deal_ii_sparse_matrix.m();
239 
240  // Init Epetra Matrix using an
241  // equidistributed map; avoid
242  // storing the nonzero
243  // elements.
244  vector_distributor = std::make_shared<Epetra_Map>
245  (static_cast<TrilinosWrappers::types::int_type>(n_rows),
246  0, communicator);
247 
248  if (trilinos_matrix.get() == nullptr)
249  trilinos_matrix = std::make_shared<SparseMatrix>();
250 
252  deal_ii_sparse_matrix, drop_tolerance, true,
253  use_this_sparsity);
254 
255  initialize (*trilinos_matrix, additional_data);
256  }
257 
258 
259 
261  {
262  ML_Epetra::MultiLevelPreconditioner *multilevel_operator =
263  dynamic_cast<ML_Epetra::MultiLevelPreconditioner *> (preconditioner.get());
264  multilevel_operator->ReComputePreconditioner();
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 
292  // explicit instantiations
293  template void PreconditionAMG::initialize (const ::SparseMatrix<double> &,
294  const AdditionalData &, const double,
295  const ::SparsityPattern *);
296  template void PreconditionAMG::initialize (const ::SparseMatrix<float> &,
297  const AdditionalData &, const double,
298  const ::SparsityPattern *);
299 
300 
301 
302 
303 
304 
305 }
306 
307 DEAL_II_NAMESPACE_CLOSE
308 
309 #endif // DEAL_II_WITH_TRILINOS
TrilinosWrappers::types::int_type global_length(const Epetra_MultiVector &vector)
std::shared_ptr< Epetra_Map > vector_distributor
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)
AdditionalData(const bool elliptic=true, const bool higher_order_elements=false, 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")
std::shared_ptr< SparseMatrix > trilinos_matrix
void initialize(const SparseMatrix &matrix, const AdditionalData &additional_data=AdditionalData())
TrilinosWrappers::types::int_type global_index(const Epetra_BlockMap &map, const ::types::global_dof_index i)
std::shared_ptr< Epetra_Operator > preconditioner