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