Ifpack2 Templated Preconditioning Package Version 1.0
Loading...
Searching...
No Matches
Ifpack2_MDF_def.hpp
1// @HEADER
2// *****************************************************************************
3// Ifpack2: Templated Object-Oriented Algebraic Preconditioner Package
4//
5// Copyright 2009 NTESS and the Ifpack2 contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef IFPACK2_MDF_DEF_HPP
11#define IFPACK2_MDF_DEF_HPP
12
13#include "Ifpack2_LocalFilter.hpp"
15#include "Tpetra_CrsMatrix.hpp"
16#include "Teuchos_StandardParameterEntryValidators.hpp"
17#include "Ifpack2_LocalSparseTriangularSolver.hpp"
18#include "Ifpack2_Details_getParamTryingTypes.hpp"
19#include "Kokkos_Core.hpp"
20#include "Kokkos_Sort.hpp"
21#include "KokkosKernels_Sorting.hpp"
22#include <exception>
23#include <type_traits>
24
25namespace Ifpack2 {
26
27namespace Details {
28
29namespace MDFImpl
30{
31
32template<class dev_view_t>
33auto copy_view(const dev_view_t & vals)
34{
35 using Kokkos::view_alloc;
36 using Kokkos::WithoutInitializing;
37 typename dev_view_t::non_const_type newvals (view_alloc (vals.label(), WithoutInitializing), vals.extent (0));
38 Kokkos::deep_copy(newvals,vals);
39 return newvals;
40}
41
42template<class array_t,class dev_view_t>
43void copy_dev_view_to_host_array(array_t & array, const dev_view_t & dev_view)
44{
45 using host_view_t = typename dev_view_t::HostMirror;
46
47 // Clear out existing and allocate
48 const auto ext = dev_view.extent(0);
49
50 TEUCHOS_TEST_FOR_EXCEPTION(
51 ext != size_t(array.size()), std::logic_error, "Ifpack2::MDF::copy_dev_view_to_host_array: "
52 "Size of permuations on host and device do not match. "
53 "Please report this bug to the Ifpack2 developers.");
54
55 //Wrap array data in view and copy
56 Kokkos::deep_copy(host_view_t(array.get(),ext),dev_view);
57}
58
59template<class scalar_type,class local_ordinal_type,class global_ordinal_type,class node_type>
60void applyReorderingPermutations(
61 const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
62 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y,
63 const Teuchos::ArrayRCP<local_ordinal_type> & perm)
64{
65 TEUCHOS_TEST_FOR_EXCEPTION(X.getNumVectors() != Y.getNumVectors(), std::runtime_error,
66 "Ifpack2::MDF::applyReorderingPermuations ERROR: X.getNumVectors() != Y.getNumVectors().");
67
68 Teuchos::ArrayRCP<Teuchos::ArrayRCP<const scalar_type> > x_ptr = X.get2dView();
69 Teuchos::ArrayRCP<Teuchos::ArrayRCP<scalar_type> > y_ptr = Y.get2dViewNonConst();
70
71 for(size_t k=0; k < X.getNumVectors(); k++)
72 for(local_ordinal_type i=0; (size_t)i< X.getLocalLength(); i++)
73 y_ptr[k][perm[i]] = x_ptr[k][i];
74}
75
76
77template<class scalar_type,class local_ordinal_type,class global_ordinal_type,class node_type>
78auto get_local_crs_row_matrix(
79 Teuchos::RCP<const Tpetra::RowMatrix<scalar_type,local_ordinal_type,global_ordinal_type,node_type>> A_local)
80{
81 using Teuchos::RCP;
82 using Teuchos::rcp;
83 using Teuchos::Array;
84 using Teuchos::rcp_const_cast;
85 using Teuchos::rcp_dynamic_cast;
86
87 using crs_matrix_type = Tpetra::CrsMatrix<scalar_type, local_ordinal_type, global_ordinal_type, node_type>;
88
89 using nonconst_local_inds_host_view_type = typename crs_matrix_type::nonconst_local_inds_host_view_type;
90 using nonconst_values_host_view_type = typename crs_matrix_type::nonconst_values_host_view_type;
91
92 RCP<const crs_matrix_type> A_local_crs = rcp_dynamic_cast<const crs_matrix_type>(A_local);
93 if (A_local_crs.is_null ()) {
94 local_ordinal_type numRows = A_local->getLocalNumRows();
95 Array<size_t> entriesPerRow(numRows);
96 for(local_ordinal_type i = 0; i < numRows; i++) {
97 entriesPerRow[i] = A_local->getNumEntriesInLocalRow(i);
98 }
99 RCP<crs_matrix_type> A_local_crs_nc =
100 rcp (new crs_matrix_type (A_local->getRowMap (),
101 A_local->getColMap (),
102 entriesPerRow()));
103 // copy entries into A_local_crs
104 nonconst_local_inds_host_view_type indices("indices",A_local->getLocalMaxNumRowEntries());
105 nonconst_values_host_view_type values("values",A_local->getLocalMaxNumRowEntries());
106 for(local_ordinal_type i = 0; i < numRows; i++) {
107 size_t numEntries = 0;
108 A_local->getLocalRowCopy(i, indices, values, numEntries);
109 A_local_crs_nc->insertLocalValues(i, numEntries, reinterpret_cast<scalar_type*>(values.data()), indices.data());
110 }
111 A_local_crs_nc->fillComplete (A_local->getDomainMap (), A_local->getRangeMap ());
112 A_local_crs = rcp_const_cast<const crs_matrix_type> (A_local_crs_nc);
113 }
114
115 return A_local_crs;
116}
117
118
119
120}
121
122}
123
124template<class MatrixType>
125MDF<MatrixType>::MDF (const Teuchos::RCP<const row_matrix_type>& Matrix_in)
126 : A_ (Matrix_in),
127 Verbosity_(0),
128 LevelOfFill_ (0),
129 Overalloc_ (2.),
130 isAllocated_ (false),
131 isInitialized_ (false),
132 isComputed_ (false),
133 numInitialize_ (0),
134 numCompute_ (0),
135 numApply_ (0),
136 initializeTime_ (0.0),
137 computeTime_ (0.0),
138 applyTime_ (0.0)
139{
140 allocateSolvers();
141 allocatePermutations();
142}
143
144
145template<class MatrixType>
146MDF<MatrixType>::MDF (const Teuchos::RCP<const crs_matrix_type>& Matrix_in)
147 : A_ (Matrix_in),
148 Verbosity_(0),
149 LevelOfFill_ (0),
150 Overalloc_ (2.),
151 isAllocated_ (false),
152 isInitialized_ (false),
153 isComputed_ (false),
154 numInitialize_ (0),
155 numCompute_ (0),
156 numApply_ (0),
157 initializeTime_ (0.0),
158 computeTime_ (0.0),
159 applyTime_ (0.0)
160{
161 allocateSolvers();
162 allocatePermutations();
163}
164
165template<class MatrixType>
166void MDF<MatrixType>::allocatePermutations (bool force)
167{
168 if (A_.is_null()) return;
169
170 // Allocate arrays as soon as size as known so their pointer is availabe
171 if (force || permutations_.is_null() || A_->getLocalNumRows() != size_t(permutations_.size()))
172 {
173 permutations_ = Teuchos::null;
174 reversePermutations_ = Teuchos::null;
175 permutations_ = permutations_type(A_->getLocalNumRows());
176 reversePermutations_ = permutations_type(A_->getLocalNumRows());
177 }
178}
179
180template<class MatrixType>
181void MDF<MatrixType>::allocateSolvers ()
182{
183 L_solver_ = Teuchos::null;
184 U_solver_ = Teuchos::null;
185 L_solver_ = Teuchos::rcp (new LocalSparseTriangularSolver<row_matrix_type> ());
186 L_solver_->setObjectLabel("lower");
187 U_solver_ = Teuchos::rcp (new LocalSparseTriangularSolver<row_matrix_type> ());
188 U_solver_->setObjectLabel("upper");
189}
190
191template<class MatrixType>
192void
193MDF<MatrixType>::setMatrix (const Teuchos::RCP<const row_matrix_type>& A)
194{
195 // It's legal for A to be null; in that case, you may not call
196 // initialize() until calling setMatrix() with a nonnull input.
197 // Regardless, setting the matrix invalidates any previous
198 // factorization.
199 if (A.getRawPtr () != A_.getRawPtr ()) {
200 isAllocated_ = false;
201 isInitialized_ = false;
202 isComputed_ = false;
203 A_local_ = Teuchos::null;
204 MDF_handle_ = Teuchos::null;
205
206 // The sparse triangular solvers get a triangular factor as their
207 // input matrix. The triangular factors L_ and U_ are getting
208 // reset, so we reset the solvers' matrices to null. Do that
209 // before setting L_ and U_ to null, so that latter step actually
210 // frees the factors.
211 if (! L_solver_.is_null ()) {
212 L_solver_->setMatrix (Teuchos::null);
213 }
214 if (! U_solver_.is_null ()) {
215 U_solver_->setMatrix (Teuchos::null);
216 }
217
218 L_ = Teuchos::null;
219 U_ = Teuchos::null;
220 A_ = A;
221
222 allocatePermutations(true);
223 }
224}
225
226
227
228template<class MatrixType>
231{
232 TEUCHOS_TEST_FOR_EXCEPTION(
233 L_.is_null (), std::runtime_error, "Ifpack2::MDF::getL: The L factor "
234 "is null. Please call initialize() and compute() "
235 "before calling this method. If the input matrix has not yet been set, "
236 "you must first call setMatrix() with a nonnull input matrix before you "
237 "may call initialize() or compute().");
238 return *L_;
239}
240
241template<class MatrixType>
242typename MDF<MatrixType>::permutations_type &
244{
245 TEUCHOS_TEST_FOR_EXCEPTION(
246 permutations_.is_null (), std::runtime_error, "Ifpack2::MDF::getPermutations: "
247 "The permulations are null. Please call initialize() and compute() "
248 "before calling this method. If the input matrix has not yet been set, "
249 "you must first call setMatrix() with a nonnull input matrix before you "
250 "may call initialize() or compute().");
251 return const_cast<permutations_type &>(permutations_);
252}
253template<class MatrixType>
254typename MDF<MatrixType>::permutations_type &
256{
257 TEUCHOS_TEST_FOR_EXCEPTION(
258 reversePermutations_.is_null (), std::runtime_error, "Ifpack2::MDF::getReversePermutations: "
259 "The permulations are null. Please call initialize() and compute() "
260 "before calling this method. If the input matrix has not yet been set, "
261 "you must first call setMatrix() with a nonnull input matrix before you "
262 "may call initialize() or compute().");
263 return const_cast<permutations_type &>(reversePermutations_);
264}
265
266template<class MatrixType>
269{
270 TEUCHOS_TEST_FOR_EXCEPTION(
271 U_.is_null (), std::runtime_error, "Ifpack2::MDF::getU: The U factor "
272 "is null. Please call initialize() and compute() "
273 "before calling this method. If the input matrix has not yet been set, "
274 "you must first call setMatrix() with a nonnull input matrix before you "
275 "may call initialize() or compute().");
276 return *U_;
277}
278
279
280template<class MatrixType>
282 TEUCHOS_TEST_FOR_EXCEPTION(
283 A_.is_null (), std::runtime_error, "Ifpack2::MDF::getNodeSmootherComplexity: "
284 "The input matrix A is null. Please call setMatrix() with a nonnull "
285 "input matrix, then call compute(), before calling this method.");
286 // MDF methods cost roughly one apply + the nnz in the upper+lower triangles
287 if(!L_.is_null() && !U_.is_null())
288 return A_->getLocalNumEntries() + L_->getLocalNumEntries() + U_->getLocalNumEntries();
289 else
290 return 0;
291}
292
293
294template<class MatrixType>
295Teuchos::RCP<const Tpetra::Map<typename MDF<MatrixType>::local_ordinal_type,
299 TEUCHOS_TEST_FOR_EXCEPTION(
300 A_.is_null (), std::runtime_error, "Ifpack2::MDF::getDomainMap: "
301 "The matrix is null. Please call setMatrix() with a nonnull input "
302 "before calling this method.");
303
304 // FIXME (mfh 25 Jan 2014) Shouldn't this just come from A_?
305 TEUCHOS_TEST_FOR_EXCEPTION(
306 L_.is_null (), std::runtime_error, "Ifpack2::MDF::getDomainMap: "
307 "The computed graph is null. Please call initialize() and compute() before calling "
308 "this method.");
309 return L_->getDomainMap ();
310}
311
312
313template<class MatrixType>
314Teuchos::RCP<const Tpetra::Map<typename MDF<MatrixType>::local_ordinal_type,
318 TEUCHOS_TEST_FOR_EXCEPTION(
319 A_.is_null (), std::runtime_error, "Ifpack2::MDF::getRangeMap: "
320 "The matrix is null. Please call setMatrix() with a nonnull input "
321 "before calling this method.");
322
323 // FIXME (mfh 25 Jan 2014) Shouldn't this just come from A_?
324 TEUCHOS_TEST_FOR_EXCEPTION(
325 L_.is_null (), std::runtime_error, "Ifpack2::MDF::getRangeMap: "
326 "The computed graph is null. Please call initialize() abd compute() before calling "
327 "this method.");
328 return L_->getRangeMap ();
329}
330
331template<class MatrixType>
332void
334setParameters (const Teuchos::ParameterList& params)
335{
336 using Teuchos::RCP;
337 using Teuchos::ParameterList;
338 using Teuchos::Array;
339 using Details::getParamTryingTypes;
340 const char prefix[] = "Ifpack2::MDF: ";
341
342 // Default values of the various parameters.
343 int fillLevel = 0;
344 double overalloc = 2.;
345 int verbosity = 0;
346
347 // "fact: mdf level-of-fill" parsing is more complicated, because
348 // we want to allow as many types as make sense. int is the native
349 // type, but we also want to accept double (for backwards
350 // compatibilty with ILUT). You can't cast arbitrary magnitude_type
351 // (e.g., Sacado::MP::Vector) to int, so we use float instead, to
352 // get coverage of the most common magnitude_type cases. Weirdly,
353 // there's an Ifpack2 test that sets the fill level as a
354 // global_ordinal_type.
355 {
356 const std::string paramName ("fact: mdf level-of-fill");
357 getParamTryingTypes<int, int, global_ordinal_type, double, float>
358 (fillLevel, params, paramName, prefix);
359
360 TEUCHOS_TEST_FOR_EXCEPTION
361 (fillLevel != 0, std::runtime_error, prefix << "MDF with level of fill != 0 is not yet implemented.");
362 }
363 {
364 const std::string paramName ("Verbosity");
365 getParamTryingTypes<int, int, global_ordinal_type, double, float>
366 (verbosity, params, paramName, prefix);
367 }
368 {
369 const std::string paramName ("fact: mdf overalloc");
370 getParamTryingTypes<double, double>
371 (overalloc, params, paramName, prefix);
372 }
373
374 // Forward to trisolvers.
375 L_solver_->setParameters(params);
376 U_solver_->setParameters(params);
377
378 // "Commit" the values only after validating all of them. This
379 // ensures that there are no side effects if this routine throws an
380 // exception.
381
382 LevelOfFill_ = fillLevel;
383 Overalloc_ = overalloc;
384 Verbosity_ = verbosity;
385}
386
387
388template<class MatrixType>
389Teuchos::RCP<const typename MDF<MatrixType>::row_matrix_type>
391 return Teuchos::rcp_implicit_cast<const row_matrix_type> (A_);
392}
393
394
395template<class MatrixType>
396Teuchos::RCP<const typename MDF<MatrixType>::crs_matrix_type>
398 return Teuchos::rcp_dynamic_cast<const crs_matrix_type> (A_, true);
399}
400
401
402template<class MatrixType>
403Teuchos::RCP<const typename MDF<MatrixType>::row_matrix_type>
404MDF<MatrixType>::makeLocalFilter (const Teuchos::RCP<const row_matrix_type>& A)
405{
406 using Teuchos::RCP;
407 using Teuchos::rcp;
408 using Teuchos::rcp_dynamic_cast;
409 using Teuchos::rcp_implicit_cast;
410
411 // If A_'s communicator only has one process, or if its column and
412 // row Maps are the same, then it is already local, so use it
413 // directly.
414 if (A->getRowMap ()->getComm ()->getSize () == 1 ||
415 A->getRowMap ()->isSameAs (* (A->getColMap ()))) {
416 return A;
417 }
418
419 // If A_ is already a LocalFilter, then use it directly. This
420 // should be the case if MDF is being used through
421 // AdditiveSchwarz, for example.
422 RCP<const LocalFilter<row_matrix_type> > A_lf_r =
423 rcp_dynamic_cast<const LocalFilter<row_matrix_type> > (A);
424 if (! A_lf_r.is_null ()) {
425 return rcp_implicit_cast<const row_matrix_type> (A_lf_r);
426 }
427 else {
428 // A_'s communicator has more than one process, its row Map and
429 // its column Map differ, and A_ is not a LocalFilter. Thus, we
430 // have to wrap it in a LocalFilter.
431 return rcp (new LocalFilter<row_matrix_type> (A));
432 }
433}
434
435
436template<class MatrixType>
438{
439 using Teuchos::RCP;
440 using Teuchos::rcp;
441 using Teuchos::rcp_const_cast;
442 using Teuchos::rcp_dynamic_cast;
443 using Teuchos::rcp_implicit_cast;
444 using Teuchos::Array;
445 using Teuchos::ArrayView;
446 const char prefix[] = "Ifpack2::MDF::initialize: ";
447
448 TEUCHOS_TEST_FOR_EXCEPTION
449 (A_.is_null (), std::runtime_error, prefix << "The matrix is null. Please "
450 "call setMatrix() with a nonnull input before calling this method.");
451 TEUCHOS_TEST_FOR_EXCEPTION
452 (! A_->isFillComplete (), std::runtime_error, prefix << "The matrix is not "
453 "fill complete. You may not invoke initialize() or compute() with this "
454 "matrix until the matrix is fill complete. If your matrix is a "
455 "Tpetra::CrsMatrix, please call fillComplete on it (with the domain and "
456 "range Maps, if appropriate) before calling this method.");
457
458 Teuchos::Time timer ("MDF::initialize");
459 double startTime = timer.wallTime();
460 { // Start timing
461 Teuchos::TimeMonitor timeMon (timer);
462
463 // Calling initialize() means that the user asserts that the graph
464 // of the sparse matrix may have changed. We must not just reuse
465 // the previous graph in that case.
466 //
467 // Regarding setting isAllocated_ to false: Eventually, we may want
468 // some kind of clever memory reuse strategy, but it's always
469 // correct just to blow everything away and start over.
470 isInitialized_ = false;
471 isAllocated_ = false;
472 isComputed_ = false;
473 MDF_handle_ = Teuchos::null;
474
475 A_local_ = makeLocalFilter (A_);
476 TEUCHOS_TEST_FOR_EXCEPTION(
477 A_local_.is_null (), std::logic_error, "Ifpack2::MDF::initialize: "
478 "makeLocalFilter returned null; it failed to compute A_local. "
479 "Please report this bug to the Ifpack2 developers.");
480
481 // FIXME (mfh 24 Jan 2014, 26 Mar 2014) It would be more efficient
482 // to rewrite MDF so that it works with any RowMatrix input, not
483 // just CrsMatrix. (That would require rewriting mdfGraph to
484 // handle a Tpetra::RowGraph.) However, to make it work for now,
485 // we just copy the input matrix if it's not a CrsMatrix.
486 {
487 RCP<const crs_matrix_type> A_local_crs = Details::MDFImpl::get_local_crs_row_matrix(A_local_);
488
489 auto A_local_device = A_local_crs->getLocalMatrixDevice();
490 MDF_handle_ = rcp( new MDF_handle_device_type(A_local_device) );
491 MDF_handle_->set_verbosity(Verbosity_);
492
493 KokkosSparse::Experimental::mdf_symbolic(A_local_device,*MDF_handle_);
494
495 isAllocated_ = true;
496 }
497
498 checkOrderingConsistency (*A_local_);
499 } // Stop timing
500
501 isInitialized_ = true;
502 ++numInitialize_;
503 initializeTime_ += (timer.wallTime() - startTime);
504}
505
506template<class MatrixType>
507void
510{
511 // First check that the local row map ordering is the same as the local portion of the column map.
512 // The extraction of the strictly lower/upper parts of A, as well as the factorization,
513 // implicitly assume that this is the case.
514 Teuchos::ArrayView<const global_ordinal_type> rowGIDs = A.getRowMap()->getLocalElementList();
515 Teuchos::ArrayView<const global_ordinal_type> colGIDs = A.getColMap()->getLocalElementList();
516 bool gidsAreConsistentlyOrdered=true;
517 global_ordinal_type indexOfInconsistentGID=0;
518 for (global_ordinal_type i=0; i<rowGIDs.size(); ++i) {
519 if (rowGIDs[i] != colGIDs[i]) {
520 gidsAreConsistentlyOrdered=false;
521 indexOfInconsistentGID=i;
522 break;
523 }
524 }
525 TEUCHOS_TEST_FOR_EXCEPTION(gidsAreConsistentlyOrdered==false, std::runtime_error,
526 "The ordering of the local GIDs in the row and column maps is not the same"
527 << std::endl << "at index " << indexOfInconsistentGID
528 << ". Consistency is required, as all calculations are done with"
529 << std::endl << "local indexing.");
530}
531
532template<class MatrixType>
534{
535 using Teuchos::RCP;
536 using Teuchos::rcp;
537 using Teuchos::rcp_const_cast;
538 using Teuchos::rcp_dynamic_cast;
539 using Teuchos::Array;
540 using Teuchos::ArrayView;
541 const char prefix[] = "Ifpack2::MDF::compute: ";
542
543 // initialize() checks this too, but it's easier for users if the
544 // error shows them the name of the method that they actually
545 // called, rather than the name of some internally called method.
546 TEUCHOS_TEST_FOR_EXCEPTION
547 (A_.is_null (), std::runtime_error, prefix << "The matrix is null. Please "
548 "call setMatrix() with a nonnull input before calling this method.");
549 TEUCHOS_TEST_FOR_EXCEPTION
550 (! A_->isFillComplete (), std::runtime_error, prefix << "The matrix is not "
551 "fill complete. You may not invoke initialize() or compute() with this "
552 "matrix until the matrix is fill complete. If your matrix is a "
553 "Tpetra::CrsMatrix, please call fillComplete on it (with the domain and "
554 "range Maps, if appropriate) before calling this method.");
555
556 if (! isInitialized ()) {
557 initialize (); // Don't count this in the compute() time
558 }
559
560 Teuchos::Time timer ("MDF::compute");
561
562 // Start timing
563 Teuchos::TimeMonitor timeMon (timer);
564 double startTime = timer.wallTime();
565
566 isComputed_ = false;
567
568 {//Make sure values in A is picked up even in case of pattern reuse
569 RCP<const crs_matrix_type> A_local_crs = Details::MDFImpl::get_local_crs_row_matrix(A_local_);
570
571 // Compute the ordering and factorize
572 auto A_local_device = A_local_crs->getLocalMatrixDevice();
573
574 KokkosSparse::Experimental::mdf_numeric(A_local_device,*MDF_handle_);
575 }
576
577 // Ordering convention for MDF impl and here are reversed. Do reverse here to avoid confusion
578 Details::MDFImpl::copy_dev_view_to_host_array(reversePermutations_, MDF_handle_->permutation);
579 Details::MDFImpl::copy_dev_view_to_host_array(permutations_, MDF_handle_->permutation_inv);
580
581 // TMR: Need to COPY the values held by the MDF handle because the CRS matrix needs to
582 // exclusively own them and the MDF_handles use_count contribution throws that off
583 {
584 auto L_mdf = MDF_handle_->getL();
585 L_ = rcp(new crs_matrix_type(
586 A_local_->getRowMap (),
587 A_local_->getColMap (),
588 Details::MDFImpl::copy_view(L_mdf.graph.row_map),
589 Details::MDFImpl::copy_view(L_mdf.graph.entries),
590 Details::MDFImpl::copy_view(L_mdf.values)
591 ));
592 }
593 {
594 auto U_mdf = MDF_handle_->getU();
595 U_ = rcp(new crs_matrix_type(
596 A_local_->getRowMap (),
597 A_local_->getColMap (),
598 Details::MDFImpl::copy_view(U_mdf.graph.row_map),
599 Details::MDFImpl::copy_view(U_mdf.graph.entries),
600 Details::MDFImpl::copy_view(U_mdf.values)
601 ));
602 }
603 L_->fillComplete ();
604 U_->fillComplete ();
605 L_solver_->setMatrix (L_);
606 L_solver_->initialize ();
607 L_solver_->compute ();
608 U_solver_->setMatrix (U_);
609 U_solver_->initialize ();
610 U_solver_->compute ();
611
612 isComputed_ = true;
613 ++numCompute_;
614 computeTime_ += (timer.wallTime() - startTime);
615}
616
617template<class MatrixType>
618void
620apply_impl (const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
621 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y,
622 Teuchos::ETransp mode,
623 scalar_type alpha,
624 scalar_type beta) const
625{
626 const scalar_type one = STS::one ();
627 const scalar_type zero = STS::zero ();
628
629 if (alpha == one && beta == zero) {
630 MV tmp (Y.getMap (), Y.getNumVectors ());
631 Details::MDFImpl::applyReorderingPermutations(X,tmp,permutations_);
632 if (mode == Teuchos::NO_TRANS) { // Solve L (D (U Y)) = X for Y.
633 // Start by solving L Y = X for Y.
634 L_solver_->apply (tmp, Y, mode);
635 U_solver_->apply (Y, tmp, mode); // Solve U Y = Y.
636 }
637 else { // Solve U^P (D^P (L^P Y)) = X for Y (where P is * or T).
638 // Start by solving U^P Y = X for Y.
639 U_solver_->apply (tmp, Y, mode);
640 L_solver_->apply (Y, tmp, mode); // Solve L^P Y = Y.
641 }
642 Details::MDFImpl::applyReorderingPermutations(tmp,Y,reversePermutations_);
643 }
644 else { // alpha != 1 or beta != 0
645 if (alpha == zero) {
646 // The special case for beta == 0 ensures that if Y contains Inf
647 // or NaN values, we replace them with 0 (following BLAS
648 // convention), rather than multiplying them by 0 to get NaN.
649 if (beta == zero) {
650 Y.putScalar (zero);
651 } else {
652 Y.scale (beta);
653 }
654 } else { // alpha != zero
655 MV Y_tmp (Y.getMap (), Y.getNumVectors ());
656 apply_impl (X, Y_tmp, mode);
657 Y.update (alpha, Y_tmp, beta);
658 }
659 }
660}
661
662template<class MatrixType>
663void
665apply (const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
666 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y,
667 Teuchos::ETransp mode,
668 scalar_type alpha,
669 scalar_type beta) const
670{
671 using Teuchos::RCP;
672 using Teuchos::rcpFromRef;
673
674 TEUCHOS_TEST_FOR_EXCEPTION(
675 A_.is_null (), std::runtime_error, "Ifpack2::MDF::apply: The matrix is "
676 "null. Please call setMatrix() with a nonnull input, then initialize() "
677 "and compute(), before calling this method.");
678 TEUCHOS_TEST_FOR_EXCEPTION(
679 ! isComputed (), std::runtime_error,
680 "Ifpack2::MDF::apply: If you have not yet called compute(), "
681 "you must call compute() before calling this method.");
682 TEUCHOS_TEST_FOR_EXCEPTION(
683 X.getNumVectors () != Y.getNumVectors (), std::invalid_argument,
684 "Ifpack2::MDF::apply: X and Y do not have the same number of columns. "
685 "X.getNumVectors() = " << X.getNumVectors ()
686 << " != Y.getNumVectors() = " << Y.getNumVectors () << ".");
687 TEUCHOS_TEST_FOR_EXCEPTION(
688 STS::isComplex && mode == Teuchos::CONJ_TRANS, std::logic_error,
689 "Ifpack2::MDF::apply: mode = Teuchos::CONJ_TRANS is not implemented for "
690 "complex Scalar type. Please talk to the Ifpack2 developers to get this "
691 "fixed. There is a FIXME in this file about this very issue.");
692#ifdef HAVE_IFPACK2_DEBUG
693 {
694 Teuchos::Array<magnitude_type> norms (X.getNumVectors ());
695 X.norm1 (norms ());
696 bool good = true;
697 for (size_t j = 0; j < X.getNumVectors (); ++j) {
698 if (STM::isnaninf (norms[j])) {
699 good = false;
700 break;
701 }
702 }
703 TEUCHOS_TEST_FOR_EXCEPTION( ! good, std::runtime_error, "Ifpack2::MDF::apply: The 1-norm of the input X is NaN or Inf.");
704 }
705#endif // HAVE_IFPACK2_DEBUG
706
707 Teuchos::Time timer ("MDF::apply");
708 double startTime = timer.wallTime();
709 { // Start timing
710 Teuchos::TimeMonitor timeMon (timer);
711 apply_impl(X,Y,mode,alpha,beta);
712 }//end timing
713
714#ifdef HAVE_IFPACK2_DEBUG
715 {
716 Teuchos::Array<magnitude_type> norms (Y.getNumVectors ());
717 Y.norm1 (norms ());
718 bool good = true;
719 for (size_t j = 0; j < Y.getNumVectors (); ++j) {
720 if (STM::isnaninf (norms[j])) {
721 good = false;
722 break;
723 }
724 }
725 TEUCHOS_TEST_FOR_EXCEPTION( ! good, std::runtime_error, "Ifpack2::MDF::apply: The 1-norm of the output Y is NaN or Inf.");
726 }
727#endif // HAVE_IFPACK2_DEBUG
728
729 ++numApply_;
730 applyTime_ += (timer.wallTime() - startTime);
731}
732
733template<class MatrixType>
735{
736 std::ostringstream os;
737
738 // Output is a valid YAML dictionary in flow style. If you don't
739 // like everything on a single line, you should call describe()
740 // instead.
741 os << "\"Ifpack2::MDF\": {";
742 os << "Initialized: " << (isInitialized () ? "true" : "false") << ", "
743 << "Computed: " << (isComputed () ? "true" : "false") << ", ";
744
745 os << "Level-of-fill: " << getLevelOfFill() << ", ";
746
747 if (A_.is_null ()) {
748 os << "Matrix: null";
749 }
750 else {
751 os << "Global matrix dimensions: ["
752 << A_->getGlobalNumRows () << ", " << A_->getGlobalNumCols () << "]"
753 << ", Global nnz: " << A_->getGlobalNumEntries();
754 }
755
756 if (! L_solver_.is_null ()) os << ", " << L_solver_->description ();
757 if (! U_solver_.is_null ()) os << ", " << U_solver_->description ();
758
759 os << "}";
760 return os.str ();
761}
762
763} // namespace Ifpack2
764
765#define IFPACK2_MDF_INSTANT(S,LO,GO,N) \
766 template class Ifpack2::MDF< Tpetra::RowMatrix<S, LO, GO, N> >;
767
768#endif /* IFPACK2_MDF_DEF_HPP */
Ifpack2::ScalingType enumerable type.
Access only local rows and columns of a sparse matrix.
Definition Ifpack2_LocalFilter_decl.hpp:131
"Preconditioner" that solves local sparse triangular systems.
Definition Ifpack2_LocalSparseTriangularSolver_decl.hpp:55
MDF (incomplete LU factorization with minimum discarded fill reordering) of a Tpetra sparse matrix.
Definition Ifpack2_MDF_decl.hpp:59
void setParameters(const Teuchos::ParameterList &params)
Definition Ifpack2_MDF_def.hpp:334
const crs_matrix_type & getL() const
Return the L factor of the MDF factorization.
Definition Ifpack2_MDF_def.hpp:230
void apply(const Tpetra::MultiVector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > &X, Tpetra::MultiVector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, scalar_type alpha=Teuchos::ScalarTraits< scalar_type >::one(), scalar_type beta=Teuchos::ScalarTraits< scalar_type >::zero()) const
Apply the (inverse of the) incomplete factorization to X, resulting in Y.
Definition Ifpack2_MDF_def.hpp:665
int getLevelOfFill() const
Get level of fill (the "k" in ILU(k)).
Definition Ifpack2_MDF_decl.hpp:335
Teuchos::RCP< LocalSparseTriangularSolver< row_matrix_type > > L_solver_
Sparse triangular solver for L.
Definition Ifpack2_MDF_decl.hpp:403
bool isComputed() const
Whether compute() has been called on this object.
Definition Ifpack2_MDF_decl.hpp:176
permutations_type & getReversePermutations() const
Return the reverse permutations of the MDF factorization.
Definition Ifpack2_MDF_def.hpp:255
permutations_type reversePermutations_
The reverse permuations from MDF factorization.
Definition Ifpack2_MDF_decl.hpp:413
const crs_matrix_type & getU() const
Return the U factor of the MDF factorization.
Definition Ifpack2_MDF_def.hpp:268
void compute()
Compute the (numeric) incomplete factorization.
Definition Ifpack2_MDF_def.hpp:533
Teuchos::RCP< crs_matrix_type > U_
The U (upper triangular) factor of ILU(k).
Definition Ifpack2_MDF_decl.hpp:405
bool isInitialized() const
Whether initialize() has been called on this object.
Definition Ifpack2_MDF_decl.hpp:172
Teuchos::RCP< const row_matrix_type > A_local_
The matrix whos numbers are used to to compute ILU(k). The graph may be computed using a crs_matrix_t...
Definition Ifpack2_MDF_decl.hpp:395
void initialize()
Initialize by computing the symbolic incomplete factorization.
Definition Ifpack2_MDF_def.hpp:437
std::string description() const
A one-line description of this object.
Definition Ifpack2_MDF_def.hpp:734
Teuchos::RCP< const Tpetra::Map< local_ordinal_type, global_ordinal_type, node_type > > getRangeMap() const
Returns the Tpetra::Map object associated with the range of this operator.
Definition Ifpack2_MDF_def.hpp:317
Tpetra::CrsMatrix< scalar_type, local_ordinal_type, global_ordinal_type, node_type > crs_matrix_type
Tpetra::CrsMatrix specialization used by this class for representing L and U.
Definition Ifpack2_MDF_decl.hpp:95
Teuchos::RCP< const Tpetra::Map< local_ordinal_type, global_ordinal_type, node_type > > getDomainMap() const
Returns the Tpetra::Map object associated with the domain of this operator.
Definition Ifpack2_MDF_def.hpp:298
Teuchos::RCP< const crs_matrix_type > getCrsMatrix() const
Return the input matrix A as a Tpetra::CrsMatrix, if possible; else throws.
Definition Ifpack2_MDF_def.hpp:397
MatrixType::node_type node_type
The Node type used by the input MatrixType.
Definition Ifpack2_MDF_decl.hpp:71
MatrixType::global_ordinal_type global_ordinal_type
The type of global indices in the input MatrixType.
Definition Ifpack2_MDF_decl.hpp:68
permutations_type permutations_
The computed permuations from MDF factorization.
Definition Ifpack2_MDF_decl.hpp:410
Teuchos::RCP< LocalSparseTriangularSolver< row_matrix_type > > U_solver_
Sparse triangular solver for U.
Definition Ifpack2_MDF_decl.hpp:407
MatrixType::scalar_type scalar_type
The type of the entries of the input MatrixType.
Definition Ifpack2_MDF_decl.hpp:62
Teuchos::RCP< const row_matrix_type > getMatrix() const
Get the input matrix.
Definition Ifpack2_MDF_def.hpp:390
virtual void setMatrix(const Teuchos::RCP< const row_matrix_type > &A)
Change the matrix to be preconditioned.
Definition Ifpack2_MDF_def.hpp:193
size_t getNodeSmootherComplexity() const
Get a rough estimate of cost per iteration.
Definition Ifpack2_MDF_def.hpp:281
Teuchos::RCP< crs_matrix_type > L_
The L (lower triangular) factor of ILU(k).
Definition Ifpack2_MDF_decl.hpp:401
Teuchos::RCP< const row_matrix_type > A_
The (original) input matrix for which to compute ILU(k).
Definition Ifpack2_MDF_decl.hpp:387
permutations_type & getPermutations() const
Return the permutations of the MDF factorization.
Definition Ifpack2_MDF_def.hpp:243
Preconditioners and smoothers for Tpetra sparse matrices.
Definition Ifpack2_AdditiveSchwarz_decl.hpp:41