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