16 #include <deal.II/lac/trilinos_precondition.h> 18 #ifdef DEAL_II_WITH_TRILINOS 19 #if DEAL_II_TRILINOS_VERSION_GTE(11,14,0) 21 # include <deal.II/lac/vector.h> 22 # include <deal.II/lac/sparse_matrix.h> 23 # include <deal.II/lac/trilinos_sparse_matrix.h> 25 DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
26 # include <Teuchos_ParameterList.hpp> 27 # include <Teuchos_RCP.hpp> 28 # include <Epetra_MultiVector.h> 29 # include <ml_include.h> 30 # include <ml_MultiLevelPreconditioner.h> 33 # include <MueLu_EpetraOperator.hpp> 34 # include <MueLu_MLParameterListInterpreter.hpp> 35 DEAL_II_ENABLE_EXTRA_DIAGNOSTICS
37 DEAL_II_NAMESPACE_OPEN
43 #ifndef DEAL_II_WITH_64BIT_INDICES 44 int n_global_rows (
const Epetra_RowMatrix &matrix)
46 return matrix.NumGlobalRows();
49 int global_length (
const Epetra_MultiVector &vector)
51 return vector.GlobalLength();
54 int gid(
const Epetra_Map &map,
unsigned int i)
59 long long int n_global_rows (
const Epetra_RowMatrix &matrix)
61 return matrix.NumGlobalRows64();
64 long long int global_length (
const Epetra_MultiVector &vector)
66 return vector.GlobalLength64();
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)
93 aggregation_threshold (aggregation_threshold),
94 constant_modes (constant_modes),
95 smoother_sweeps (smoother_sweeps),
96 smoother_overlap (smoother_overlap),
97 output_details (output_details),
98 smoother_type (smoother_type),
99 coarse_type (coarse_type)
106 #ifdef DEAL_II_WITH_64BIT_INDICES 125 initialize(matrix.trilinos_matrix(), additional_data);
135 Teuchos::ParameterList parameter_list;
137 if (additional_data.
elliptic ==
true)
138 ML_Epetra::SetDefaults(
"SA",parameter_list);
141 ML_Epetra::SetDefaults(
"NSSA",parameter_list);
142 parameter_list.set(
"aggregation: block scaling",
true);
147 parameter_list.set(
"aggregation: type",
"Uncoupled");
149 parameter_list.set(
"smoother: type", additional_data.
smoother_type);
150 parameter_list.set(
"coarse: type", additional_data.
coarse_type);
152 parameter_list.set(
"smoother: sweeps",
154 parameter_list.set(
"cycle applications",
155 static_cast<int>(additional_data.
n_cycles));
156 if (additional_data.
w_cycle ==
true)
157 parameter_list.set(
"prec type",
"MGW");
159 parameter_list.set(
"prec type",
"MGV");
161 parameter_list.set(
"smoother: Chebyshev alpha",10.);
162 parameter_list.set(
"smoother: ifpack overlap",
164 parameter_list.set(
"aggregation: threshold",
166 parameter_list.set(
"coarse: max size", 2000);
169 parameter_list.set(
"ML output", 10);
171 parameter_list.set(
"ML output", 0);
173 const Epetra_Map &domain_map = matrix.OperatorDomainMap();
175 const size_type constant_modes_dimension =
177 Epetra_MultiVector distributed_constant_modes (domain_map,
178 constant_modes_dimension > 0 ?
179 constant_modes_dimension : 1);
180 std::vector<double> dummy (constant_modes_dimension);
182 if (constant_modes_dimension > 0)
184 const size_type n_rows = n_global_rows(matrix);
185 const bool constant_modes_are_global =
188 constant_modes_are_global ? n_rows : additional_data.
constant_modes[0].size();
189 const size_type my_size = domain_map.NumMyElements();
190 if (constant_modes_are_global ==
false)
191 Assert (n_relevant_rows == my_size,
194 static_cast<size_type>(global_length(distributed_constant_modes)),
196 global_length(distributed_constant_modes)));
198 (void)n_relevant_rows;
203 for (
size_type d=0; d<constant_modes_dimension; ++d)
204 for (
size_type row=0; row<my_size; ++row)
206 TrilinosWrappers::types::int_type global_row_id =
207 constant_modes_are_global ? gid(domain_map,row) : row;
208 distributed_constant_modes[d][row] =
212 parameter_list.set(
"null space: type",
"pre-computed");
213 parameter_list.set(
"null space: dimension",
214 distributed_constant_modes.NumVectors());
216 parameter_list.set(
"null space: vectors",
217 distributed_constant_modes.Values());
221 parameter_list.set(
"null space: vectors",
232 Teuchos::ParameterList &muelu_parameters)
234 initialize(matrix.trilinos_matrix(), muelu_parameters);
241 Teuchos::ParameterList &muelu_parameters)
248 typedef KokkosClassic::DefaultNode::DefaultNodeType node;
253 Teuchos::RCP<Epetra_CrsMatrix> rcp_matrix = Teuchos::rcpFromRef(
254 *(const_cast<Epetra_CrsMatrix *>(&matrix)));
255 Teuchos::RCP<Xpetra::CrsMatrix<double,int,int,node> > muelu_crs_matrix =
256 Teuchos::rcp(
new Xpetra::EpetraCrsMatrix (rcp_matrix));
257 Teuchos::RCP<Xpetra::Matrix<double,int,int,node> > muelu_matrix =
258 Teuchos::rcp(
new Xpetra::CrsMatrixWrap<double,int,int,node> (muelu_crs_matrix));
261 Teuchos::RCP<MueLu::HierarchyManager<double,int,int,node> > hierarchy_factory;
262 hierarchy_factory = Teuchos::rcp(
263 new MueLu::MLParameterListInterpreter<double,int,int,node> (muelu_parameters));
264 Teuchos::RCP<MueLu::Hierarchy<double,int,int,node> > hierarchy = hierarchy_factory->CreateHierarchy();
265 hierarchy->GetLevel(0)->Set(
"A",muelu_matrix);
266 hierarchy_factory->SetupHierarchy(*hierarchy);
275 template <
typename number>
278 initialize (const ::SparseMatrix<number> &deal_ii_sparse_matrix,
280 const double drop_tolerance,
281 const ::SparsityPattern *use_this_sparsity)
284 const size_type n_rows = deal_ii_sparse_matrix.m();
290 vector_distributor.reset (
new Epetra_Map(static_cast<TrilinosWrappers::types::int_type>(n_rows),
297 deal_ii_sparse_matrix, drop_tolerance,
true,
316 unsigned int memory =
sizeof(*this);
329 const AdditionalData &,
const double,
330 const ::SparsityPattern *);
332 const AdditionalData &,
const double,
333 const ::SparsityPattern *);
337 DEAL_II_NAMESPACE_CLOSE
339 #endif // DEAL_II_TRILINOS_VERSION_GTE(11,14,0) 340 #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")
const char * smoother_type
#define AssertThrow(cond, exc)
size_type memory_consumption() const
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
double aggregation_threshold
#define Assert(cond, exc)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
std_cxx11::shared_ptr< SparseMatrix > trilinos_matrix
std::vector< std::vector< bool > > constant_modes
unsigned int smoother_overlap
::types::global_dof_index size_type
Epetra_MpiComm communicator