Reference documentation for deal.II version GIT 29f9da0a34 2023-12-07 10:00:01+00:00
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
trilinos_sparse_matrix.cc
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2008 - 2023 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.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
18 
19 #ifdef DEAL_II_WITH_TRILINOS
20 
22 # include <deal.II/base/utilities.h>
23 
30 
31 # include <boost/container/small_vector.hpp>
32 
33 # ifdef DEAL_II_TRILINOS_WITH_EPETRAEXT
34 # include <EpetraExt_MatrixMatrix.h>
35 # endif
36 # include <Epetra_Export.h>
37 # include <Teuchos_RCP.hpp>
38 # include <ml_epetra_utils.h>
39 # include <ml_struct.h>
40 
41 # include <memory>
42 
44 
45 namespace TrilinosWrappers
46 {
47  namespace internal
48  {
49  template <typename VectorType>
50  typename VectorType::value_type *
51  begin(VectorType &V)
52  {
53  return V.begin();
54  }
55 
56  template <typename VectorType>
57  const typename VectorType::value_type *
58  begin(const VectorType &V)
59  {
60  return V.begin();
61  }
62 
63  template <typename VectorType>
64  typename VectorType::value_type *
65  end(VectorType &V)
66  {
67  return V.end();
68  }
69 
70  template <typename VectorType>
71  const typename VectorType::value_type *
72  end(const VectorType &V)
73  {
74  return V.end();
75  }
76 
77  template <>
78  double *
80  {
81  return V.trilinos_vector()[0];
82  }
83 
84  template <>
85  const double *
87  {
88  return V.trilinos_vector()[0];
89  }
90 
91  template <>
92  double *
94  {
95  return V.trilinos_vector()[0] + V.trilinos_vector().MyLength();
96  }
97 
98  template <>
99  const double *
101  {
102  return V.trilinos_vector()[0] + V.trilinos_vector().MyLength();
103  }
104 
105 # ifdef DEAL_II_TRILINOS_WITH_TPETRA
106  template <typename Number>
107  Number *
109  {
110  return V.trilinos_vector().getDataNonConst().get();
111  }
112 
113  template <typename Number>
114  const Number *
116  {
117  return V.trilinos_vector().getData().get();
118  }
119 
120  template <typename Number>
121  Number *
123  {
124  return V.trilinos_vector().getDataNonConst().get() +
125  V.trilinos_vector().getLocalLength();
126  }
127 
128  template <typename Number>
129  const Number *
131  {
132  return V.trilinos_vector().getData().get() +
133  V.trilinos_vector().getLocalLength();
134  }
135 # endif
136  } // namespace internal
137 
138 
139  namespace SparseMatrixIterators
140  {
141  void
143  {
144  // if we are asked to visit the past-the-end line, then simply
145  // release all our caches and go on with life.
146  //
147  // do the same if the row we're supposed to visit is not locally
148  // owned. this is simply going to make non-locally owned rows
149  // look like they're empty
150  if ((this->a_row == matrix->m()) ||
151  (matrix->in_local_range(this->a_row) == false))
152  {
153  colnum_cache.reset();
154  value_cache.reset();
155 
156  return;
157  }
158 
159  // get a representation of the present row
160  int ncols;
162  if (value_cache.get() == nullptr)
163  {
164  value_cache =
165  std::make_shared<std::vector<TrilinosScalar>>(matrix->n());
166  colnum_cache = std::make_shared<std::vector<size_type>>(matrix->n());
167  }
168  else
169  {
170  value_cache->resize(matrix->n());
171  colnum_cache->resize(matrix->n());
172  }
173 
174  int ierr = matrix->trilinos_matrix().ExtractGlobalRowCopy(
175  this->a_row,
176  colnums,
177  ncols,
178  value_cache->data(),
179  reinterpret_cast<TrilinosWrappers::types::int_type *>(
180  colnum_cache->data()));
181  value_cache->resize(ncols);
182  colnum_cache->resize(ncols);
183  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
184 
185  // copy it into our caches if the
186  // line isn't empty. if it is, then
187  // we've done something wrong, since
188  // we shouldn't have initialized an
189  // iterator for an empty line (what
190  // would it point to?)
191  }
192  } // namespace SparseMatrixIterators
193 
194 
195  // The constructor is actually the
196  // only point where we have to check
197  // whether we build a serial or a
198  // parallel Trilinos matrix.
199  // Actually, it does not even matter
200  // how many threads there are, but
201  // only if we use an MPI compiler or
202  // a standard compiler. So, even one
203  // thread on a configuration with
204  // MPI will still get a parallel
205  // interface.
207  : column_space_map(new Epetra_Map(0, 0, Utilities::Trilinos::comm_self()))
208  , matrix(
209  new Epetra_FECrsMatrix(View, *column_space_map, *column_space_map, 0))
210  , last_action(Zero)
211  , compressed(true)
212  {
213  matrix->FillComplete();
214  }
215 
216 
217 
219  const size_type n,
220  const unsigned int n_max_entries_per_row)
221  : column_space_map(
222  new Epetra_Map(static_cast<TrilinosWrappers::types::int_type>(n),
223  0,
224  Utilities::Trilinos::comm_self()))
225  ,
226 
227  // on one processor only, we know how the
228  // columns of the matrix will be
229  // distributed (everything on one
230  // processor), so we can hand in this
231  // information to the constructor. we
232  // can't do so in parallel, where the
233  // information from columns is only
234  // available when entries have been added
235  matrix(new Epetra_FECrsMatrix(
236  Copy,
237  Epetra_Map(static_cast<TrilinosWrappers::types::int_type>(m),
238  0,
239  Utilities::Trilinos::comm_self()),
240  *column_space_map,
241  n_max_entries_per_row,
242  false))
243  , last_action(Zero)
244  , compressed(false)
245  {}
246 
247 
248 
250  const size_type n,
251  const std::vector<unsigned int> &n_entries_per_row)
252  : column_space_map(
253  new Epetra_Map(static_cast<TrilinosWrappers::types::int_type>(n),
254  0,
255  Utilities::Trilinos::comm_self()))
256  , matrix(new Epetra_FECrsMatrix(
257  Copy,
258  Epetra_Map(static_cast<TrilinosWrappers::types::int_type>(m),
259  0,
260  Utilities::Trilinos::comm_self()),
261  *column_space_map,
262  reinterpret_cast<int *>(
263  const_cast<unsigned int *>(n_entries_per_row.data())),
264  false))
265  , last_action(Zero)
266  , compressed(false)
267  {}
268 
269 
270 
271  SparseMatrix::SparseMatrix(const IndexSet &parallel_partitioning,
272  const MPI_Comm communicator,
273  const unsigned int n_max_entries_per_row)
274  : column_space_map(new Epetra_Map(
275  parallel_partitioning.make_trilinos_map(communicator, false)))
276  , matrix(new Epetra_FECrsMatrix(Copy,
277  *column_space_map,
278  n_max_entries_per_row,
279  false))
280  , last_action(Zero)
281  , compressed(false)
282  {}
283 
284 
285 
286  SparseMatrix::SparseMatrix(const IndexSet &parallel_partitioning,
287  const MPI_Comm communicator,
288  const std::vector<unsigned int> &n_entries_per_row)
289  : column_space_map(new Epetra_Map(
290  parallel_partitioning.make_trilinos_map(communicator, false)))
291  , matrix(new Epetra_FECrsMatrix(Copy,
292  *column_space_map,
293  reinterpret_cast<int *>(
294  const_cast<unsigned int *>(
295  n_entries_per_row.data())),
296  false))
297  , last_action(Zero)
298  , compressed(false)
299  {}
300 
301 
302 
303  SparseMatrix::SparseMatrix(const IndexSet &row_parallel_partitioning,
304  const IndexSet &col_parallel_partitioning,
305  const MPI_Comm communicator,
306  const size_type n_max_entries_per_row)
307  : column_space_map(new Epetra_Map(
308  col_parallel_partitioning.make_trilinos_map(communicator, false)))
309  , matrix(new Epetra_FECrsMatrix(
310  Copy,
311  row_parallel_partitioning.make_trilinos_map(communicator, false),
312  n_max_entries_per_row,
313  false))
314  , last_action(Zero)
315  , compressed(false)
316  {}
317 
318 
319 
320  SparseMatrix::SparseMatrix(const IndexSet &row_parallel_partitioning,
321  const IndexSet &col_parallel_partitioning,
322  const MPI_Comm communicator,
323  const std::vector<unsigned int> &n_entries_per_row)
324  : column_space_map(new Epetra_Map(
325  col_parallel_partitioning.make_trilinos_map(communicator, false)))
326  , matrix(new Epetra_FECrsMatrix(
327  Copy,
328  row_parallel_partitioning.make_trilinos_map(communicator, false),
329  reinterpret_cast<int *>(
330  const_cast<unsigned int *>(n_entries_per_row.data())),
331  false))
332  , last_action(Zero)
333  , compressed(false)
334  {}
335 
336 
337 
339  : column_space_map(new Epetra_Map(sparsity_pattern.domain_partitioner()))
340  , matrix(
341  new Epetra_FECrsMatrix(Copy,
342  sparsity_pattern.trilinos_sparsity_pattern(),
343  false))
344  , last_action(Zero)
345  , compressed(true)
346  {
347  Assert(sparsity_pattern.trilinos_sparsity_pattern().Filled() == true,
348  ExcMessage(
349  "The Trilinos sparsity pattern has not been compressed."));
351  }
352 
353 
354 
356  : column_space_map(std::move(other.column_space_map))
357  , matrix(std::move(other.matrix))
358  , nonlocal_matrix(std::move(other.nonlocal_matrix))
359  , nonlocal_matrix_exporter(std::move(other.nonlocal_matrix_exporter))
360  , last_action(other.last_action)
361  , compressed(other.compressed)
362  {
363  other.last_action = Zero;
364  other.compressed = false;
365  }
366 
367 
368 
369  void
371  {
372  if (this == &rhs)
373  return;
374 
375  nonlocal_matrix.reset();
376  nonlocal_matrix_exporter.reset();
377 
378  // check whether we need to update the whole matrix layout (we have
379  // different maps or if we detect a row where the columns of the two
380  // matrices do not match)
381  bool needs_deep_copy =
382  !matrix->RowMap().SameAs(rhs.matrix->RowMap()) ||
383  !matrix->ColMap().SameAs(rhs.matrix->ColMap()) ||
384  !matrix->DomainMap().SameAs(rhs.matrix->DomainMap()) ||
386  if (!needs_deep_copy)
387  {
388  // Try to copy all the rows of the matrix one by one. In case of error
389  // (i.e., the column indices are different), we need to abort and blow
390  // away the matrix.
391  for (const auto row : locally_owned_range_indices())
392  {
393  const int row_local = matrix->RowMap().LID(
394  static_cast<TrilinosWrappers::types::int_type>(row));
395  Assert((row_local >= 0), ExcAccessToNonlocalRow(row));
396 
397  int n_entries, rhs_n_entries;
398  TrilinosScalar *value_ptr, *rhs_value_ptr;
399  int *index_ptr, *rhs_index_ptr;
400  int ierr = rhs.matrix->ExtractMyRowView(row_local,
401  rhs_n_entries,
402  rhs_value_ptr,
403  rhs_index_ptr);
404  (void)ierr;
405  Assert(ierr == 0, ExcTrilinosError(ierr));
406 
407  ierr = matrix->ExtractMyRowView(row_local,
408  n_entries,
409  value_ptr,
410  index_ptr);
411  Assert(ierr == 0, ExcTrilinosError(ierr));
412 
413  if (n_entries != rhs_n_entries ||
414  std::memcmp(static_cast<void *>(index_ptr),
415  static_cast<void *>(rhs_index_ptr),
416  sizeof(int) * n_entries) != 0)
417  {
418  needs_deep_copy = true;
419  break;
420  }
421 
422  for (int i = 0; i < n_entries; ++i)
423  value_ptr[i] = rhs_value_ptr[i];
424  }
425  }
426 
427  if (needs_deep_copy)
428  {
430  std::make_unique<Epetra_Map>(rhs.trilinos_matrix().DomainMap());
431 
432  // release memory before reallocation
433  matrix = std::make_unique<Epetra_FECrsMatrix>(*rhs.matrix);
434 
435  matrix->FillComplete(*column_space_map, matrix->RowMap());
436  }
437 
438  if (rhs.nonlocal_matrix.get() != nullptr)
440  std::make_unique<Epetra_CrsMatrix>(Copy, rhs.nonlocal_matrix->Graph());
441  }
442 
443 
444 
445  namespace
446  {
448 
449  template <typename SparsityPatternType>
450  void
451  reinit_matrix(const IndexSet &row_parallel_partitioning,
452  const IndexSet &column_parallel_partitioning,
453  const SparsityPatternType &sparsity_pattern,
454  const bool exchange_data,
455  const MPI_Comm communicator,
456  std::unique_ptr<Epetra_Map> &column_space_map,
457  std::unique_ptr<Epetra_FECrsMatrix> &matrix,
458  std::unique_ptr<Epetra_CrsMatrix> &nonlocal_matrix,
459  std::unique_ptr<Epetra_Export> &nonlocal_matrix_exporter)
460  {
461  // release memory before reallocation
462  matrix.reset();
463  nonlocal_matrix.reset();
464  nonlocal_matrix_exporter.reset();
465 
466  column_space_map = std::make_unique<Epetra_Map>(
467  column_parallel_partitioning.make_trilinos_map(communicator, false));
468 
469  if (column_space_map->Comm().MyPID() == 0)
470  {
471  AssertDimension(sparsity_pattern.n_rows(),
472  row_parallel_partitioning.size());
473  AssertDimension(sparsity_pattern.n_cols(),
474  column_parallel_partitioning.size());
475  }
476 
477  Epetra_Map row_space_map =
478  row_parallel_partitioning.make_trilinos_map(communicator, false);
479 
480  // if we want to exchange data, build a usual Trilinos sparsity pattern
481  // and let that handle the exchange. otherwise, manually create a
482  // CrsGraph, which consumes considerably less memory because it can set
483  // correct number of indices right from the start
484  if (exchange_data)
485  {
486  SparsityPattern trilinos_sparsity;
487  trilinos_sparsity.reinit(row_parallel_partitioning,
488  column_parallel_partitioning,
489  sparsity_pattern,
490  communicator,
491  exchange_data);
492  matrix = std::make_unique<Epetra_FECrsMatrix>(
493  Copy, trilinos_sparsity.trilinos_sparsity_pattern(), false);
494 
495  return;
496  }
497 
498  const size_type first_row = TrilinosWrappers::min_my_gid(row_space_map),
499  last_row =
500  TrilinosWrappers::max_my_gid(row_space_map) + 1;
501  std::vector<int> n_entries_per_row(last_row - first_row);
502 
503  for (size_type row = first_row; row < last_row; ++row)
504  n_entries_per_row[row - first_row] = sparsity_pattern.row_length(row);
505 
506  // The deal.II notation of a Sparsity pattern corresponds to the Epetra
507  // concept of a Graph. Hence, we generate a graph by copying the
508  // sparsity pattern into it, and then build up the matrix from the
509  // graph. This is considerable faster than directly filling elements
510  // into the matrix. Moreover, it consumes less memory, since the
511  // internal reordering is done on ints only, and we can leave the
512  // doubles aside.
513 
514  // for more than one processor, need to specify only row map first and
515  // let the matrix entries decide about the column map (which says which
516  // columns are present in the matrix, not to be confused with the
517  // col_map that tells how the domain dofs of the matrix will be
518  // distributed). for only one processor, we can directly assign the
519  // columns as well. Compare this with bug # 4123 in the Sandia Bugzilla.
520  std::unique_ptr<Epetra_CrsGraph> graph;
521  if (row_space_map.Comm().NumProc() > 1)
522  graph = std::make_unique<Epetra_CrsGraph>(Copy,
523  row_space_map,
524  n_entries_per_row.data(),
525  true);
526  else
527  graph = std::make_unique<Epetra_CrsGraph>(Copy,
528  row_space_map,
529  *column_space_map,
530  n_entries_per_row.data(),
531  true);
532 
533  // This functions assumes that the sparsity pattern sits on all
534  // processors (completely). The parallel version uses an Epetra graph
535  // that is already distributed.
536 
537  // now insert the indices
538  std::vector<TrilinosWrappers::types::int_type> row_indices;
539 
540  for (size_type row = first_row; row < last_row; ++row)
541  {
542  const int row_length = sparsity_pattern.row_length(row);
543  if (row_length == 0)
544  continue;
545 
546  row_indices.resize(row_length, -1);
547  {
548  typename SparsityPatternType::iterator p =
549  sparsity_pattern.begin(row);
550  for (size_type col = 0; p != sparsity_pattern.end(row); ++p, ++col)
551  row_indices[col] = p->column();
552  }
553  graph->Epetra_CrsGraph::InsertGlobalIndices(row,
554  row_length,
555  row_indices.data());
556  }
557 
558  // Eventually, optimize the graph structure (sort indices, make memory
559  // contiguous, etc). note that the documentation of the function indeed
560  // states that we first need to provide the column (domain) map and then
561  // the row (range) map
562  graph->FillComplete(*column_space_map, row_space_map);
563  graph->OptimizeStorage();
564 
565  // check whether we got the number of columns right.
566  AssertDimension(sparsity_pattern.n_cols(),
568  (void)n_global_cols;
569 
570  // And now finally generate the matrix.
571  matrix = std::make_unique<Epetra_FECrsMatrix>(Copy, *graph, false);
572  }
573 
574 
575 
576  // for the non-local graph, we need to circumvent the problem that some
577  // processors will not add into the non-local graph at all: We do not want
578  // to insert dummy elements on >5000 processors because that gets very
579  // slow. Thus, we set a flag in Epetra_CrsGraph that sets the correct
580  // flag. Since it is protected, we need to expose this information by
581  // deriving a class from Epetra_CrsGraph for the purpose of creating the
582  // data structure
583  class Epetra_CrsGraphMod : public Epetra_CrsGraph
584  {
585  public:
586  Epetra_CrsGraphMod(const Epetra_Map &row_map,
587  const int *n_entries_per_row)
588  : Epetra_CrsGraph(Copy, row_map, n_entries_per_row, true)
589  {}
590 
591  void
592  SetIndicesAreGlobal()
593  {
594  this->Epetra_CrsGraph::SetIndicesAreGlobal(true);
595  }
596  };
597 
598 
599 
600  // specialization for DynamicSparsityPattern which can provide us with
601  // more information about the non-locally owned rows
602  template <>
603  void
604  reinit_matrix(const IndexSet &row_parallel_partitioning,
605  const IndexSet &column_parallel_partitioning,
606  const DynamicSparsityPattern &sparsity_pattern,
607  const bool exchange_data,
608  const MPI_Comm communicator,
609  std::unique_ptr<Epetra_Map> &column_space_map,
610  std::unique_ptr<Epetra_FECrsMatrix> &matrix,
611  std::unique_ptr<Epetra_CrsMatrix> &nonlocal_matrix,
612  std::unique_ptr<Epetra_Export> &nonlocal_matrix_exporter)
613  {
614  matrix.reset();
615  nonlocal_matrix.reset();
616  nonlocal_matrix_exporter.reset();
617 
618  column_space_map = std::make_unique<Epetra_Map>(
619  column_parallel_partitioning.make_trilinos_map(communicator, false));
620 
621  AssertDimension(sparsity_pattern.n_rows(),
622  row_parallel_partitioning.size());
623  AssertDimension(sparsity_pattern.n_cols(),
624  column_parallel_partitioning.size());
625 
626  Epetra_Map row_space_map =
627  row_parallel_partitioning.make_trilinos_map(communicator, false);
628 
629  IndexSet relevant_rows(sparsity_pattern.row_index_set());
630  // serial case
631  if (relevant_rows.size() == 0)
632  {
633  relevant_rows.set_size(
634  TrilinosWrappers::n_global_elements(row_space_map));
635  relevant_rows.add_range(
636  0, TrilinosWrappers::n_global_elements(row_space_map));
637  }
638  relevant_rows.compress();
639  Assert(relevant_rows.n_elements() >=
640  static_cast<unsigned int>(row_space_map.NumMyElements()),
641  ExcMessage(
642  "Locally relevant rows of sparsity pattern must contain "
643  "all locally owned rows"));
644 
645  // check whether the relevant rows correspond to exactly the same map as
646  // the owned rows. In that case, do not create the nonlocal graph and
647  // fill the columns by demand
648  const bool have_ghost_rows = [&]() {
649  const std::vector<::types::global_dof_index> indices =
650  relevant_rows.get_index_vector();
651  Epetra_Map relevant_map(
653  TrilinosWrappers::types::int_type(relevant_rows.n_elements()),
654  (indices.empty() ?
655  nullptr :
656  reinterpret_cast<const TrilinosWrappers::types::int_type *>(
657  indices.data())),
658  0,
659  row_space_map.Comm());
660  return !relevant_map.SameAs(row_space_map);
661  }();
662 
663  const unsigned int n_rows = relevant_rows.n_elements();
664  std::vector<TrilinosWrappers::types::int_type> ghost_rows;
665  std::vector<int> n_entries_per_row(row_space_map.NumMyElements());
666  std::vector<int> n_entries_per_ghost_row;
667  for (unsigned int i = 0, own = 0; i < n_rows; ++i)
668  {
669  const TrilinosWrappers::types::int_type global_row =
670  relevant_rows.nth_index_in_set(i);
671  if (row_space_map.MyGID(global_row))
672  n_entries_per_row[own++] = sparsity_pattern.row_length(global_row);
673  else if (sparsity_pattern.row_length(global_row) > 0)
674  {
675  ghost_rows.push_back(global_row);
676  n_entries_per_ghost_row.push_back(
677  sparsity_pattern.row_length(global_row));
678  }
679  }
680 
681  Epetra_Map off_processor_map(-1,
682  ghost_rows.size(),
683  (ghost_rows.size() > 0) ?
684  (ghost_rows.data()) :
685  nullptr,
686  0,
687  row_space_map.Comm());
688 
689  std::unique_ptr<Epetra_CrsGraph> graph;
690  std::unique_ptr<Epetra_CrsGraphMod> nonlocal_graph;
691  if (row_space_map.Comm().NumProc() > 1)
692  {
693  graph =
694  std::make_unique<Epetra_CrsGraph>(Copy,
695  row_space_map,
696  (n_entries_per_row.size() > 0) ?
697  (n_entries_per_row.data()) :
698  nullptr,
699  exchange_data ? false : true);
700  if (have_ghost_rows == true)
701  nonlocal_graph = std::make_unique<Epetra_CrsGraphMod>(
702  off_processor_map, n_entries_per_ghost_row.data());
703  }
704  else
705  graph =
706  std::make_unique<Epetra_CrsGraph>(Copy,
707  row_space_map,
708  *column_space_map,
709  (n_entries_per_row.size() > 0) ?
710  (n_entries_per_row.data()) :
711  nullptr,
712  true);
713 
714  // now insert the indices, select between the right matrix
715  std::vector<TrilinosWrappers::types::int_type> row_indices;
716 
717  for (unsigned int i = 0; i < n_rows; ++i)
718  {
719  const TrilinosWrappers::types::int_type global_row =
720  relevant_rows.nth_index_in_set(i);
721  const int row_length = sparsity_pattern.row_length(global_row);
722  if (row_length == 0)
723  continue;
724 
725  row_indices.resize(row_length, -1);
726  for (int col = 0; col < row_length; ++col)
727  row_indices[col] = sparsity_pattern.column_number(global_row, col);
728 
729  if (row_space_map.MyGID(global_row))
730  graph->InsertGlobalIndices(global_row,
731  row_length,
732  row_indices.data());
733  else
734  {
735  Assert(nonlocal_graph.get() != nullptr, ExcInternalError());
736  nonlocal_graph->InsertGlobalIndices(global_row,
737  row_length,
738  row_indices.data());
739  }
740  }
741 
742  // finalize nonlocal graph and create nonlocal matrix
743  if (nonlocal_graph.get() != nullptr)
744  {
745  // must make sure the IndicesAreGlobal flag is set on all processors
746  // because some processors might not call InsertGlobalIndices (and
747  // we do not want to insert dummy indices on all processors for
748  // large-scale simulations due to the bad impact on performance)
749  nonlocal_graph->SetIndicesAreGlobal();
750  Assert(nonlocal_graph->IndicesAreGlobal() == true,
751  ExcInternalError());
752  nonlocal_graph->FillComplete(*column_space_map, row_space_map);
753  nonlocal_graph->OptimizeStorage();
754 
755  // insert data from nonlocal graph into the final sparsity pattern
756  if (exchange_data)
757  {
758  Epetra_Export exporter(nonlocal_graph->RowMap(), row_space_map);
759  int ierr = graph->Export(*nonlocal_graph, exporter, Add);
760  (void)ierr;
761  Assert(ierr == 0, ExcTrilinosError(ierr));
762  }
763 
764  nonlocal_matrix =
765  std::make_unique<Epetra_CrsMatrix>(Copy, *nonlocal_graph);
766  }
767 
768  graph->FillComplete(*column_space_map, row_space_map);
769  graph->OptimizeStorage();
770 
771  AssertDimension(sparsity_pattern.n_cols(),
773 
774  matrix = std::make_unique<Epetra_FECrsMatrix>(Copy, *graph, false);
775  }
776  } // namespace
777 
778 
779 
780  template <typename SparsityPatternType>
781  void
782  SparseMatrix::reinit(const SparsityPatternType &sparsity_pattern)
783  {
784  reinit_matrix(complete_index_set(sparsity_pattern.n_rows()),
785  complete_index_set(sparsity_pattern.n_cols()),
786  sparsity_pattern,
787  false,
788  MPI_COMM_SELF,
790  matrix,
793  }
794 
795 
796 
797  template <typename SparsityPatternType>
798  inline std::enable_if_t<
799  !std::is_same_v<SparsityPatternType, ::SparseMatrix<double>>>
800  SparseMatrix::reinit(const IndexSet &row_parallel_partitioning,
801  const IndexSet &col_parallel_partitioning,
802  const SparsityPatternType &sparsity_pattern,
803  const MPI_Comm communicator,
804  const bool exchange_data)
805  {
806  reinit_matrix(row_parallel_partitioning,
807  col_parallel_partitioning,
808  sparsity_pattern,
809  exchange_data,
810  communicator,
812  matrix,
815 
816  // In the end, the matrix needs to be compressed in order to be really
817  // ready.
818  last_action = Zero;
820  }
821 
822 
823 
824  void
825  SparseMatrix::reinit(const SparsityPattern &sparsity_pattern)
826  {
827  matrix.reset();
828  nonlocal_matrix_exporter.reset();
829 
830  // reinit with a (parallel) Trilinos sparsity pattern.
832  std::make_unique<Epetra_Map>(sparsity_pattern.domain_partitioner());
833  matrix = std::make_unique<Epetra_FECrsMatrix>(
834  Copy, sparsity_pattern.trilinos_sparsity_pattern(), false);
835 
836  if (sparsity_pattern.nonlocal_graph.get() != nullptr)
838  std::make_unique<Epetra_CrsMatrix>(Copy,
839  *sparsity_pattern.nonlocal_graph);
840  else
841  nonlocal_matrix.reset();
842 
843  last_action = Zero;
845  }
846 
847 
848 
849  void
850  SparseMatrix::reinit(const SparseMatrix &sparse_matrix)
851  {
852  if (this == &sparse_matrix)
853  return;
854 
856  std::make_unique<Epetra_Map>(sparse_matrix.trilinos_matrix().DomainMap());
857  matrix.reset();
858  nonlocal_matrix_exporter.reset();
859  matrix = std::make_unique<Epetra_FECrsMatrix>(
860  Copy, sparse_matrix.trilinos_sparsity_pattern(), false);
861 
862  if (sparse_matrix.nonlocal_matrix != nullptr)
863  nonlocal_matrix = std::make_unique<Epetra_CrsMatrix>(
864  Copy, sparse_matrix.nonlocal_matrix->Graph());
865  else
866  nonlocal_matrix.reset();
867 
868  last_action = Zero;
870  }
871 
872 
873 
874  template <typename number>
875  inline void
877  const IndexSet &row_parallel_partitioning,
878  const IndexSet &col_parallel_partitioning,
879  const ::SparseMatrix<number> &dealii_sparse_matrix,
880  const MPI_Comm communicator,
881  const double drop_tolerance,
882  const bool copy_values,
883  const ::SparsityPattern *use_this_sparsity)
884  {
885  if (copy_values == false)
886  {
887  // in case we do not copy values, just
888  // call the other function.
889  if (use_this_sparsity == nullptr)
890  reinit(row_parallel_partitioning,
891  col_parallel_partitioning,
892  dealii_sparse_matrix.get_sparsity_pattern(),
893  communicator,
894  false);
895  else
896  reinit(row_parallel_partitioning,
897  col_parallel_partitioning,
898  *use_this_sparsity,
899  communicator,
900  false);
901  return;
902  }
903 
904  const size_type n_rows = dealii_sparse_matrix.m();
905 
906  AssertDimension(row_parallel_partitioning.size(), n_rows);
907  AssertDimension(col_parallel_partitioning.size(), dealii_sparse_matrix.n());
908 
909  const ::SparsityPattern &sparsity_pattern =
910  (use_this_sparsity != nullptr) ?
911  *use_this_sparsity :
912  dealii_sparse_matrix.get_sparsity_pattern();
913 
914  if (matrix.get() == nullptr || m() != n_rows ||
915  n_nonzero_elements() != sparsity_pattern.n_nonzero_elements())
916  {
917  reinit(row_parallel_partitioning,
918  col_parallel_partitioning,
919  sparsity_pattern,
920  communicator,
921  false);
922  }
923 
924  // fill the values. the same as above: go through all rows of the
925  // matrix, and then all columns. since the sparsity patterns of the
926  // input matrix and the specified sparsity pattern might be different,
927  // need to go through the row for both these sparsity structures
928  // simultaneously in order to really set the correct values.
929  size_type maximum_row_length = matrix->MaxNumEntries();
930  std::vector<size_type> row_indices(maximum_row_length);
931  std::vector<TrilinosScalar> values(maximum_row_length);
932 
933  for (size_type row = 0; row < n_rows; ++row)
934  // see if the row is locally stored on this processor
935  if (row_parallel_partitioning.is_element(row) == true)
936  {
937  ::SparsityPattern::iterator select_index =
938  sparsity_pattern.begin(row);
939  typename ::SparseMatrix<number>::const_iterator it =
940  dealii_sparse_matrix.begin(row);
941  size_type col = 0;
942  if (sparsity_pattern.n_rows() == sparsity_pattern.n_cols())
943  {
944  // optimized diagonal
945  AssertDimension(it->column(), row);
946  if (std::fabs(it->value()) > drop_tolerance)
947  {
948  values[col] = it->value();
949  row_indices[col++] = it->column();
950  }
951  ++select_index;
952  ++it;
953  }
954 
955  while (it != dealii_sparse_matrix.end(row) &&
956  select_index != sparsity_pattern.end(row))
957  {
958  while (select_index->column() < it->column() &&
959  select_index != sparsity_pattern.end(row))
960  ++select_index;
961  while (it->column() < select_index->column() &&
962  it != dealii_sparse_matrix.end(row))
963  ++it;
964 
965  if (it == dealii_sparse_matrix.end(row))
966  break;
967  if (std::fabs(it->value()) > drop_tolerance)
968  {
969  values[col] = it->value();
970  row_indices[col++] = it->column();
971  }
972  ++select_index;
973  ++it;
974  }
975  set(row,
976  col,
977  reinterpret_cast<size_type *>(row_indices.data()),
978  values.data(),
979  false);
980  }
982  }
983 
984 
985 
986  template <typename number>
987  void
989  const ::SparseMatrix<number> &dealii_sparse_matrix,
990  const double drop_tolerance,
991  const bool copy_values,
992  const ::SparsityPattern *use_this_sparsity)
993  {
994  reinit(complete_index_set(dealii_sparse_matrix.m()),
995  complete_index_set(dealii_sparse_matrix.n()),
996  dealii_sparse_matrix,
997  MPI_COMM_SELF,
998  drop_tolerance,
999  copy_values,
1000  use_this_sparsity);
1001  }
1002 
1003 
1004 
1005  void
1006  SparseMatrix::reinit(const Epetra_CrsMatrix &input_matrix,
1007  const bool copy_values)
1008  {
1009  Assert(input_matrix.Filled() == true,
1010  ExcMessage("Input CrsMatrix has not called FillComplete()!"));
1011 
1012  column_space_map = std::make_unique<Epetra_Map>(input_matrix.DomainMap());
1013 
1014  const Epetra_CrsGraph *graph = &input_matrix.Graph();
1015 
1016  nonlocal_matrix.reset();
1017  nonlocal_matrix_exporter.reset();
1018  matrix.reset();
1019  matrix = std::make_unique<Epetra_FECrsMatrix>(Copy, *graph, false);
1020 
1021  matrix->FillComplete(*column_space_map, input_matrix.RangeMap(), true);
1022 
1023  if (copy_values == true)
1024  {
1025  // point to the first data entry in the two
1026  // matrices and copy the content
1027  const TrilinosScalar *in_values = input_matrix[0];
1028  TrilinosScalar *values = (*matrix)[0];
1029  const size_type my_nonzeros = input_matrix.NumMyNonzeros();
1030  std::memcpy(values, in_values, my_nonzeros * sizeof(TrilinosScalar));
1031  }
1032 
1033  last_action = Zero;
1035  }
1036 
1037 
1038 
1039  void
1041  {
1042  Epetra_CombineMode mode = last_action;
1043  if (last_action == Zero)
1044  {
1045  if ((operation == VectorOperation::add) ||
1046  (operation == VectorOperation::unknown))
1047  mode = Add;
1048  else if (operation == VectorOperation::insert)
1049  mode = Insert;
1050  else
1051  Assert(
1052  false,
1053  ExcMessage(
1054  "compress() can only be called with VectorOperation add, insert, or unknown"));
1055  }
1056  else
1057  {
1058  Assert(
1059  ((last_action == Add) && (operation != VectorOperation::insert)) ||
1060  ((last_action == Insert) && (operation != VectorOperation::add)),
1061  ExcMessage("Operation and argument to compress() do not match"));
1062  }
1063 
1064  // flush buffers
1065  int ierr;
1066  if (nonlocal_matrix.get() != nullptr && mode == Add)
1067  {
1068  // do only export in case of an add() operation, otherwise the owning
1069  // processor must have set the correct entry
1070  nonlocal_matrix->FillComplete(*column_space_map, matrix->RowMap());
1071  if (nonlocal_matrix_exporter.get() == nullptr)
1073  std::make_unique<Epetra_Export>(nonlocal_matrix->RowMap(),
1074  matrix->RowMap());
1075  ierr =
1077  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
1078  ierr = matrix->FillComplete(*column_space_map, matrix->RowMap());
1079  nonlocal_matrix->PutScalar(0);
1080  }
1081  else
1082  ierr =
1083  matrix->GlobalAssemble(*column_space_map, matrix->RowMap(), true, mode);
1084 
1085  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
1086 
1087  ierr = matrix->OptimizeStorage();
1088  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
1089 
1090  last_action = Zero;
1091 
1092  compressed = true;
1093  }
1094 
1095 
1096 
1097  void
1099  {
1100  // When we clear the matrix, reset
1101  // the pointer and generate an
1102  // empty matrix.
1104  std::make_unique<Epetra_Map>(0, 0, Utilities::Trilinos::comm_self());
1105  matrix = std::make_unique<Epetra_FECrsMatrix>(View, *column_space_map, 0);
1106  nonlocal_matrix.reset();
1107  nonlocal_matrix_exporter.reset();
1108 
1109  matrix->FillComplete();
1110 
1111  compressed = true;
1112  }
1113 
1114 
1115 
1116  void
1118  const TrilinosScalar new_diag_value)
1119  {
1120  Assert(matrix->Filled() == true, ExcMatrixNotCompressed());
1121 
1122  // Only do this on the rows owned
1123  // locally on this processor.
1124  int local_row =
1125  matrix->LRID(static_cast<TrilinosWrappers::types::int_type>(row));
1126  if (local_row >= 0)
1127  {
1129  int *col_indices;
1130  int num_entries;
1131  const int ierr =
1132  matrix->ExtractMyRowView(local_row, num_entries, values, col_indices);
1133  (void)ierr;
1134 
1135  Assert(ierr == 0, ExcTrilinosError(ierr));
1136 
1137  const std::ptrdiff_t diag_index =
1138  std::find(col_indices, col_indices + num_entries, local_row) -
1139  col_indices;
1140 
1141  for (TrilinosWrappers::types::int_type j = 0; j < num_entries; ++j)
1142  if (diag_index != j || new_diag_value == 0)
1143  values[j] = 0.;
1144 
1145  if (diag_index != num_entries)
1146  values[diag_index] = new_diag_value;
1147  }
1148  }
1149 
1150 
1151 
1152  void
1153  SparseMatrix::clear_rows(const std::vector<size_type> &rows,
1154  const TrilinosScalar new_diag_value)
1155  {
1156  for (const auto row : rows)
1157  clear_row(row, new_diag_value);
1158  }
1159 
1160 
1161 
1164  {
1165  // Extract local indices in
1166  // the matrix.
1167  int trilinos_i =
1168  matrix->LRID(static_cast<TrilinosWrappers::types::int_type>(i)),
1169  trilinos_j =
1170  matrix->LCID(static_cast<TrilinosWrappers::types::int_type>(j));
1171  TrilinosScalar value = 0.;
1172 
1173  // If the data is not on the
1174  // present processor, we throw
1175  // an exception. This is one of
1176  // the two tiny differences to
1177  // the el(i,j) call, which does
1178  // not throw any assertions.
1179  if (trilinos_i == -1)
1180  {
1181  Assert(false,
1183  i, j, local_range().first, local_range().second - 1));
1184  }
1185  else
1186  {
1187  // Check whether the matrix has
1188  // already been transformed to local
1189  // indices.
1190  Assert(matrix->Filled(), ExcMatrixNotCompressed());
1191 
1192  // Prepare pointers for extraction
1193  // of a view of the row.
1194  int nnz_present = matrix->NumMyEntries(trilinos_i);
1195  int nnz_extracted;
1196  int *col_indices;
1198 
1199  // Generate the view and make
1200  // sure that we have not generated
1201  // an error.
1202  // TODO Check that col_indices are int and not long long
1203  int ierr = matrix->ExtractMyRowView(trilinos_i,
1204  nnz_extracted,
1205  values,
1206  col_indices);
1207  (void)ierr;
1208  Assert(ierr == 0, ExcTrilinosError(ierr));
1209 
1210  Assert(nnz_present == nnz_extracted,
1211  ExcDimensionMismatch(nnz_present, nnz_extracted));
1212 
1213  // Search the index where we
1214  // look for the value, and then
1215  // finally get it.
1216  const std::ptrdiff_t local_col_index =
1217  std::find(col_indices, col_indices + nnz_present, trilinos_j) -
1218  col_indices;
1219 
1220  // This is actually the only
1221  // difference to the el(i,j)
1222  // function, which means that
1223  // we throw an exception in
1224  // this case instead of just
1225  // returning zero for an
1226  // element that is not present
1227  // in the sparsity pattern.
1228  if (local_col_index == nnz_present)
1229  {
1230  Assert(false, ExcInvalidIndex(i, j));
1231  }
1232  else
1233  value = values[local_col_index];
1234  }
1235 
1236  return value;
1237  }
1238 
1239 
1240 
1242  SparseMatrix::el(const size_type i, const size_type j) const
1243  {
1244  // Extract local indices in
1245  // the matrix.
1246  int trilinos_i =
1247  matrix->LRID(static_cast<TrilinosWrappers::types::int_type>(i)),
1248  trilinos_j =
1249  matrix->LCID(static_cast<TrilinosWrappers::types::int_type>(j));
1250  TrilinosScalar value = 0.;
1251 
1252  // If the data is not on the
1253  // present processor, we can't
1254  // continue. Just print out zero
1255  // as discussed in the
1256  // documentation of this
1257  // function. if you want error
1258  // checking, use operator().
1259  if ((trilinos_i == -1) || (trilinos_j == -1))
1260  return 0.;
1261  else
1262  {
1263  // Check whether the matrix
1264  // already is transformed to
1265  // local indices.
1266  Assert(matrix->Filled(), ExcMatrixNotCompressed());
1267 
1268  // Prepare pointers for extraction
1269  // of a view of the row.
1270  int nnz_present = matrix->NumMyEntries(trilinos_i);
1271  int nnz_extracted;
1272  int *col_indices;
1274 
1275  // Generate the view and make
1276  // sure that we have not generated
1277  // an error.
1278  int ierr = matrix->ExtractMyRowView(trilinos_i,
1279  nnz_extracted,
1280  values,
1281  col_indices);
1282  (void)ierr;
1283  Assert(ierr == 0, ExcTrilinosError(ierr));
1284 
1285  Assert(nnz_present == nnz_extracted,
1286  ExcDimensionMismatch(nnz_present, nnz_extracted));
1287 
1288  // Search the index where we
1289  // look for the value, and then
1290  // finally get it.
1291  const std::ptrdiff_t local_col_index =
1292  std::find(col_indices, col_indices + nnz_present, trilinos_j) -
1293  col_indices;
1294 
1295  // This is actually the only
1296  // difference to the () function
1297  // querying (i,j), where we throw an
1298  // exception instead of just
1299  // returning zero for an element
1300  // that is not present in the
1301  // sparsity pattern.
1302  if (local_col_index == nnz_present)
1303  value = 0;
1304  else
1305  value = values[local_col_index];
1306  }
1307 
1308  return value;
1309  }
1310 
1311 
1312 
1315  {
1316  Assert(m() == n(), ExcNotQuadratic());
1317 
1318 # ifdef DEBUG
1319  // use operator() in debug mode because
1320  // it checks if this is a valid element
1321  // (in parallel)
1322  return operator()(i, i);
1323 # else
1324  // Trilinos doesn't seem to have a
1325  // more efficient way to access the
1326  // diagonal than by just using the
1327  // standard el(i,j) function.
1328  return el(i, i);
1329 # endif
1330  }
1331 
1332 
1333 
1334  unsigned int
1336  {
1337  Assert(row < m(), ExcInternalError());
1338 
1339  // get a representation of the
1340  // present row
1341  int ncols = -1;
1342  int local_row =
1343  matrix->LRID(static_cast<TrilinosWrappers::types::int_type>(row));
1344  Assert((local_row >= 0), ExcAccessToNonlocalRow(row));
1345 
1346  // on the processor who owns this
1347  // row, we'll have a non-negative
1348  // value.
1349  if (local_row >= 0)
1350  {
1351  int ierr = matrix->NumMyRowEntries(local_row, ncols);
1352  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
1353  }
1354 
1355  return static_cast<unsigned int>(ncols);
1356  }
1357 
1358 
1359 
1360  void
1361  SparseMatrix::set(const std::vector<size_type> &row_indices,
1362  const std::vector<size_type> &col_indices,
1364  const bool elide_zero_values)
1365  {
1366  Assert(row_indices.size() == values.m(),
1367  ExcDimensionMismatch(row_indices.size(), values.m()));
1368  Assert(col_indices.size() == values.n(),
1369  ExcDimensionMismatch(col_indices.size(), values.n()));
1370 
1371  for (size_type i = 0; i < row_indices.size(); ++i)
1372  set(row_indices[i],
1373  col_indices.size(),
1374  col_indices.data(),
1375  &values(i, 0),
1376  elide_zero_values);
1377  }
1378 
1379 
1380 
1381  void
1383  const std::vector<size_type> &col_indices,
1384  const std::vector<TrilinosScalar> &values,
1385  const bool elide_zero_values)
1386  {
1387  Assert(col_indices.size() == values.size(),
1388  ExcDimensionMismatch(col_indices.size(), values.size()));
1389 
1390  set(row,
1391  col_indices.size(),
1392  col_indices.data(),
1393  values.data(),
1394  elide_zero_values);
1395  }
1396 
1397 
1398 
1399  template <>
1400  void
1401  SparseMatrix::set<TrilinosScalar>(const size_type row,
1402  const size_type n_cols,
1403  const size_type *col_indices,
1404  const TrilinosScalar *values,
1405  const bool elide_zero_values)
1406  {
1407  AssertIndexRange(row, this->m());
1408 
1409  int ierr;
1410  if (last_action == Add)
1411  {
1412  ierr =
1413  matrix->GlobalAssemble(*column_space_map, matrix->RowMap(), true);
1414 
1415  Assert(ierr == 0, ExcTrilinosError(ierr));
1416  }
1417 
1418  last_action = Insert;
1419 
1420  const TrilinosWrappers::types::int_type *col_index_ptr;
1421  const TrilinosScalar *col_value_ptr;
1422  const TrilinosWrappers::types::int_type trilinos_row = row;
1424 
1425  boost::container::small_vector<TrilinosScalar, 200> local_value_array(
1426  n_cols);
1427  boost::container::small_vector<TrilinosWrappers::types::int_type, 200>
1428  local_index_array(n_cols);
1429 
1430  // If we don't elide zeros, the pointers are already available... need to
1431  // cast to non-const pointers as that is the format taken by Trilinos (but
1432  // we will not modify const data)
1433  if (elide_zero_values == false)
1434  {
1435  col_index_ptr =
1436  reinterpret_cast<const TrilinosWrappers::types::int_type *>(
1437  col_indices);
1438  col_value_ptr = values;
1439  n_columns = n_cols;
1440  }
1441  else
1442  {
1443  // Otherwise, extract nonzero values in each row and get the
1444  // respective indices.
1445  col_index_ptr = local_index_array.data();
1446  col_value_ptr = local_value_array.data();
1447 
1448  n_columns = 0;
1449  for (size_type j = 0; j < n_cols; ++j)
1450  {
1451  const double value = values[j];
1452  AssertIsFinite(value);
1453  if (value != 0)
1454  {
1455  local_index_array[n_columns] = col_indices[j];
1456  local_value_array[n_columns] = value;
1457  n_columns++;
1458  }
1459  }
1460 
1461  AssertIndexRange(n_columns, n_cols + 1);
1462  }
1463 
1464 
1465  // If the calling matrix owns the row to which we want to insert values,
1466  // we can directly call the Epetra_CrsMatrix input function, which is much
1467  // faster than the Epetra_FECrsMatrix function. We distinguish between two
1468  // cases: the first one is when the matrix is not filled (i.e., it is
1469  // possible to add new elements to the sparsity pattern), and the second
1470  // one is when the pattern is already fixed. In the former case, we add
1471  // the possibility to insert new values, and in the second we just replace
1472  // data.
1473  if (matrix->RowMap().MyGID(
1474  static_cast<TrilinosWrappers::types::int_type>(row)) == true)
1475  {
1476  if (matrix->Filled() == false)
1477  {
1478  ierr = matrix->Epetra_CrsMatrix::InsertGlobalValues(
1479  row, static_cast<int>(n_columns), col_value_ptr, col_index_ptr);
1480 
1481  // When inserting elements, we do not want to create exceptions in
1482  // the case when inserting non-local data (since that's what we
1483  // want to do right now).
1484  if (ierr > 0)
1485  ierr = 0;
1486  }
1487  else
1488  ierr = matrix->Epetra_CrsMatrix::ReplaceGlobalValues(row,
1489  n_columns,
1490  col_value_ptr,
1491  col_index_ptr);
1492  }
1493  else
1494  {
1495  // When we're at off-processor data, we have to stick with the
1496  // standard Insert/ReplaceGlobalValues function. Nevertheless, the way
1497  // we call it is the fastest one (any other will lead to repeated
1498  // allocation and deallocation of memory in order to call the function
1499  // we already use, which is very inefficient if writing one element at
1500  // a time).
1501  compressed = false;
1502 
1503  if (matrix->Filled() == false)
1504  {
1505  ierr = matrix->InsertGlobalValues(1,
1506  &trilinos_row,
1507  n_columns,
1508  col_index_ptr,
1509  &col_value_ptr,
1510  Epetra_FECrsMatrix::ROW_MAJOR);
1511  if (ierr > 0)
1512  ierr = 0;
1513  }
1514  else
1515  ierr = matrix->ReplaceGlobalValues(1,
1516  &trilinos_row,
1517  n_columns,
1518  col_index_ptr,
1519  &col_value_ptr,
1520  Epetra_FECrsMatrix::ROW_MAJOR);
1521  // use the FECrsMatrix facilities for set even in the case when we
1522  // have explicitly set the off-processor rows because that only works
1523  // properly when adding elements, not when setting them (since we want
1524  // to only touch elements that have been set explicitly, and there is
1525  // no way on the receiving processor to identify them otherwise)
1526  }
1527 
1528  Assert(ierr <= 0, ExcAccessToNonPresentElement(row, col_index_ptr[0]));
1529  AssertThrow(ierr >= 0, ExcTrilinosError(ierr));
1530  }
1531 
1532 
1533 
1534  void
1535  SparseMatrix::add(const std::vector<size_type> &indices,
1537  const bool elide_zero_values)
1538  {
1539  Assert(indices.size() == values.m(),
1540  ExcDimensionMismatch(indices.size(), values.m()));
1541  Assert(values.m() == values.n(), ExcNotQuadratic());
1542 
1543  for (size_type i = 0; i < indices.size(); ++i)
1544  add(indices[i],
1545  indices.size(),
1546  indices.data(),
1547  &values(i, 0),
1548  elide_zero_values);
1549  }
1550 
1551 
1552 
1553  void
1554  SparseMatrix::add(const std::vector<size_type> &row_indices,
1555  const std::vector<size_type> &col_indices,
1557  const bool elide_zero_values)
1558  {
1559  Assert(row_indices.size() == values.m(),
1560  ExcDimensionMismatch(row_indices.size(), values.m()));
1561  Assert(col_indices.size() == values.n(),
1562  ExcDimensionMismatch(col_indices.size(), values.n()));
1563 
1564  for (size_type i = 0; i < row_indices.size(); ++i)
1565  add(row_indices[i],
1566  col_indices.size(),
1567  col_indices.data(),
1568  &values(i, 0),
1569  elide_zero_values);
1570  }
1571 
1572 
1573 
1574  void
1576  const std::vector<size_type> &col_indices,
1577  const std::vector<TrilinosScalar> &values,
1578  const bool elide_zero_values)
1579  {
1580  Assert(col_indices.size() == values.size(),
1581  ExcDimensionMismatch(col_indices.size(), values.size()));
1582 
1583  add(row,
1584  col_indices.size(),
1585  col_indices.data(),
1586  values.data(),
1587  elide_zero_values);
1588  }
1589 
1590 
1591 
1592  void
1594  const size_type n_cols,
1595  const size_type *col_indices,
1596  const TrilinosScalar *values,
1597  const bool elide_zero_values,
1598  const bool /*col_indices_are_sorted*/)
1599  {
1600  AssertIndexRange(row, this->m());
1601  int ierr;
1602  if (last_action == Insert)
1603  {
1604  // TODO: this could lead to a dead lock when only one processor
1605  // calls GlobalAssemble.
1606  ierr =
1607  matrix->GlobalAssemble(*column_space_map, matrix->RowMap(), false);
1608 
1609  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
1610  }
1611 
1612  last_action = Add;
1613 
1614  const TrilinosWrappers::types::int_type *col_index_ptr;
1615  const TrilinosScalar *col_value_ptr;
1616  const TrilinosWrappers::types::int_type trilinos_row = row;
1618 
1619  boost::container::small_vector<TrilinosScalar, 100> local_value_array(
1620  n_cols);
1621  boost::container::small_vector<TrilinosWrappers::types::int_type, 100>
1622  local_index_array(n_cols);
1623 
1624  // If we don't elide zeros, the pointers are already available... need to
1625  // cast to non-const pointers as that is the format taken by Trilinos (but
1626  // we will not modify const data)
1627  if (elide_zero_values == false)
1628  {
1629  col_index_ptr =
1630  reinterpret_cast<const TrilinosWrappers::types::int_type *>(
1631  col_indices);
1632  col_value_ptr = values;
1633  n_columns = n_cols;
1634 # ifdef DEBUG
1635  for (size_type j = 0; j < n_cols; ++j)
1636  AssertIsFinite(values[j]);
1637 # endif
1638  }
1639  else
1640  {
1641  // Otherwise, extract nonzero values in each row and the corresponding
1642  // index.
1643  col_index_ptr = local_index_array.data();
1644  col_value_ptr = local_value_array.data();
1645 
1646  n_columns = 0;
1647  for (size_type j = 0; j < n_cols; ++j)
1648  {
1649  const double value = values[j];
1650 
1651  AssertIsFinite(value);
1652  if (value != 0)
1653  {
1654  local_index_array[n_columns] = col_indices[j];
1655  local_value_array[n_columns] = value;
1656  ++n_columns;
1657  }
1658  }
1659 
1660  AssertIndexRange(n_columns, n_cols + 1);
1661  }
1662  // Exit early if there is nothing to do
1663  if (n_columns == 0)
1664  {
1665  return;
1666  }
1667 
1668  // If the calling processor owns the row to which we want to add values, we
1669  // can directly call the Epetra_CrsMatrix input function, which is much
1670  // faster than the Epetra_FECrsMatrix function.
1671  if (matrix->RowMap().MyGID(
1672  static_cast<TrilinosWrappers::types::int_type>(row)) == true)
1673  {
1674  ierr = matrix->Epetra_CrsMatrix::SumIntoGlobalValues(row,
1675  n_columns,
1676  col_value_ptr,
1677  col_index_ptr);
1678  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
1679  }
1680  else if (nonlocal_matrix.get() != nullptr)
1681  {
1682  compressed = false;
1683  // this is the case when we have explicitly set the off-processor rows
1684  // and want to create a separate matrix object for them (to retain
1685  // thread-safety)
1686  Assert(nonlocal_matrix->RowMap().LID(
1687  static_cast<TrilinosWrappers::types::int_type>(row)) != -1,
1688  ExcMessage("Attempted to write into off-processor matrix row "
1689  "that has not be specified as being writable upon "
1690  "initialization"));
1691  ierr = nonlocal_matrix->SumIntoGlobalValues(row,
1692  n_columns,
1693  col_value_ptr,
1694  col_index_ptr);
1695  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
1696  }
1697  else
1698  {
1699  // When we're at off-processor data, we have to stick with the
1700  // standard SumIntoGlobalValues function. Nevertheless, the way we
1701  // call it is the fastest one (any other will lead to repeated
1702  // allocation and deallocation of memory in order to call the function
1703  // we already use, which is very inefficient if writing one element at
1704  // a time).
1705  compressed = false;
1706 
1707  ierr = matrix->SumIntoGlobalValues(1,
1708  &trilinos_row,
1709  n_columns,
1710  col_index_ptr,
1711  &col_value_ptr,
1712  Epetra_FECrsMatrix::ROW_MAJOR);
1713  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
1714  }
1715 
1716 # ifdef DEBUG
1717  if (ierr > 0)
1718  {
1719  std::cout << "------------------------------------------" << std::endl;
1720  std::cout << "Got error " << ierr << " in row " << row << " of proc "
1721  << matrix->RowMap().Comm().MyPID()
1722  << " when trying to add the columns:" << std::endl;
1723  for (TrilinosWrappers::types::int_type i = 0; i < n_columns; ++i)
1724  std::cout << col_index_ptr[i] << " ";
1725  std::cout << std::endl << std::endl;
1726  std::cout << "Matrix row "
1727  << (matrix->RowMap().MyGID(
1728  static_cast<TrilinosWrappers::types::int_type>(row)) ==
1729  false ?
1730  "(nonlocal part)" :
1731  "")
1732  << " has the following indices:" << std::endl;
1733  std::vector<TrilinosWrappers::types::int_type> indices;
1734  const Epetra_CrsGraph *graph =
1735  (nonlocal_matrix.get() != nullptr &&
1736  matrix->RowMap().MyGID(
1737  static_cast<TrilinosWrappers::types::int_type>(row)) == false) ?
1738  &nonlocal_matrix->Graph() :
1739  &matrix->Graph();
1740 
1741  indices.resize(graph->NumGlobalIndices(row));
1742  int n_indices = 0;
1743  graph->ExtractGlobalRowCopy(row,
1744  indices.size(),
1745  n_indices,
1746  indices.data());
1747  AssertDimension(n_indices, indices.size());
1748 
1749  for (TrilinosWrappers::types::int_type i = 0; i < n_indices; ++i)
1750  std::cout << indices[i] << " ";
1751  std::cout << std::endl << std::endl;
1752  Assert(ierr <= 0, ExcAccessToNonPresentElement(row, col_index_ptr[0]));
1753  }
1754 # endif
1755  Assert(ierr >= 0, ExcTrilinosError(ierr));
1756  }
1757 
1758 
1759 
1760  SparseMatrix &
1762  {
1764  compress(VectorOperation::unknown); // TODO: why do we do this? Should we
1765  // not check for is_compressed?
1766 
1767  const int ierr = matrix->PutScalar(d);
1768  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
1769  if (nonlocal_matrix.get() != nullptr)
1770  nonlocal_matrix->PutScalar(d);
1771 
1772  return *this;
1773  }
1774 
1775 
1776 
1777  void
1779  {
1780  AssertDimension(rhs.m(), m());
1781  AssertDimension(rhs.n(), n());
1782  AssertDimension(rhs.local_range().first, local_range().first);
1783  AssertDimension(rhs.local_range().second, local_range().second);
1784  Assert(matrix->RowMap().SameAs(rhs.matrix->RowMap()),
1785  ExcMessage("Can only add matrices with same distribution of rows"));
1786  Assert(matrix->Filled() && rhs.matrix->Filled(),
1787  ExcMessage("Addition of matrices only allowed if matrices are "
1788  "filled, i.e., compress() has been called"));
1789 
1790  const bool same_col_map = matrix->ColMap().SameAs(rhs.matrix->ColMap());
1791 
1792  for (const auto row : locally_owned_range_indices())
1793  {
1794  const int row_local = matrix->RowMap().LID(
1795  static_cast<TrilinosWrappers::types::int_type>(row));
1796  Assert((row_local >= 0), ExcAccessToNonlocalRow(row));
1797 
1798  // First get a view to the matrix columns of both matrices. Note that
1799  // the data is in local index spaces so we need to be careful not only
1800  // to compare column indices in case they are derived from the same
1801  // map.
1802  int n_entries, rhs_n_entries;
1803  TrilinosScalar *value_ptr, *rhs_value_ptr;
1804  int *index_ptr, *rhs_index_ptr;
1805  int ierr = rhs.matrix->ExtractMyRowView(row_local,
1806  rhs_n_entries,
1807  rhs_value_ptr,
1808  rhs_index_ptr);
1809  (void)ierr;
1810  Assert(ierr == 0, ExcTrilinosError(ierr));
1811 
1812  ierr =
1813  matrix->ExtractMyRowView(row_local, n_entries, value_ptr, index_ptr);
1814  Assert(ierr == 0, ExcTrilinosError(ierr));
1815  bool expensive_checks = (n_entries != rhs_n_entries || !same_col_map);
1816  if (!expensive_checks)
1817  {
1818  // check if the column indices are the same. If yes, can simply
1819  // copy over the data.
1820  expensive_checks = std::memcmp(static_cast<void *>(index_ptr),
1821  static_cast<void *>(rhs_index_ptr),
1822  sizeof(int) * n_entries) != 0;
1823  if (!expensive_checks)
1824  for (int i = 0; i < n_entries; ++i)
1825  value_ptr[i] += rhs_value_ptr[i] * factor;
1826  }
1827  // Now to the expensive case where we need to check all column indices
1828  // against each other (transformed into global index space) and where
1829  // we need to make sure that all entries we are about to add into the
1830  // lhs matrix actually exist
1831  if (expensive_checks)
1832  {
1833  for (int i = 0; i < rhs_n_entries; ++i)
1834  {
1835  if (rhs_value_ptr[i] == 0.)
1836  continue;
1837  const TrilinosWrappers::types::int_type rhs_global_col =
1838  global_column_index(*rhs.matrix, rhs_index_ptr[i]);
1839  int local_col = matrix->ColMap().LID(rhs_global_col);
1840  int *local_index = Utilities::lower_bound(index_ptr,
1841  index_ptr + n_entries,
1842  local_col);
1843  Assert(local_index != index_ptr + n_entries &&
1844  *local_index == local_col,
1845  ExcMessage(
1846  "Adding the entries from the other matrix "
1847  "failed, because the sparsity pattern "
1848  "of that matrix includes more elements than the "
1849  "calling matrix, which is not allowed."));
1850  value_ptr[local_index - index_ptr] += factor * rhs_value_ptr[i];
1851  }
1852  }
1853  }
1854  }
1855 
1856 
1857 
1858  void
1860  {
1861  // This only flips a flag that tells
1862  // Trilinos that any vmult operation
1863  // should be done with the
1864  // transpose. However, the matrix
1865  // structure is not reset.
1866  int ierr;
1867 
1868  if (!matrix->UseTranspose())
1869  {
1870  ierr = matrix->SetUseTranspose(true);
1871  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
1872  }
1873  else
1874  {
1875  ierr = matrix->SetUseTranspose(false);
1876  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
1877  }
1878  }
1879 
1880 
1881 
1882  SparseMatrix &
1884  {
1885  const int ierr = matrix->Scale(a);
1886  Assert(ierr == 0, ExcTrilinosError(ierr));
1887  (void)ierr; // removes -Wunused-variable in optimized mode
1888 
1889  return *this;
1890  }
1891 
1892 
1893 
1894  SparseMatrix &
1896  {
1897  Assert(a != 0, ExcDivideByZero());
1898 
1899  const TrilinosScalar factor = 1. / a;
1900 
1901  const int ierr = matrix->Scale(factor);
1902  Assert(ierr == 0, ExcTrilinosError(ierr));
1903  (void)ierr; // removes -Wunused-variable in optimized mode
1904 
1905  return *this;
1906  }
1907 
1908 
1909 
1912  {
1913  Assert(matrix->Filled(), ExcMatrixNotCompressed());
1914  return matrix->NormOne();
1915  }
1916 
1917 
1918 
1921  {
1922  Assert(matrix->Filled(), ExcMatrixNotCompressed());
1923  return matrix->NormInf();
1924  }
1925 
1926 
1927 
1930  {
1931  Assert(matrix->Filled(), ExcMatrixNotCompressed());
1932  return matrix->NormFrobenius();
1933  }
1934 
1935 
1936 
1937  namespace internal
1938  {
1939  namespace SparseMatrixImplementation
1940  {
1941  template <typename VectorType>
1942  inline void
1943  check_vector_map_equality(const Epetra_CrsMatrix &,
1944  const VectorType &,
1945  const VectorType &)
1946  {}
1947 
1948  inline void
1949  check_vector_map_equality(const Epetra_CrsMatrix &m,
1951  const TrilinosWrappers::MPI::Vector &out)
1952  {
1953  Assert(in.trilinos_partitioner().SameAs(m.DomainMap()) == true,
1954  ExcMessage("The column partitioning of a matrix does not match "
1955  "the partitioning of a vector you are trying to "
1956  "multiply it with. Are you multiplying the "
1957  "matrix with a vector that has ghost elements?"));
1958  Assert(out.trilinos_partitioner().SameAs(m.RangeMap()) == true,
1959  ExcMessage("The row partitioning of a matrix does not match "
1960  "the partitioning of a vector you are trying to "
1961  "put the result of a matrix-vector product in. "
1962  "Are you trying to put the product of the "
1963  "matrix with a vector into a vector that has "
1964  "ghost elements?"));
1965  (void)m;
1966  (void)in;
1967  (void)out;
1968  }
1969  } // namespace SparseMatrixImplementation
1970  } // namespace internal
1971 
1972 
1973  template <typename VectorType>
1974  std::enable_if_t<
1975  std::is_same_v<typename VectorType::value_type, TrilinosScalar>>
1976  SparseMatrix::vmult(VectorType &dst, const VectorType &src) const
1977  {
1978  Assert(&src != &dst, ExcSourceEqualsDestination());
1979  Assert(matrix->Filled(), ExcMatrixNotCompressed());
1980  (void)src;
1981  (void)dst;
1982 
1984  src,
1985  dst);
1986  const size_type dst_local_size = internal::end(dst) - internal::begin(dst);
1987  AssertDimension(dst_local_size, matrix->RangeMap().NumMyPoints());
1988  const size_type src_local_size = internal::end(src) - internal::begin(src);
1989  AssertDimension(src_local_size, matrix->DomainMap().NumMyPoints());
1990 
1991  Epetra_MultiVector tril_dst(
1992  View, matrix->RangeMap(), internal::begin(dst), dst_local_size, 1);
1993  Epetra_MultiVector tril_src(View,
1994  matrix->DomainMap(),
1995  const_cast<TrilinosScalar *>(
1996  internal::begin(src)),
1997  src_local_size,
1998  1);
1999 
2000  const int ierr = matrix->Multiply(false, tril_src, tril_dst);
2001  Assert(ierr == 0, ExcTrilinosError(ierr));
2002  (void)ierr; // removes -Wunused-variable in optimized mode
2003  }
2004 
2005 
2006 
2007  template <typename VectorType>
2008  std::enable_if_t<
2009  !std::is_same_v<typename VectorType::value_type, TrilinosScalar>>
2010  SparseMatrix::vmult(VectorType & /*dst*/, const VectorType & /*src*/) const
2011  {
2012  AssertThrow(false, ExcNotImplemented());
2013  }
2014 
2015 
2016 
2017  template <typename VectorType>
2018  std::enable_if_t<
2019  std::is_same_v<typename VectorType::value_type, TrilinosScalar>>
2020  SparseMatrix::Tvmult(VectorType &dst, const VectorType &src) const
2021  {
2022  Assert(&src != &dst, ExcSourceEqualsDestination());
2023  Assert(matrix->Filled(), ExcMatrixNotCompressed());
2024 
2026  dst,
2027  src);
2028  const size_type dst_local_size = internal::end(dst) - internal::begin(dst);
2029  AssertDimension(dst_local_size, matrix->DomainMap().NumMyPoints());
2030  const size_type src_local_size = internal::end(src) - internal::begin(src);
2031  AssertDimension(src_local_size, matrix->RangeMap().NumMyPoints());
2032 
2033  Epetra_MultiVector tril_dst(
2034  View, matrix->DomainMap(), internal::begin(dst), dst_local_size, 1);
2035  Epetra_MultiVector tril_src(View,
2036  matrix->RangeMap(),
2037  const_cast<double *>(internal::begin(src)),
2038  src_local_size,
2039  1);
2040 
2041  const int ierr = matrix->Multiply(true, tril_src, tril_dst);
2042  Assert(ierr == 0, ExcTrilinosError(ierr));
2043  (void)ierr; // removes -Wunused-variable in optimized mode
2044  }
2045 
2046 
2047 
2048  template <typename VectorType>
2049  std::enable_if_t<
2050  !std::is_same_v<typename VectorType::value_type, TrilinosScalar>>
2051  SparseMatrix::Tvmult(VectorType & /*dst*/, const VectorType & /*src*/) const
2052  {
2053  AssertThrow(false, ExcNotImplemented());
2054  }
2055 
2056 
2057 
2058  template <typename VectorType>
2059  void
2060  SparseMatrix::vmult_add(VectorType &dst, const VectorType &src) const
2061  {
2062  Assert(&src != &dst, ExcSourceEqualsDestination());
2063 
2064  // Reinit a temporary vector with fast argument set, which does not
2065  // overwrite the content (to save time).
2066  VectorType tmp_vector;
2067  tmp_vector.reinit(dst, true);
2068  vmult(tmp_vector, src);
2069  dst += tmp_vector;
2070  }
2071 
2072 
2073 
2074  template <typename VectorType>
2075  void
2076  SparseMatrix::Tvmult_add(VectorType &dst, const VectorType &src) const
2077  {
2078  Assert(&src != &dst, ExcSourceEqualsDestination());
2079 
2080  // Reinit a temporary vector with fast argument set, which does not
2081  // overwrite the content (to save time).
2082  VectorType tmp_vector;
2083  tmp_vector.reinit(dst, true);
2084  Tvmult(tmp_vector, src);
2085  dst += tmp_vector;
2086  }
2087 
2088 
2089 
2092  {
2093  Assert(matrix->RowMap().SameAs(matrix->DomainMap()), ExcNotQuadratic());
2094 
2095  MPI::Vector temp_vector;
2096  temp_vector.reinit(v, true);
2097 
2098  vmult(temp_vector, v);
2099  return temp_vector * v;
2100  }
2101 
2102 
2103 
2106  const MPI::Vector &v) const
2107  {
2108  Assert(matrix->RowMap().SameAs(matrix->DomainMap()), ExcNotQuadratic());
2109 
2110  MPI::Vector temp_vector;
2111  temp_vector.reinit(v, true);
2112 
2113  vmult(temp_vector, v);
2114  return u * temp_vector;
2115  }
2116 
2117 
2118 
2121  const MPI::Vector &x,
2122  const MPI::Vector &b) const
2123  {
2124  vmult(dst, x);
2125  dst -= b;
2126  dst *= -1.;
2127 
2128  return dst.l2_norm();
2129  }
2130 
2131 
2132 
2133  namespace internals
2134  {
2136 
2137  void
2138  perform_mmult(const SparseMatrix &inputleft,
2139  const SparseMatrix &inputright,
2140  SparseMatrix &result,
2141  const MPI::Vector &V,
2142  const bool transpose_left)
2143  {
2144  const bool use_vector = (V.size() == inputright.m() ? true : false);
2145  if (transpose_left == false)
2146  {
2147  Assert(inputleft.n() == inputright.m(),
2148  ExcDimensionMismatch(inputleft.n(), inputright.m()));
2149  Assert(inputleft.trilinos_matrix().DomainMap().SameAs(
2150  inputright.trilinos_matrix().RangeMap()),
2151  ExcMessage("Parallel partitioning of A and B does not fit."));
2152  }
2153  else
2154  {
2155  Assert(inputleft.m() == inputright.m(),
2156  ExcDimensionMismatch(inputleft.m(), inputright.m()));
2157  Assert(inputleft.trilinos_matrix().RangeMap().SameAs(
2158  inputright.trilinos_matrix().RangeMap()),
2159  ExcMessage("Parallel partitioning of A and B does not fit."));
2160  }
2161 
2162  result.clear();
2163 
2164  // create a suitable operator B: in case
2165  // we do not use a vector, all we need to
2166  // do is to set the pointer. Otherwise,
2167  // we insert the data from B, but
2168  // multiply each row with the respective
2169  // vector element.
2170  Teuchos::RCP<Epetra_CrsMatrix> mod_B;
2171  if (use_vector == false)
2172  {
2173  mod_B = Teuchos::rcp(const_cast<Epetra_CrsMatrix *>(
2174  &inputright.trilinos_matrix()),
2175  false);
2176  }
2177  else
2178  {
2179  mod_B = Teuchos::rcp(
2180  new Epetra_CrsMatrix(Copy, inputright.trilinos_sparsity_pattern()),
2181  true);
2182  mod_B->FillComplete(inputright.trilinos_matrix().DomainMap(),
2183  inputright.trilinos_matrix().RangeMap());
2184  Assert(inputright.local_range() == V.local_range(),
2185  ExcMessage("Parallel distribution of matrix B and vector V "
2186  "does not match."));
2187 
2188  const int local_N = inputright.local_size();
2189  for (int i = 0; i < local_N; ++i)
2190  {
2191  int N_entries = -1;
2192  double *new_data, *B_data;
2193  mod_B->ExtractMyRowView(i, N_entries, new_data);
2194  inputright.trilinos_matrix().ExtractMyRowView(i,
2195  N_entries,
2196  B_data);
2197  double value = V.trilinos_vector()[0][i];
2198  for (TrilinosWrappers::types::int_type j = 0; j < N_entries; ++j)
2199  new_data[j] = value * B_data[j];
2200  }
2201  }
2202 
2203 
2204  SparseMatrix tmp_result(transpose_left ?
2205  inputleft.locally_owned_domain_indices() :
2206  inputleft.locally_owned_range_indices(),
2207  inputright.locally_owned_domain_indices(),
2208  inputleft.get_mpi_communicator());
2209 
2210 # ifdef DEAL_II_TRILINOS_WITH_EPETRAEXT
2211  EpetraExt::MatrixMatrix::Multiply(inputleft.trilinos_matrix(),
2212  transpose_left,
2213  *mod_B,
2214  false,
2215  const_cast<Epetra_CrsMatrix &>(
2216  tmp_result.trilinos_matrix()));
2217 # else
2218  Assert(false,
2219  ExcMessage("This function requires that the Trilinos "
2220  "installation found while running the deal.II "
2221  "CMake scripts contains the optional Trilinos "
2222  "package 'EpetraExt'. However, this optional "
2223  "part of Trilinos was not found."));
2224 # endif
2225  result.reinit(tmp_result.trilinos_matrix());
2226  }
2227  } // namespace internals
2228 
2229 
2230  void
2232  const SparseMatrix &B,
2233  const MPI::Vector &V) const
2234  {
2235  internals::perform_mmult(*this, B, C, V, false);
2236  }
2237 
2238 
2239 
2240  void
2242  const SparseMatrix &B,
2243  const MPI::Vector &V) const
2244  {
2245  internals::perform_mmult(*this, B, C, V, true);
2246  }
2247 
2248 
2249 
2250  void
2252  {
2253  Assert(false, ExcNotImplemented());
2254  }
2255 
2256 
2257 
2258  // As of now, no particularly neat
2259  // output is generated in case of
2260  // multiple processors.
2261  void
2262  SparseMatrix::print(std::ostream &out,
2263  const bool print_detailed_trilinos_information) const
2264  {
2265  if (print_detailed_trilinos_information == true)
2266  out << *matrix;
2267  else
2268  {
2269  double *values;
2270  int *indices;
2271  int num_entries;
2272 
2273  for (int i = 0; i < matrix->NumMyRows(); ++i)
2274  {
2275  const int ierr =
2276  matrix->ExtractMyRowView(i, num_entries, values, indices);
2277  (void)ierr;
2278  Assert(ierr == 0, ExcTrilinosError(ierr));
2279 
2280  for (TrilinosWrappers::types::int_type j = 0; j < num_entries; ++j)
2281  out << "(" << TrilinosWrappers::global_row_index(*matrix, i)
2282  << ","
2284  << ") " << values[j] << std::endl;
2285  }
2286  }
2287 
2288  AssertThrow(out.fail() == false, ExcIO());
2289  }
2290 
2291 
2292 
2295  {
2296  size_type static_memory =
2297  sizeof(*this) + sizeof(*matrix) + sizeof(*matrix->Graph().DataPtr());
2298  return (
2300  matrix->NumMyNonzeros() +
2301  sizeof(int) * local_size() + static_memory);
2302  }
2303 
2304 
2305 
2306  MPI_Comm
2308  {
2309  const Epetra_MpiComm *mpi_comm =
2310  dynamic_cast<const Epetra_MpiComm *>(&matrix->RangeMap().Comm());
2311  Assert(mpi_comm != nullptr, ExcInternalError());
2312  return mpi_comm->Comm();
2313  }
2314 } // namespace TrilinosWrappers
2315 
2316 
2317 namespace TrilinosWrappers
2318 {
2319  namespace internal
2320  {
2321  namespace LinearOperatorImplementation
2322  {
2324  : use_transpose(false)
2325  , communicator(MPI_COMM_SELF)
2326  , domain_map(IndexSet().make_trilinos_map(communicator.Comm()))
2327  , range_map(IndexSet().make_trilinos_map(communicator.Comm()))
2328  {
2329  vmult = [](Range &, const Domain &) {
2330  Assert(false,
2331  ExcMessage("Uninitialized TrilinosPayload::vmult called "
2332  "(Default constructor)"));
2333  };
2334 
2335  Tvmult = [](Domain &, const Range &) {
2336  Assert(false,
2337  ExcMessage("Uninitialized TrilinosPayload::Tvmult called "
2338  "(Default constructor)"));
2339  };
2340 
2341  inv_vmult = [](Domain &, const Range &) {
2342  Assert(false,
2343  ExcMessage("Uninitialized TrilinosPayload::inv_vmult called "
2344  "(Default constructor)"));
2345  };
2346 
2347  inv_Tvmult = [](Range &, const Domain &) {
2348  Assert(false,
2349  ExcMessage("Uninitialized TrilinosPayload::inv_Tvmult called "
2350  "(Default constructor)"));
2351  };
2352  }
2353 
2354 
2355 
2357  const TrilinosWrappers::SparseMatrix &matrix_exemplar,
2359  : TrilinosPayload(const_cast<Epetra_CrsMatrix &>(
2360  matrix.trilinos_matrix()),
2361  /*op_supports_inverse_operations = */ false,
2362  matrix_exemplar.trilinos_matrix().UseTranspose(),
2363  matrix_exemplar.get_mpi_communicator(),
2364  matrix_exemplar.locally_owned_domain_indices(),
2365  matrix_exemplar.locally_owned_range_indices())
2366  {}
2367 
2368 
2369 
2371  const TrilinosPayload &payload_exemplar,
2373 
2374  : TrilinosPayload(const_cast<Epetra_CrsMatrix &>(
2375  matrix.trilinos_matrix()),
2376  /*op_supports_inverse_operations = */ false,
2377  payload_exemplar.UseTranspose(),
2378  payload_exemplar.get_mpi_communicator(),
2379  payload_exemplar.locally_owned_domain_indices(),
2380  payload_exemplar.locally_owned_range_indices())
2381  {}
2382 
2383 
2384 
2386  const TrilinosWrappers::SparseMatrix &matrix_exemplar,
2387  const TrilinosWrappers::PreconditionBase &preconditioner)
2388  : TrilinosPayload(preconditioner.trilinos_operator(),
2389  /*op_supports_inverse_operations = */ true,
2390  matrix_exemplar.trilinos_matrix().UseTranspose(),
2391  matrix_exemplar.get_mpi_communicator(),
2392  matrix_exemplar.locally_owned_domain_indices(),
2393  matrix_exemplar.locally_owned_range_indices())
2394  {}
2395 
2396 
2397 
2399  const TrilinosWrappers::PreconditionBase &preconditioner_exemplar,
2400  const TrilinosWrappers::PreconditionBase &preconditioner)
2401  : TrilinosPayload(
2402  preconditioner.trilinos_operator(),
2403  /*op_supports_inverse_operations = */ true,
2404  preconditioner_exemplar.trilinos_operator().UseTranspose(),
2405  preconditioner_exemplar.get_mpi_communicator(),
2406  preconditioner_exemplar.locally_owned_domain_indices(),
2407  preconditioner_exemplar.locally_owned_range_indices())
2408  {}
2409 
2410 
2411 
2413  const TrilinosPayload &payload_exemplar,
2414  const TrilinosWrappers::PreconditionBase &preconditioner)
2415  : TrilinosPayload(preconditioner.trilinos_operator(),
2416  /*op_supports_inverse_operations = */ true,
2417  payload_exemplar.UseTranspose(),
2418  payload_exemplar.get_mpi_communicator(),
2419  payload_exemplar.locally_owned_domain_indices(),
2420  payload_exemplar.locally_owned_range_indices())
2421  {}
2422 
2423 
2424 
2426  : vmult(payload.vmult)
2427  , Tvmult(payload.Tvmult)
2428  , inv_vmult(payload.inv_vmult)
2429  , inv_Tvmult(payload.inv_Tvmult)
2430  , use_transpose(payload.use_transpose)
2431  , communicator(payload.communicator)
2432  , domain_map(payload.domain_map)
2433  , range_map(payload.range_map)
2434  {}
2435 
2436 
2437 
2438  // Composite copy constructor
2439  // This is required for PackagedOperations
2441  const TrilinosPayload &second_op)
2442  : use_transpose(false)
2443  , // The combination of operators provides the exact
2444  // definition of the operation
2445  communicator(first_op.communicator)
2446  , domain_map(second_op.domain_map)
2447  , range_map(first_op.range_map)
2448  {}
2449 
2450 
2451 
2454  {
2455  TrilinosPayload return_op(*this);
2456 
2457  return_op.vmult = [](Range &tril_dst, const Range &tril_src) {
2458  tril_dst = tril_src;
2459  };
2460 
2461  return_op.Tvmult = [](Range &tril_dst, const Range &tril_src) {
2462  tril_dst = tril_src;
2463  };
2464 
2465  return_op.inv_vmult = [](Range &tril_dst, const Range &tril_src) {
2466  tril_dst = tril_src;
2467  };
2468 
2469  return_op.inv_Tvmult = [](Range &tril_dst, const Range &tril_src) {
2470  tril_dst = tril_src;
2471  };
2472 
2473  return return_op;
2474  }
2475 
2476 
2477 
2480  {
2481  TrilinosPayload return_op(*this);
2482 
2483  return_op.vmult = [](Range &tril_dst, const Domain &) {
2484  const int ierr = tril_dst.PutScalar(0.0);
2485 
2486  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
2487  };
2488 
2489  return_op.Tvmult = [](Domain &tril_dst, const Range &) {
2490  const int ierr = tril_dst.PutScalar(0.0);
2491 
2492  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
2493  };
2494 
2495  return_op.inv_vmult = [](Domain &tril_dst, const Range &) {
2496  AssertThrow(false,
2497  ExcMessage("Cannot compute inverse of null operator"));
2498 
2499  const int ierr = tril_dst.PutScalar(0.0);
2500 
2501  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
2502  };
2503 
2504  return_op.inv_Tvmult = [](Range &tril_dst, const Domain &) {
2505  AssertThrow(false,
2506  ExcMessage("Cannot compute inverse of null operator"));
2507 
2508  const int ierr = tril_dst.PutScalar(0.0);
2509 
2510  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
2511  };
2512 
2513  return return_op;
2514  }
2515 
2516 
2517 
2520  {
2521  TrilinosPayload return_op(*this);
2522  return_op.transpose();
2523  return return_op;
2524  }
2525 
2526 
2527 
2528  IndexSet
2530  {
2531  return IndexSet(domain_map);
2532  }
2533 
2534 
2535 
2536  IndexSet
2538  {
2539  return IndexSet(range_map);
2540  }
2541 
2542 
2543 
2544  MPI_Comm
2546  {
2547  return communicator.Comm();
2548  }
2549 
2550 
2551 
2552  void
2554  {
2556  }
2557 
2558 
2559 
2560  bool
2562  {
2563  return use_transpose;
2564  }
2565 
2566 
2567 
2568  int
2570  {
2571  if (use_transpose != UseTranspose)
2572  {
2577  }
2578  return 0;
2579  }
2580 
2581 
2582 
2583  int
2585  {
2586  // The transposedness of the operations is taken care of
2587  // when we hit the transpose flag.
2588  vmult(Y, X);
2589  return 0;
2590  }
2591 
2592 
2593 
2594  int
2596  {
2597  // The transposedness of the operations is taken care of
2598  // when we hit the transpose flag.
2599  inv_vmult(X, Y);
2600  return 0;
2601  }
2602 
2603 
2604 
2605  const char *
2607  {
2608  return "TrilinosPayload";
2609  }
2610 
2611 
2612 
2613  const Epetra_Comm &
2615  {
2616  return communicator;
2617  }
2618 
2619 
2620 
2621  const Epetra_Map &
2623  {
2624  return domain_map;
2625  }
2626 
2627 
2628 
2629  const Epetra_Map &
2631  {
2632  return range_map;
2633  }
2634 
2635 
2636 
2637  bool
2639  {
2640  return false;
2641  }
2642 
2643 
2644 
2645  double
2647  {
2648  AssertThrow(false, ExcNotImplemented());
2649  return 0.0;
2650  }
2651 
2652 
2653 
2655  operator+(const TrilinosPayload &first_op,
2656  const TrilinosPayload &second_op)
2657  {
2658  using Domain = typename TrilinosPayload::Domain;
2659  using Range = typename TrilinosPayload::Range;
2660  using Intermediate = typename TrilinosPayload::VectorType;
2661  using GVMVectorType = TrilinosWrappers::MPI::Vector;
2662 
2663  Assert(first_op.locally_owned_domain_indices() ==
2664  second_op.locally_owned_domain_indices(),
2665  ExcMessage(
2666  "Operators are set to work on incompatible IndexSets."));
2667  Assert(first_op.locally_owned_range_indices() ==
2668  second_op.locally_owned_range_indices(),
2669  ExcMessage(
2670  "Operators are set to work on incompatible IndexSets."));
2671 
2672  TrilinosPayload return_op(first_op, second_op);
2673 
2674  // Capture by copy so the payloads are always valid
2675  return_op.vmult = [first_op, second_op](Range &tril_dst,
2676  const Domain &tril_src) {
2677  // Duplicated from LinearOperator::operator*
2678  // TODO: Template the constructor on GrowingVectorMemory vector type?
2679  GrowingVectorMemory<GVMVectorType> vector_memory;
2680  VectorMemory<GVMVectorType>::Pointer i(vector_memory);
2681 
2682  // Initialize intermediate vector
2683  const Epetra_Map &first_op_init_map = first_op.OperatorDomainMap();
2684  i->reinit(IndexSet(first_op_init_map),
2685  first_op.get_mpi_communicator(),
2686  /*bool omit_zeroing_entries =*/true);
2687 
2688  // Duplicated from TrilinosWrappers::SparseMatrix::vmult
2689  const size_type i_local_size = i->end() - i->begin();
2690  AssertDimension(i_local_size, first_op_init_map.NumMyPoints());
2691  const Epetra_Map &second_op_init_map = second_op.OperatorDomainMap();
2692  AssertDimension(i_local_size, second_op_init_map.NumMyPoints());
2693  (void)second_op_init_map;
2694  Intermediate tril_int(View,
2695  first_op_init_map,
2696  const_cast<TrilinosScalar *>(i->begin()),
2697  i_local_size,
2698  1);
2699 
2700  // These operators may themselves be transposed or not, so we let them
2701  // decide what the intended outcome is
2702  second_op.Apply(tril_src, tril_int);
2703  first_op.Apply(tril_src, tril_dst);
2704  const int ierr = tril_dst.Update(1.0, tril_int, 1.0);
2705  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
2706  };
2707 
2708  return_op.Tvmult = [first_op, second_op](Domain &tril_dst,
2709  const Range &tril_src) {
2710  // Duplicated from LinearOperator::operator*
2711  // TODO: Template the constructor on GrowingVectorMemory vector type?
2712  GrowingVectorMemory<GVMVectorType> vector_memory;
2713  VectorMemory<GVMVectorType>::Pointer i(vector_memory);
2714 
2715  // These operators may themselves be transposed or not, so we let them
2716  // decide what the intended outcome is
2717  // We must first transpose the operators to get the right IndexSets
2718  // for the input, intermediate and result vectors
2719  const_cast<TrilinosPayload &>(first_op).transpose();
2720  const_cast<TrilinosPayload &>(second_op).transpose();
2721 
2722  // Initialize intermediate vector
2723  const Epetra_Map &first_op_init_map = first_op.OperatorRangeMap();
2724  i->reinit(IndexSet(first_op_init_map),
2725  first_op.get_mpi_communicator(),
2726  /*bool omit_zeroing_entries =*/true);
2727 
2728  // Duplicated from TrilinosWrappers::SparseMatrix::vmult
2729  const size_type i_local_size = i->end() - i->begin();
2730  AssertDimension(i_local_size, first_op_init_map.NumMyPoints());
2731  const Epetra_Map &second_op_init_map = second_op.OperatorRangeMap();
2732  AssertDimension(i_local_size, second_op_init_map.NumMyPoints());
2733  (void)second_op_init_map;
2734  Intermediate tril_int(View,
2735  first_op_init_map,
2736  const_cast<TrilinosScalar *>(i->begin()),
2737  i_local_size,
2738  1);
2739 
2740  // These operators may themselves be transposed or not, so we let them
2741  // decide what the intended outcome is
2742  second_op.Apply(tril_src, tril_int);
2743  first_op.Apply(tril_src, tril_dst);
2744  const int ierr = tril_dst.Update(1.0, tril_int, 1.0);
2745  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
2746 
2747  // Reset transpose flag
2748  const_cast<TrilinosPayload &>(first_op).transpose();
2749  const_cast<TrilinosPayload &>(second_op).transpose();
2750  };
2751 
2752  return_op.inv_vmult = [first_op, second_op](Domain &tril_dst,
2753  const Range &tril_src) {
2754  // Duplicated from LinearOperator::operator*
2755  // TODO: Template the constructor on GrowingVectorMemory vector type?
2756  GrowingVectorMemory<GVMVectorType> vector_memory;
2757  VectorMemory<GVMVectorType>::Pointer i(vector_memory);
2758 
2759  // Initialize intermediate vector
2760  const Epetra_Map &first_op_init_map = first_op.OperatorRangeMap();
2761  i->reinit(IndexSet(first_op_init_map),
2762  first_op.get_mpi_communicator(),
2763  /*bool omit_zeroing_entries =*/true);
2764 
2765  // Duplicated from TrilinosWrappers::SparseMatrix::vmult
2766  const size_type i_local_size = i->end() - i->begin();
2767  AssertDimension(i_local_size, first_op_init_map.NumMyPoints());
2768  const Epetra_Map &second_op_init_map = second_op.OperatorRangeMap();
2769  AssertDimension(i_local_size, second_op_init_map.NumMyPoints());
2770  (void)second_op_init_map;
2771  Intermediate tril_int(View,
2772  first_op_init_map,
2773  const_cast<TrilinosScalar *>(i->begin()),
2774  i_local_size,
2775  1);
2776 
2777  // These operators may themselves be transposed or not, so we let them
2778  // decide what the intended outcome is
2779  second_op.ApplyInverse(tril_src, tril_int);
2780  first_op.ApplyInverse(tril_src, tril_dst);
2781  const int ierr = tril_dst.Update(1.0, tril_int, 1.0);
2782  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
2783  };
2784 
2785  return_op.inv_Tvmult = [first_op, second_op](Range &tril_dst,
2786  const Domain &tril_src) {
2787  // Duplicated from LinearOperator::operator*
2788  // TODO: Template the constructor on GrowingVectorMemory vector type?
2789  GrowingVectorMemory<GVMVectorType> vector_memory;
2790  VectorMemory<GVMVectorType>::Pointer i(vector_memory);
2791 
2792  // These operators may themselves be transposed or not, so we let them
2793  // decide what the intended outcome is
2794  // We must first transpose the operators to get the right IndexSets
2795  // for the input, intermediate and result vectors
2796  const_cast<TrilinosPayload &>(first_op).transpose();
2797  const_cast<TrilinosPayload &>(second_op).transpose();
2798 
2799  // Initialize intermediate vector
2800  const Epetra_Map &first_op_init_map = first_op.OperatorDomainMap();
2801  i->reinit(IndexSet(first_op_init_map),
2802  first_op.get_mpi_communicator(),
2803  /*bool omit_zeroing_entries =*/true);
2804 
2805  // Duplicated from TrilinosWrappers::SparseMatrix::vmult
2806  const size_type i_local_size = i->end() - i->begin();
2807  AssertDimension(i_local_size, first_op_init_map.NumMyPoints());
2808  const Epetra_Map &second_op_init_map = second_op.OperatorDomainMap();
2809  AssertDimension(i_local_size, second_op_init_map.NumMyPoints());
2810  (void)second_op_init_map;
2811  Intermediate tril_int(View,
2812  first_op_init_map,
2813  const_cast<TrilinosScalar *>(i->begin()),
2814  i_local_size,
2815  1);
2816 
2817  // These operators may themselves be transposed or not, so we let them
2818  // decide what the intended outcome is
2819  second_op.ApplyInverse(tril_src, tril_int);
2820  first_op.ApplyInverse(tril_src, tril_dst);
2821  const int ierr = tril_dst.Update(1.0, tril_int, 1.0);
2822  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
2823 
2824  // Reset transpose flag
2825  const_cast<TrilinosPayload &>(first_op).transpose();
2826  const_cast<TrilinosPayload &>(second_op).transpose();
2827  };
2828 
2829  return return_op;
2830  }
2831 
2832 
2833 
2834  TrilinosPayload
2835  operator*(const TrilinosPayload &first_op,
2836  const TrilinosPayload &second_op)
2837  {
2838  using Domain = typename TrilinosPayload::Domain;
2839  using Range = typename TrilinosPayload::Range;
2840  using Intermediate = typename TrilinosPayload::VectorType;
2841  using GVMVectorType = TrilinosWrappers::MPI::Vector;
2842 
2844  second_op.locally_owned_range_indices(),
2845  ExcMessage(
2846  "Operators are set to work on incompatible IndexSets."));
2847 
2848  TrilinosPayload return_op(first_op, second_op);
2849 
2850  // Capture by copy so the payloads are always valid
2851  return_op.vmult = [first_op, second_op](Range &tril_dst,
2852  const Domain &tril_src) {
2853  // Duplicated from LinearOperator::operator*
2854  // TODO: Template the constructor on GrowingVectorMemory vector type?
2855  GrowingVectorMemory<GVMVectorType> vector_memory;
2856  VectorMemory<GVMVectorType>::Pointer i(vector_memory);
2857 
2858  // Initialize intermediate vector
2859  const Epetra_Map &first_op_init_map = first_op.OperatorDomainMap();
2860  i->reinit(IndexSet(first_op_init_map),
2861  first_op.get_mpi_communicator(),
2862  /*bool omit_zeroing_entries =*/true);
2863 
2864  // Duplicated from TrilinosWrappers::SparseMatrix::vmult
2865  const size_type i_local_size = i->end() - i->begin();
2866  AssertDimension(i_local_size, first_op_init_map.NumMyPoints());
2867  const Epetra_Map &second_op_init_map = second_op.OperatorRangeMap();
2868  AssertDimension(i_local_size, second_op_init_map.NumMyPoints());
2869  (void)second_op_init_map;
2870  Intermediate tril_int(View,
2871  first_op_init_map,
2872  const_cast<TrilinosScalar *>(i->begin()),
2873  i_local_size,
2874  1);
2875 
2876  // These operators may themselves be transposed or not, so we let them
2877  // decide what the intended outcome is
2878  second_op.Apply(tril_src, tril_int);
2879  first_op.Apply(tril_int, tril_dst);
2880  };
2881 
2882  return_op.Tvmult = [first_op, second_op](Domain &tril_dst,
2883  const Range &tril_src) {
2884  // Duplicated from LinearOperator::operator*
2885  // TODO: Template the constructor on GrowingVectorMemory vector type?
2886  GrowingVectorMemory<GVMVectorType> vector_memory;
2887  VectorMemory<GVMVectorType>::Pointer i(vector_memory);
2888 
2889  // These operators may themselves be transposed or not, so we let them
2890  // decide what the intended outcome is
2891  // We must first transpose the operators to get the right IndexSets
2892  // for the input, intermediate and result vectors
2893  const_cast<TrilinosPayload &>(first_op).transpose();
2894  const_cast<TrilinosPayload &>(second_op).transpose();
2895 
2896  // Initialize intermediate vector
2897  const Epetra_Map &first_op_init_map = first_op.OperatorRangeMap();
2898  i->reinit(IndexSet(first_op_init_map),
2899  first_op.get_mpi_communicator(),
2900  /*bool omit_zeroing_entries =*/true);
2901 
2902  // Duplicated from TrilinosWrappers::SparseMatrix::vmult
2903  const size_type i_local_size = i->end() - i->begin();
2904  AssertDimension(i_local_size, first_op_init_map.NumMyPoints());
2905  const Epetra_Map &second_op_init_map = second_op.OperatorDomainMap();
2906  AssertDimension(i_local_size, second_op_init_map.NumMyPoints());
2907  (void)second_op_init_map;
2908  Intermediate tril_int(View,
2909  first_op_init_map,
2910  const_cast<TrilinosScalar *>(i->begin()),
2911  i_local_size,
2912  1);
2913 
2914  // Apply the operators in the reverse order to vmult
2915  first_op.Apply(tril_src, tril_int);
2916  second_op.Apply(tril_int, tril_dst);
2917 
2918  // Reset transpose flag
2919  const_cast<TrilinosPayload &>(first_op).transpose();
2920  const_cast<TrilinosPayload &>(second_op).transpose();
2921  };
2922 
2923  return_op.inv_vmult = [first_op, second_op](Domain &tril_dst,
2924  const Range &tril_src) {
2925  // Duplicated from LinearOperator::operator*
2926  // TODO: Template the constructor on GrowingVectorMemory vector type?
2927  GrowingVectorMemory<GVMVectorType> vector_memory;
2928  VectorMemory<GVMVectorType>::Pointer i(vector_memory);
2929 
2930  // Initialize intermediate vector
2931  const Epetra_Map &first_op_init_map = first_op.OperatorRangeMap();
2932  i->reinit(IndexSet(first_op_init_map),
2933  first_op.get_mpi_communicator(),
2934  /*bool omit_zeroing_entries =*/true);
2935 
2936  // Duplicated from TrilinosWrappers::SparseMatrix::vmult
2937  const size_type i_local_size = i->end() - i->begin();
2938  AssertDimension(i_local_size, first_op_init_map.NumMyPoints());
2939  const Epetra_Map &second_op_init_map = second_op.OperatorDomainMap();
2940  AssertDimension(i_local_size, second_op_init_map.NumMyPoints());
2941  (void)second_op_init_map;
2942  Intermediate tril_int(View,
2943  first_op_init_map,
2944  const_cast<TrilinosScalar *>(i->begin()),
2945  i_local_size,
2946  1);
2947 
2948  // Apply the operators in the reverse order to vmult
2949  // and the same order as Tvmult
2950  first_op.ApplyInverse(tril_src, tril_int);
2951  second_op.ApplyInverse(tril_int, tril_dst);
2952  };
2953 
2954  return_op.inv_Tvmult = [first_op, second_op](Range &tril_dst,
2955  const Domain &tril_src) {
2956  // Duplicated from LinearOperator::operator*
2957  // TODO: Template the constructor on GrowingVectorMemory vector type?
2958  GrowingVectorMemory<GVMVectorType> vector_memory;
2959  VectorMemory<GVMVectorType>::Pointer i(vector_memory);
2960 
2961  // These operators may themselves be transposed or not, so we let them
2962  // decide what the intended outcome is
2963  // We must first transpose the operators to get the right IndexSets
2964  // for the input, intermediate and result vectors
2965  const_cast<TrilinosPayload &>(first_op).transpose();
2966  const_cast<TrilinosPayload &>(second_op).transpose();
2967 
2968  // Initialize intermediate vector
2969  const Epetra_Map &first_op_init_map = first_op.OperatorDomainMap();
2970  i->reinit(IndexSet(first_op_init_map),
2971  first_op.get_mpi_communicator(),
2972  /*bool omit_zeroing_entries =*/true);
2973 
2974  // Duplicated from TrilinosWrappers::SparseMatrix::vmult
2975  const size_type i_local_size = i->end() - i->begin();
2976  AssertDimension(i_local_size, first_op_init_map.NumMyPoints());
2977  const Epetra_Map &second_op_init_map = second_op.OperatorRangeMap();
2978  AssertDimension(i_local_size, second_op_init_map.NumMyPoints());
2979  (void)second_op_init_map;
2980  Intermediate tril_int(View,
2981  first_op_init_map,
2982  const_cast<TrilinosScalar *>(i->begin()),
2983  i_local_size,
2984  1);
2985 
2986  // These operators may themselves be transposed or not, so we let them
2987  // decide what the intended outcome is
2988  // Apply the operators in the reverse order to Tvmult
2989  // and the same order as vmult
2990  second_op.ApplyInverse(tril_src, tril_int);
2991  first_op.ApplyInverse(tril_int, tril_dst);
2992 
2993  // Reset transpose flag
2994  const_cast<TrilinosPayload &>(first_op).transpose();
2995  const_cast<TrilinosPayload &>(second_op).transpose();
2996  };
2997 
2998  return return_op;
2999  }
3000 
3001  } // namespace LinearOperatorImplementation
3002  } /* namespace internal */
3003 } /* namespace TrilinosWrappers */
3004 
3005 
3006 
3007 // explicit instantiations
3008 # include "trilinos_sparse_matrix.inst"
3009 
3010 # ifndef DOXYGEN
3011 // TODO: put these instantiations into generic file
3012 namespace TrilinosWrappers
3013 {
3014  template void
3015  SparseMatrix::reinit(const ::SparsityPattern &);
3016 
3017  template void
3019 
3020  template void
3022  const IndexSet &,
3023  const ::SparsityPattern &,
3024  const MPI_Comm,
3025  const bool);
3026 
3027  template void
3029  const IndexSet &,
3030  const DynamicSparsityPattern &,
3031  const MPI_Comm,
3032  const bool);
3033 
3034  template void
3035  SparseMatrix::vmult(MPI::Vector &, const MPI::Vector &) const;
3036 
3037  template void
3039  const ::Vector<double> &) const;
3040 
3041  template void
3044  const ::LinearAlgebra::distributed::Vector<double> &) const;
3045 
3046 # ifdef DEAL_II_TRILINOS_WITH_TPETRA
3047 # if defined(HAVE_TPETRA_INST_DOUBLE)
3048  template void
3051  const ::LinearAlgebra::TpetraWrappers::Vector<double> &) const;
3052 # endif
3053 
3054 # if defined(HAVE_TPETRA_INST_FLOAT)
3055  template void
3058  const ::LinearAlgebra::TpetraWrappers::Vector<float> &) const;
3059 # endif
3060 # endif
3061 
3062  template void
3066 
3067  template void
3068  SparseMatrix::Tvmult(MPI::Vector &, const MPI::Vector &) const;
3069 
3070  template void
3072  const ::Vector<double> &) const;
3073 
3074  template void
3077  const ::LinearAlgebra::distributed::Vector<double> &) const;
3078 
3079 # ifdef DEAL_II_TRILINOS_WITH_TPETRA
3080 # if defined(HAVE_TPETRA_INST_DOUBLE)
3081  template void
3084  const ::LinearAlgebra::TpetraWrappers::Vector<double> &) const;
3085 # endif
3086 
3087 # if defined(HAVE_TPETRA_INST_FLOAT)
3088  template void
3091  const ::LinearAlgebra::TpetraWrappers::Vector<float> &) const;
3092 # endif
3093 # endif
3094 
3095  template void
3099 
3100  template void
3102 
3103  template void
3105  const ::Vector<double> &) const;
3106 
3107  template void
3110  const ::LinearAlgebra::distributed::Vector<double> &) const;
3111 
3112 # ifdef DEAL_II_TRILINOS_WITH_TPETRA
3113 # if defined(HAVE_TPETRA_INST_DOUBLE)
3114  template void
3117  const ::LinearAlgebra::TpetraWrappers::Vector<double> &) const;
3118 # endif
3119 
3120 # if defined(HAVE_TPETRA_INST_FLOAT)
3121  template void
3124  const ::LinearAlgebra::TpetraWrappers::Vector<float> &) const;
3125 # endif
3126 # endif
3127 
3128  template void
3132 
3133  template void
3135 
3136  template void
3138  const ::Vector<double> &) const;
3139 
3140  template void
3143  const ::LinearAlgebra::distributed::Vector<double> &) const;
3144 
3145 # ifdef DEAL_II_TRILINOS_WITH_TPETRA
3146 # if defined(HAVE_TPETRA_INST_DOUBLE)
3147  template void
3150  const ::LinearAlgebra::TpetraWrappers::Vector<double> &) const;
3151 # endif
3152 
3153 # if defined(HAVE_TPETRA_INST_FLOAT)
3154  template void
3157  const ::LinearAlgebra::TpetraWrappers::Vector<float> &) const;
3158 # endif
3159 # endif
3160 
3161  template void
3165 } // namespace TrilinosWrappers
3166 # endif // DOXYGEN
3167 
3169 
3170 #endif // DEAL_II_WITH_TRILINOS
const IndexSet & row_index_set() const
size_type row_length(const size_type row) const
size_type column_number(const size_type row, const size_type index) const
size_type size() const
Definition: index_set.h:1761
bool is_element(const size_type index) const
Definition: index_set.h:1879
Epetra_Map make_trilinos_map(const MPI_Comm communicator=MPI_COMM_WORLD, const bool overlapping=false) const
Definition: index_set.cc:1037
size_type n_rows() const
size_type n_cols() const
void reinit(const size_type m, const size_type n, const ArrayView< const unsigned int > &row_lengths)
const Epetra_BlockMap & trilinos_partitioner() const
void reinit(const Vector &v, const bool omit_zeroing_entries=false, const bool allow_different_maps=false)
std::shared_ptr< std::vector< TrilinosScalar > > value_cache
std::shared_ptr< std::vector< size_type > > colnum_cache
void set(const size_type i, const size_type j, const TrilinosScalar value)
std::unique_ptr< Epetra_Map > column_space_map
SparseMatrix & operator*=(const TrilinosScalar factor)
void mmult(SparseMatrix &C, const SparseMatrix &B, const MPI::Vector &V=MPI::Vector()) const
std::unique_ptr< Epetra_Export > nonlocal_matrix_exporter
void compress(VectorOperation::values operation)
std::unique_ptr< Epetra_FECrsMatrix > matrix
const Epetra_CrsMatrix & trilinos_matrix() const
::types::global_dof_index size_type
TrilinosScalar matrix_norm_square(const MPI::Vector &v) const
IndexSet locally_owned_range_indices() const
void print(std::ostream &out, const bool write_extended_trilinos_info=false) const
void Tmmult(SparseMatrix &C, const SparseMatrix &B, const MPI::Vector &V=MPI::Vector()) const
std::enable_if_t< std::is_same_v< typename VectorType::value_type, TrilinosScalar > > Tvmult(VectorType &dst, const VectorType &src) const
void clear_row(const size_type row, const TrilinosScalar new_diag_value=0)
std::pair< size_type, size_type > local_range() const
const Epetra_CrsGraph & trilinos_sparsity_pattern() const
void reinit(const SparsityPatternType &sparsity_pattern)
void vmult_add(VectorType &dst, const VectorType &src) const
SparseMatrix & operator/=(const TrilinosScalar factor)
void Tvmult_add(VectorType &dst, const VectorType &src) const
void clear_rows(const std::vector< size_type > &rows, const TrilinosScalar new_diag_value=0)
TrilinosScalar el(const size_type i, const size_type j) const
IndexSet locally_owned_domain_indices() const
TrilinosScalar residual(MPI::Vector &dst, const MPI::Vector &x, const MPI::Vector &b) const
bool in_local_range(const size_type index) const
void copy_from(const SparseMatrix &source)
unsigned int row_length(const size_type row) const
std::uint64_t n_nonzero_elements() const
std::enable_if_t< std::is_same_v< typename VectorType::value_type, TrilinosScalar > > vmult(VectorType &dst, const VectorType &src) const
TrilinosScalar diag_element(const size_type i) const
void add(const size_type i, const size_type j, const TrilinosScalar value)
unsigned int local_size() const
TrilinosScalar operator()(const size_type i, const size_type j) const
TrilinosScalar matrix_scalar_product(const MPI::Vector &u, const MPI::Vector &v) const
SparseMatrix & operator=(const SparseMatrix &)=delete
std::unique_ptr< Epetra_CrsMatrix > nonlocal_matrix
std::unique_ptr< Epetra_CrsGraph > nonlocal_graph
const Epetra_Map & domain_partitioner() const
const Epetra_FECrsGraph & trilinos_sparsity_pattern() const
std::function< void(VectorType &, const VectorType &)> inv_Tvmult
virtual int ApplyInverse(const VectorType &Y, VectorType &X) const override
virtual int Apply(const VectorType &X, VectorType &Y) const override
std::function< void(VectorType &, const VectorType &)> inv_vmult
Definition: vector.h:110
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:477
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:478
DerivativeForm< 1, spacedim, dim, Number > transpose(const DerivativeForm< 1, dim, spacedim, Number > &DF)
Point< 2 > second
Definition: grid_out.cc:4615
Point< 2 > first
Definition: grid_out.cc:4614
static ::ExceptionBase & ExcTrilinosError(int arg1)
static ::ExceptionBase & ExcAccessToNonlocalRow(std::size_t arg1)
static ::ExceptionBase & ExcDivideByZero()
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcInvalidIndex(size_type arg1, size_type arg2)
static ::ExceptionBase & ExcSourceEqualsDestination()
static ::ExceptionBase & ExcScalarAssignmentOnlyForZeroValue()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
#define Assert(cond, exc)
Definition: exceptions.h:1631
static ::ExceptionBase & ExcNotImplemented()
#define AssertIsFinite(number)
Definition: exceptions.h:1915
static ::ExceptionBase & ExcTrilinosError(int arg1)
static ::ExceptionBase & ExcAccessToNonLocalElement(size_type arg1, size_type arg2, size_type arg3, size_type arg4)
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1820
#define AssertIndexRange(index, range)
Definition: exceptions.h:1888
static ::ExceptionBase & ExcIO()
static ::ExceptionBase & ExcMatrixNotCompressed()
static ::ExceptionBase & ExcAccessToNonPresentElement(size_type arg1, size_type arg2)
static ::ExceptionBase & ExcMessage(std::string arg1)
static ::ExceptionBase & ExcNotQuadratic()
#define AssertThrow(cond, exc)
Definition: exceptions.h:1732
IndexSet complete_index_set(const IndexSet::size_type N)
Definition: index_set.h:1189
Expression fabs(const Expression &x)
typename ::internal::FEInterfaceViews::ViewType< dim, spacedim, Extractor >::type View
@ matrix
Contents is actually a matrix.
static const char V
void swap(MemorySpaceData< T, MemorySpace > &u, MemorySpaceData< T, MemorySpace > &v)
SymmetricTensor< 2, dim, Number > C(const Tensor< 2, dim, Number > &F)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
TrilinosPayload operator+(const TrilinosPayload &first_op, const TrilinosPayload &second_op)
TrilinosPayload operator*(const TrilinosPayload &first_op, const TrilinosPayload &second_op)
void check_vector_map_equality(const Epetra_CrsMatrix &, const VectorType &, const VectorType &)
VectorType::value_type * begin(VectorType &V)
VectorType::value_type * end(VectorType &V)
::types::global_dof_index size_type
void perform_mmult(const SparseMatrix &inputleft, const SparseMatrix &inputright, SparseMatrix &result, const MPI::Vector &V, const bool transpose_left)
TrilinosWrappers::types::int_type min_my_gid(const Epetra_BlockMap &map)
TrilinosWrappers::types::int64_type n_global_elements(const Epetra_BlockMap &map)
TrilinosWrappers::types::int_type global_column_index(const Epetra_CrsMatrix &matrix, const ::types::global_dof_index i)
TrilinosWrappers::types::int_type max_my_gid(const Epetra_BlockMap &map)
TrilinosWrappers::types::int_type global_row_index(const Epetra_CrsMatrix &matrix, const ::types::global_dof_index i)
TrilinosWrappers::types::int_type n_global_cols(const Epetra_CrsGraph &graph)
const Epetra_Comm & comm_self()
Iterator lower_bound(Iterator first, Iterator last, const T &val)
Definition: utilities.h:1008
Definition: types.h:33
unsigned int global_dof_index
Definition: types.h:82