Reference documentation for deal.II version 8.5.1
trilinos_precondition_ml.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2008 - 2017 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_precondition.h>
17 
18 #ifdef DEAL_II_WITH_TRILINOS
19 
20 # include <deal.II/lac/vector.h>
21 # include <deal.II/lac/sparse_matrix.h>
22 # include <deal.II/lac/trilinos_sparse_matrix.h>
23 
24 DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
25 # include <Ifpack.h>
26 # include <Ifpack_Chebyshev.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>
32 DEAL_II_ENABLE_EXTRA_DIAGNOSTICS
33 
34 DEAL_II_NAMESPACE_OPEN
35 
36 namespace TrilinosWrappers
37 {
38  namespace
39  {
40 #ifndef DEAL_II_WITH_64BIT_INDICES
41  int n_global_rows (const Epetra_RowMatrix &matrix)
42  {
43  return matrix.NumGlobalRows();
44  }
45 
46  int global_length (const Epetra_MultiVector &vector)
47  {
48  return vector.GlobalLength();
49  }
50 
51  int gid(const Epetra_Map &map, unsigned int i)
52  {
53  return map.GID(i);
54  }
55 #else
56  long long int n_global_rows (const Epetra_RowMatrix &matrix)
57  {
58  return matrix.NumGlobalRows64();
59  }
60 
61  long long int global_length (const Epetra_MultiVector &vector)
62  {
63  return vector.GlobalLength64();
64  }
65 
66  long long int gid(const Epetra_Map &map, ::types::global_dof_index i)
67  {
68  return map.GID64(i);
69  }
70 #endif
71  }
72 
73 
74 
75  /* -------------------------- PreconditionAMG -------------------------- */
76 
78  AdditionalData (const bool elliptic,
79  const bool higher_order_elements,
80  const unsigned int n_cycles,
81  const bool w_cycle,
82  const double aggregation_threshold,
83  const std::vector<std::vector<bool> > &constant_modes,
84  const unsigned int smoother_sweeps,
85  const unsigned int smoother_overlap,
86  const bool output_details,
87  const char *smoother_type,
88  const char *coarse_type)
89  :
90  elliptic (elliptic),
91  higher_order_elements (higher_order_elements),
92  n_cycles (n_cycles),
93  w_cycle (w_cycle),
94  aggregation_threshold (aggregation_threshold),
95  constant_modes (constant_modes),
96  smoother_sweeps (smoother_sweeps),
97  smoother_overlap (smoother_overlap),
98  output_details (output_details),
99  smoother_type (smoother_type),
100  coarse_type (coarse_type)
101  {}
102 
103 
105  {
106  preconditioner.reset();
107  trilinos_matrix.reset();
108  }
109 
110 
111 
112  void
114  const AdditionalData &additional_data)
115  {
116  initialize(matrix.trilinos_matrix(), additional_data);
117  }
118 
119 
120 
121  void
122  PreconditionAMG::initialize (const Epetra_RowMatrix &matrix,
123  const AdditionalData &additional_data)
124  {
125  // Build the AMG preconditioner.
126  Teuchos::ParameterList parameter_list;
127 
128  if (additional_data.elliptic == true)
129  {
130  ML_Epetra::SetDefaults("SA",parameter_list);
131 
132  // uncoupled mode can give a lot of warnings or even fail when there
133  // are too many entries per row and aggreggation gets complicated, but
134  // MIS does not work if too few elements are located on one
135  // processor. work around these warnings by choosing the different
136  // strategies in different situations: for low order, always use the
137  // standard choice uncoupled. if higher order, right now we also just
138  // use Uncoupled, but we should be aware that maybe MIS might be
139  // needed
140  if (additional_data.higher_order_elements)
141  parameter_list.set("aggregation: type", "Uncoupled");
142  }
143  else
144  {
145  ML_Epetra::SetDefaults("NSSA",parameter_list);
146  parameter_list.set("aggregation: type", "Uncoupled");
147  parameter_list.set("aggregation: block scaling", true);
148  }
149 
150  parameter_list.set("smoother: type", additional_data.smoother_type);
151  parameter_list.set("coarse: type", additional_data.coarse_type);
152 
153  // Force re-initialization of the random seed to make ML deterministic
154  // (only supported in trilinos >12.2):
155 #if DEAL_II_TRILINOS_VERSION_GTE(12,4,0)
156  parameter_list.set("initialize random seed", true);
157 #endif
158 
159  parameter_list.set("smoother: sweeps",
160  static_cast<int>(additional_data.smoother_sweeps));
161  parameter_list.set("cycle applications",
162  static_cast<int>(additional_data.n_cycles));
163  if (additional_data.w_cycle == true)
164  parameter_list.set("prec type", "MGW");
165  else
166  parameter_list.set("prec type", "MGV");
167 
168  parameter_list.set("smoother: Chebyshev alpha",10.);
169  parameter_list.set("smoother: ifpack overlap",
170  static_cast<int>(additional_data.smoother_overlap));
171  parameter_list.set("aggregation: threshold",
172  additional_data.aggregation_threshold);
173  parameter_list.set("coarse: max size", 2000);
174 
175  if (additional_data.output_details)
176  parameter_list.set("ML output", 10);
177  else
178  parameter_list.set("ML output", 0);
179 
180  const Epetra_Map &domain_map = matrix.OperatorDomainMap();
181 
182  const size_type constant_modes_dimension =
183  additional_data.constant_modes.size();
184  Epetra_MultiVector distributed_constant_modes (domain_map,
185  constant_modes_dimension > 0 ?
186  constant_modes_dimension : 1);
187  std::vector<double> dummy (constant_modes_dimension);
188 
189  if (constant_modes_dimension > 0)
190  {
191  const size_type global_size = n_global_rows(matrix);
192  (void)global_length; // work around compiler warning about unused function in release mode
193  Assert (global_size ==
194  static_cast<size_type>(global_length(distributed_constant_modes)),
195  ExcDimensionMismatch(global_size,
196  global_length(distributed_constant_modes)));
197  const bool constant_modes_are_global
198  = additional_data.constant_modes[0].size() == global_size;
199  const size_type my_size = domain_map.NumMyElements();
200 
201  // Reshape null space as a contiguous vector of doubles so that
202  // Trilinos can read from it.
203  const size_type expected_mode_size =
204  constant_modes_are_global ? global_size : my_size;
205  for (size_type d=0; d<constant_modes_dimension; ++d)
206  {
207  Assert (additional_data.constant_modes[d].size() == expected_mode_size,
208  ExcDimensionMismatch(additional_data.constant_modes[d].size(), expected_mode_size));
209  for (size_type row=0; row<my_size; ++row)
210  {
211  const TrilinosWrappers::types::int_type mode_index =
212  constant_modes_are_global ? gid(domain_map,row) : row;
213  distributed_constant_modes[d][row] =
214  additional_data.constant_modes[d][mode_index];
215  }
216  }
217  (void)expected_mode_size;
218 
219  parameter_list.set("null space: type", "pre-computed");
220  parameter_list.set("null space: dimension",
221  distributed_constant_modes.NumVectors());
222  if (my_size > 0)
223  parameter_list.set("null space: vectors",
224  distributed_constant_modes.Values());
225  // We need to set a valid pointer to data even if there is no data on
226  // the current processor. Therefore, pass a dummy in that case
227  else
228  parameter_list.set("null space: vectors",
229  &dummy[0]);
230  }
231 
232  initialize (matrix, parameter_list);
233 
234  if (additional_data.output_details)
235  {
236  ML_Epetra::MultiLevelPreconditioner *multilevel_operator =
237  dynamic_cast<ML_Epetra::MultiLevelPreconditioner *> (preconditioner.get());
238  Assert (multilevel_operator != 0,
239  ExcMessage ("Preconditioner setup failed."));
240  multilevel_operator->PrintUnused(0);
241  }
242  }
243 
244 
245 
246  void
248  const Teuchos::ParameterList &ml_parameters)
249  {
250  initialize(matrix.trilinos_matrix(), ml_parameters);
251  }
252 
253 
254 
255  void
256  PreconditionAMG::initialize (const Epetra_RowMatrix &matrix,
257  const Teuchos::ParameterList &ml_parameters)
258  {
259  preconditioner.reset ();
260  preconditioner.reset (new ML_Epetra::MultiLevelPreconditioner
261  (matrix, ml_parameters));
262  }
263 
264 
265 
266  template <typename number>
267  void
269  initialize (const ::SparseMatrix<number> &deal_ii_sparse_matrix,
270  const AdditionalData &additional_data,
271  const double drop_tolerance,
272  const ::SparsityPattern *use_this_sparsity)
273  {
274  preconditioner.reset();
275  const size_type n_rows = deal_ii_sparse_matrix.m();
276 
277  // Init Epetra Matrix using an
278  // equidistributed map; avoid
279  // storing the nonzero
280  // elements.
281  vector_distributor.reset (new Epetra_Map(static_cast<TrilinosWrappers::types::int_type>(n_rows),
282  0, communicator));
283 
284  if (trilinos_matrix.get() == 0)
285  trilinos_matrix.reset (new SparseMatrix());
286 
288  deal_ii_sparse_matrix, drop_tolerance, true,
289  use_this_sparsity);
290 
291  initialize (*trilinos_matrix, additional_data);
292  }
293 
294 
295 
297  {
298  ML_Epetra::MultiLevelPreconditioner *multilevel_operator =
299  dynamic_cast<ML_Epetra::MultiLevelPreconditioner *> (preconditioner.get());
300  multilevel_operator->ReComputePreconditioner();
301  }
302 
303 
304 
306  {
308  trilinos_matrix.reset();
309  }
310 
311 
312 
315  {
316  unsigned int memory = sizeof(*this);
317 
318  // todo: find a way to read out ML's data
319  // sizes
320  if (trilinos_matrix.get() != 0)
321  memory += trilinos_matrix->memory_consumption();
322  return memory;
323  }
324 
325 
326 
327 
328  // explicit instantiations
329  template void PreconditionAMG::initialize (const ::SparseMatrix<double> &,
330  const AdditionalData &, const double,
331  const ::SparsityPattern *);
332  template void PreconditionAMG::initialize (const ::SparseMatrix<float> &,
333  const AdditionalData &, const double,
334  const ::SparsityPattern *);
335 
336 
337 
338 
339 
340 
341 }
342 
343 DEAL_II_NAMESPACE_CLOSE
344 
345 #endif // DEAL_II_WITH_TRILINOS
std_cxx11::shared_ptr< SparseMatrix > trilinos_matrix
std_cxx11::shared_ptr< Epetra_Map > vector_distributor
static ::ExceptionBase & ExcMessage(std::string arg1)
std_cxx11::shared_ptr< Epetra_Operator > preconditioner
unsigned int global_dof_index
Definition: types.h:88
#define Assert(cond, exc)
Definition: exceptions.h:313
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")
void initialize(const SparseMatrix &matrix, const AdditionalData &additional_data=AdditionalData())