16 #include <deal.II/lac/trilinos_precondition.h> 18 #ifdef DEAL_II_WITH_TRILINOS 20 # include <deal.II/lac/vector.h> 21 # include <deal.II/lac/sparse_matrix.h> 22 # include <deal.II/lac/trilinos_sparse_matrix.h> 24 DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
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
34 DEAL_II_NAMESPACE_OPEN
40 #ifndef DEAL_II_WITH_64BIT_INDICES 41 int n_global_rows (
const Epetra_RowMatrix &matrix)
43 return matrix.NumGlobalRows();
46 int global_length (
const Epetra_MultiVector &vector)
48 return vector.GlobalLength();
51 int gid(
const Epetra_Map &map,
unsigned int i)
56 long long int n_global_rows (
const Epetra_RowMatrix &matrix)
58 return matrix.NumGlobalRows64();
61 long long int global_length (
const Epetra_MultiVector &vector)
63 return vector.GlobalLength64();
79 const bool higher_order_elements,
80 const unsigned int n_cycles,
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)
91 higher_order_elements (higher_order_elements),
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)
116 initialize(matrix.trilinos_matrix(), additional_data);
126 Teuchos::ParameterList parameter_list;
128 if (additional_data.
elliptic ==
true)
130 ML_Epetra::SetDefaults(
"SA",parameter_list);
141 parameter_list.set(
"aggregation: type",
"Uncoupled");
145 ML_Epetra::SetDefaults(
"NSSA",parameter_list);
146 parameter_list.set(
"aggregation: type",
"Uncoupled");
147 parameter_list.set(
"aggregation: block scaling",
true);
150 parameter_list.set(
"smoother: type", additional_data.
smoother_type);
151 parameter_list.set(
"coarse: type", additional_data.
coarse_type);
155 #if DEAL_II_TRILINOS_VERSION_GTE(12,4,0) 156 parameter_list.set(
"initialize random seed",
true);
159 parameter_list.set(
"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");
166 parameter_list.set(
"prec type",
"MGV");
168 parameter_list.set(
"smoother: Chebyshev alpha",10.);
169 parameter_list.set(
"smoother: ifpack overlap",
171 parameter_list.set(
"aggregation: threshold",
173 parameter_list.set(
"coarse: max size", 2000);
176 parameter_list.set(
"ML output", 10);
178 parameter_list.set(
"ML output", 0);
180 const Epetra_Map &domain_map = matrix.OperatorDomainMap();
182 const size_type constant_modes_dimension =
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);
189 if (constant_modes_dimension > 0)
191 const size_type global_size = n_global_rows(matrix);
194 static_cast<size_type>(global_length(distributed_constant_modes)),
196 global_length(distributed_constant_modes)));
197 const bool constant_modes_are_global
199 const size_type my_size = domain_map.NumMyElements();
204 constant_modes_are_global ? global_size : my_size;
205 for (
size_type d=0; d<constant_modes_dimension; ++d)
209 for (
size_type row=0; row<my_size; ++row)
211 const TrilinosWrappers::types::int_type mode_index =
212 constant_modes_are_global ? gid(domain_map,row) : row;
213 distributed_constant_modes[d][row] =
217 (void)expected_mode_size;
219 parameter_list.set(
"null space: type",
"pre-computed");
220 parameter_list.set(
"null space: dimension",
221 distributed_constant_modes.NumVectors());
223 parameter_list.set(
"null space: vectors",
224 distributed_constant_modes.Values());
228 parameter_list.set(
"null space: vectors",
236 ML_Epetra::MultiLevelPreconditioner *multilevel_operator =
237 dynamic_cast<ML_Epetra::MultiLevelPreconditioner *
> (
preconditioner.get());
238 Assert (multilevel_operator != 0,
240 multilevel_operator->PrintUnused(0);
248 const Teuchos::ParameterList &ml_parameters)
250 initialize(matrix.trilinos_matrix(), ml_parameters);
257 const Teuchos::ParameterList &ml_parameters)
261 (matrix, ml_parameters));
266 template <
typename number>
269 initialize (const ::SparseMatrix<number> &deal_ii_sparse_matrix,
271 const double drop_tolerance,
272 const ::SparsityPattern *use_this_sparsity)
275 const size_type n_rows = deal_ii_sparse_matrix.m();
281 vector_distributor.reset (
new Epetra_Map(static_cast<TrilinosWrappers::types::int_type>(n_rows),
288 deal_ii_sparse_matrix, drop_tolerance,
true,
298 ML_Epetra::MultiLevelPreconditioner *multilevel_operator =
299 dynamic_cast<ML_Epetra::MultiLevelPreconditioner *
> (
preconditioner.get());
300 multilevel_operator->ReComputePreconditioner();
316 unsigned int memory =
sizeof(*this);
330 const AdditionalData &,
const double,
331 const ::SparsityPattern *);
333 const AdditionalData &,
const double,
334 const ::SparsityPattern *);
343 DEAL_II_NAMESPACE_CLOSE
345 #endif // DEAL_II_WITH_TRILINOS
double aggregation_threshold
std::vector< std::vector< bool > > constant_modes
unsigned int smoother_overlap
unsigned int smoother_sweeps
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
const char * smoother_type
#define Assert(cond, exc)
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")
size_type memory_consumption() const
void initialize(const SparseMatrix &matrix, const AdditionalData &additional_data=AdditionalData())
::types::global_dof_index size_type
Epetra_MpiComm communicator
bool higher_order_elements