Reference documentation for deal.II version 8.4.2
trilinos_sparse_matrix.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2008 - 2015 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 at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 #include <deal.II/lac/trilinos_sparse_matrix.h>
17 
18 #ifdef DEAL_II_WITH_TRILINOS
19 
20 # include <deal.II/base/utilities.h>
21 # include <deal.II/lac/sparse_matrix.h>
22 # include <deal.II/lac/trilinos_sparsity_pattern.h>
23 # include <deal.II/lac/sparsity_pattern.h>
24 # include <deal.II/lac/dynamic_sparsity_pattern.h>
25 # include <deal.II/lac/sparsity_tools.h>
26 # include <deal.II/lac/parallel_vector.h>
27 
28 DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
29 # include <Epetra_Export.h>
30 # include <ml_epetra_utils.h>
31 # include <ml_struct.h>
32 # include <Teuchos_RCP.hpp>
34 
35 DEAL_II_NAMESPACE_OPEN
36 
37 namespace TrilinosWrappers
38 {
39  namespace
40  {
41 #ifndef DEAL_II_WITH_64BIT_INDICES
42  // define a helper function that queries the size of an Epetra_Map object
43  // by calling either the 32- or 64-bit function necessary, and returns the
44  // result in the correct data type so that we can use it in calling other
45  // Epetra member functions that are overloaded by index type
46  int n_global_elements (const Epetra_BlockMap &map)
47  {
48  return map.NumGlobalElements();
49  }
50 
51  int min_my_gid(const Epetra_BlockMap &map)
52  {
53  return map.MinMyGID();
54  }
55 
56  int max_my_gid(const Epetra_BlockMap &map)
57  {
58  return map.MaxMyGID();
59  }
60 
61  int n_global_cols(const Epetra_CrsGraph &graph)
62  {
63  return graph.NumGlobalCols();
64  }
65 
66  int global_column_index(const Epetra_CrsMatrix &matrix, int i)
67  {
68  return matrix.GCID(i);
69  }
70 
71  int global_row_index(const Epetra_CrsMatrix &matrix, int i)
72  {
73  return matrix.GRID(i);
74  }
75 #else
76  // define a helper function that queries the size of an Epetra_Map object
77  // by calling either the 32- or 64-bit function necessary, and returns the
78  // result in the correct data type so that we can use it in calling other
79  // Epetra member functions that are overloaded by index type
80  long long int n_global_elements (const Epetra_BlockMap &map)
81  {
82  return map.NumGlobalElements64();
83  }
84 
85  long long int min_my_gid(const Epetra_BlockMap &map)
86  {
87  return map.MinMyGID64();
88  }
89 
90  long long int max_my_gid(const Epetra_BlockMap &map)
91  {
92  return map.MaxMyGID64();
93  }
94 
95  long long int n_global_cols(const Epetra_CrsGraph &graph)
96  {
97  return graph.NumGlobalCols64();
98  }
99 
100  long long int global_column_index(const Epetra_CrsMatrix &matrix, int i)
101  {
102  return matrix.GCID64(i);
103  }
104 
105  long long int global_row_index(const Epetra_CrsMatrix &matrix, int i)
106  {
107  return matrix.GRID64(i);
108  }
109 #endif
110  }
111 
112 
113  namespace SparseMatrixIterators
114  {
115  void
117  {
118  // if we are asked to visit the past-the-end line, then simply
119  // release all our caches and go on with life.
120  //
121  // do the same if the row we're supposed to visit is not locally
122  // owned. this is simply going to make non-locally owned rows
123  // look like they're empty
124  if ((this->a_row == matrix->m())
125  ||
126  (matrix->in_local_range (this->a_row) == false))
127  {
128  colnum_cache.reset ();
129  value_cache.reset ();
130 
131  return;
132  }
133 
134  // get a representation of the present row
135  int ncols;
136  TrilinosWrappers::types::int_type colnums = matrix->n();
137  if (value_cache.get() == 0)
138  {
139  value_cache.reset (new std::vector<TrilinosScalar> (matrix->n()));
140  colnum_cache.reset (new std::vector<size_type> (matrix->n()));
141  }
142  else
143  {
144  value_cache->resize (matrix->n());
145  colnum_cache->resize (matrix->n());
146  }
147 
148  int ierr = matrix->trilinos_matrix().
149  ExtractGlobalRowCopy((TrilinosWrappers::types::int_type)this->a_row,
150  colnums,
151  ncols, &((*value_cache)[0]),
152  reinterpret_cast<TrilinosWrappers::types::int_type *>(&((*colnum_cache)[0])));
153  value_cache->resize (ncols);
154  colnum_cache->resize (ncols);
155  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
156 
157  // copy it into our caches if the
158  // line isn't empty. if it is, then
159  // we've done something wrong, since
160  // we shouldn't have initialized an
161  // iterator for an empty line (what
162  // would it point to?)
163  }
164  }
165 
166 
167  // The constructor is actually the
168  // only point where we have to check
169  // whether we build a serial or a
170  // parallel Trilinos matrix.
171  // Actually, it does not even matter
172  // how many threads there are, but
173  // only if we use an MPI compiler or
174  // a standard compiler. So, even one
175  // thread on a configuration with
176  // MPI will still get a parallel
177  // interface.
179  :
180  column_space_map (new Epetra_Map (0, 0,
181  Utilities::Trilinos::comm_self())),
182  matrix (new Epetra_FECrsMatrix(View, *column_space_map,
183  *column_space_map, 0)),
184  last_action (Zero),
185  compressed (true)
186  {
187  matrix->FillComplete();
188  }
189 
190 
191 
192  SparseMatrix::SparseMatrix (const Epetra_Map &input_map,
193  const size_type n_max_entries_per_row)
194  :
195  column_space_map (new Epetra_Map (input_map)),
196  matrix (new Epetra_FECrsMatrix(Copy, *column_space_map,
197  TrilinosWrappers::types::int_type(n_max_entries_per_row), false)),
198  last_action (Zero),
199  compressed (false)
200  {}
201 
202 
203 
204  SparseMatrix::SparseMatrix (const Epetra_Map &input_map,
205  const std::vector<unsigned int> &n_entries_per_row)
206  :
207  column_space_map (new Epetra_Map (input_map)),
208  matrix (new Epetra_FECrsMatrix
209  (Copy, *column_space_map,
210  (int *)const_cast<unsigned int *>(&(n_entries_per_row[0])),
211  false)),
212  last_action (Zero),
213  compressed (false)
214  {}
215 
216 
217 
218  SparseMatrix::SparseMatrix (const Epetra_Map &input_row_map,
219  const Epetra_Map &input_col_map,
220  const size_type n_max_entries_per_row)
221  :
222  column_space_map (new Epetra_Map (input_col_map)),
223  matrix (new Epetra_FECrsMatrix(Copy, input_row_map,
224  TrilinosWrappers::types::int_type(n_max_entries_per_row), false)),
225  last_action (Zero),
226  compressed (false)
227  {}
228 
229 
230 
231  SparseMatrix::SparseMatrix (const Epetra_Map &input_row_map,
232  const Epetra_Map &input_col_map,
233  const std::vector<unsigned int> &n_entries_per_row)
234  :
235  column_space_map (new Epetra_Map (input_col_map)),
236  matrix (new Epetra_FECrsMatrix(Copy, input_row_map,
237  (int *)const_cast<unsigned int *>(&(n_entries_per_row[0])),
238  false)),
239  last_action (Zero),
240  compressed (false)
241  {}
242 
243 
244 
246  const size_type n,
247  const unsigned int n_max_entries_per_row)
248  :
249  column_space_map (new Epetra_Map (static_cast<TrilinosWrappers::types::int_type>(n), 0,
250  Utilities::Trilinos::comm_self())),
251 
252  // on one processor only, we know how the
253  // columns of the matrix will be
254  // distributed (everything on one
255  // processor), so we can hand in this
256  // information to the constructor. we
257  // can't do so in parallel, where the
258  // information from columns is only
259  // available when entries have been added
260  matrix (new Epetra_FECrsMatrix(Copy,
261  Epetra_Map (static_cast<TrilinosWrappers::types::int_type>(m), 0,
262  Utilities::Trilinos::comm_self()),
264  n_max_entries_per_row,
265  false)),
266  last_action (Zero),
267  compressed (false)
268  {}
269 
270 
271 
273  const size_type n,
274  const std::vector<unsigned int> &n_entries_per_row)
275  :
276  column_space_map (new Epetra_Map (static_cast<TrilinosWrappers::types::int_type>(n), 0,
277  Utilities::Trilinos::comm_self())),
278  matrix (new Epetra_FECrsMatrix(Copy,
279  Epetra_Map (static_cast<TrilinosWrappers::types::int_type>(m), 0,
280  Utilities::Trilinos::comm_self()),
282  (int *)const_cast<unsigned int *>(&(n_entries_per_row[0])),
283  false)),
284  last_action (Zero),
285  compressed (false)
286  {}
287 
288 
289 
290  SparseMatrix::SparseMatrix (const IndexSet &parallel_partitioning,
291  const MPI_Comm &communicator,
292  const unsigned int n_max_entries_per_row)
293  :
294  column_space_map (new Epetra_Map(parallel_partitioning.
295  make_trilinos_map(communicator, false))),
296  matrix (new Epetra_FECrsMatrix(Copy,
298  n_max_entries_per_row,
299  false)),
300  last_action (Zero),
301  compressed (false)
302  {}
303 
304 
305 
306  SparseMatrix::SparseMatrix (const IndexSet &parallel_partitioning,
307  const MPI_Comm &communicator,
308  const std::vector<unsigned int> &n_entries_per_row)
309  :
310  column_space_map (new Epetra_Map(parallel_partitioning.
311  make_trilinos_map(communicator, false))),
312  matrix (new Epetra_FECrsMatrix(Copy,
314  (int *)const_cast<unsigned int *>(&(n_entries_per_row[0])),
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 size_type n_max_entries_per_row)
326  :
327  column_space_map (new Epetra_Map(col_parallel_partitioning.
328  make_trilinos_map(communicator, false))),
329  matrix (new Epetra_FECrsMatrix(Copy,
330  row_parallel_partitioning.
331  make_trilinos_map(communicator, false),
332  n_max_entries_per_row,
333  false)),
334  last_action (Zero),
335  compressed (false)
336  {}
337 
338 
339 
340  SparseMatrix::SparseMatrix (const IndexSet &row_parallel_partitioning,
341  const IndexSet &col_parallel_partitioning,
342  const MPI_Comm &communicator,
343  const std::vector<unsigned int> &n_entries_per_row)
344  :
345  column_space_map (new Epetra_Map(col_parallel_partitioning.
346  make_trilinos_map(communicator, false))),
347  matrix (new Epetra_FECrsMatrix(Copy,
348  row_parallel_partitioning.
349  make_trilinos_map(communicator, false),
350  (int *)const_cast<unsigned int *>(&(n_entries_per_row[0])),
351  false)),
352  last_action (Zero),
353  compressed (false)
354  {}
355 
356 
357 
358  SparseMatrix::SparseMatrix (const SparsityPattern &sparsity_pattern)
359  :
360  column_space_map (new Epetra_Map (sparsity_pattern.domain_partitioner())),
361  matrix (new Epetra_FECrsMatrix(Copy,
362  sparsity_pattern.trilinos_sparsity_pattern(),
363  false)),
364  last_action (Zero),
365  compressed (true)
366  {
367  Assert(sparsity_pattern.trilinos_sparsity_pattern().Filled() == true,
368  ExcMessage("The Trilinos sparsity pattern has not been compressed."));
369  compress(VectorOperation::insert);
370  }
371 
372 
373 
375  {}
376 
377 
378 
379  void
381  {
382  if (this == &rhs)
383  return;
384 
385  nonlocal_matrix.reset();
386  nonlocal_matrix_exporter.reset();
387 
388  // check whether we need to update the whole matrix layout (we have
389  // different maps or if we detect a row where the columns of the two
390  // matrices do not match)
391  bool needs_deep_copy =
392  !matrix->RowMap().SameAs(rhs.matrix->RowMap()) ||
393  !matrix->ColMap().SameAs(rhs.matrix->ColMap()) ||
394  !matrix->DomainMap().SameAs(rhs.matrix->DomainMap()) ||
396  if (!needs_deep_copy)
397  {
398  const std::pair<size_type, size_type>
399  local_range = rhs.local_range();
400 
401  int ierr;
402  // Try to copy all the rows of the matrix one by one. In case of error
403  // (i.e., the column indices are different), we need to abort and blow
404  // away the matrix.
405  for (size_type row=local_range.first; row < local_range.second; ++row)
406  {
407  const int row_local =
408  matrix->RowMap().LID(static_cast<TrilinosWrappers::types::int_type>(row));
409 
410  int n_entries, rhs_n_entries;
411  TrilinosScalar *value_ptr, *rhs_value_ptr;
412  int *index_ptr, *rhs_index_ptr;
413  ierr = rhs.matrix->ExtractMyRowView (row_local, rhs_n_entries,
414  rhs_value_ptr, rhs_index_ptr);
415  (void)ierr;
416  Assert (ierr == 0, ExcTrilinosError(ierr));
417 
418  ierr = matrix->ExtractMyRowView (row_local, n_entries, value_ptr,
419  index_ptr);
420  Assert (ierr == 0, ExcTrilinosError(ierr));
421 
422  if (n_entries != rhs_n_entries ||
423  std::memcmp(static_cast<void *>(index_ptr),
424  static_cast<void *>(rhs_index_ptr),
425  sizeof(int)*n_entries) != 0)
426  {
427  needs_deep_copy = true;
428  break;
429  }
430 
431  for (int i=0; i<n_entries; ++i)
432  value_ptr[i] = rhs_value_ptr[i];
433  }
434  }
435 
436  if (needs_deep_copy)
437  {
438  column_space_map.reset (new Epetra_Map (rhs.domain_partitioner()));
439 
440  // release memory before reallocation
441  matrix.reset ();
442  matrix.reset (new Epetra_FECrsMatrix(*rhs.matrix));
443 
444  matrix->FillComplete(*column_space_map, matrix->RowMap());
445  }
446 
447  if (rhs.nonlocal_matrix.get() != 0)
448  nonlocal_matrix.reset(new Epetra_CrsMatrix(Copy, rhs.nonlocal_matrix->Graph()));
449  }
450 
451 
452 
453  namespace
454  {
456 
457  template <typename SparsityPatternType>
458  void
459  reinit_matrix (const Epetra_Map &input_row_map,
460  const Epetra_Map &input_col_map,
461  const SparsityPatternType &sparsity_pattern,
462  const bool exchange_data,
463  std_cxx11::shared_ptr<Epetra_Map> &column_space_map,
464  std_cxx11::shared_ptr<Epetra_FECrsMatrix> &matrix,
465  std_cxx11::shared_ptr<Epetra_CrsMatrix> &nonlocal_matrix,
466  std_cxx11::shared_ptr<Epetra_Export> &nonlocal_matrix_exporter)
467  {
468  // release memory before reallocation
469  matrix.reset();
470  nonlocal_matrix.reset();
471  nonlocal_matrix_exporter.reset();
472 
473  if (input_row_map.Comm().MyPID() == 0)
474  {
475  AssertDimension (sparsity_pattern.n_rows(),
476  static_cast<size_type>(n_global_elements(input_row_map)));
477  AssertDimension (sparsity_pattern.n_cols(),
478  static_cast<size_type>(n_global_elements(input_col_map)));
479  }
480 
481  column_space_map.reset (new Epetra_Map (input_col_map));
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 (input_row_map, input_col_map,
491  sparsity_pattern, exchange_data);
492  matrix.reset (new Epetra_FECrsMatrix
493  (Copy, trilinos_sparsity.trilinos_sparsity_pattern(), false));
494 
495  return;
496  }
497 
498  const size_type first_row = min_my_gid(input_row_map),
499  last_row = max_my_gid(input_row_map)+1;
500  std::vector<int> n_entries_per_row(last_row-first_row);
501 
502  for (size_type row=first_row; row<last_row; ++row)
503  n_entries_per_row[row-first_row] = sparsity_pattern.row_length(row);
504 
505  // The deal.II notation of a Sparsity pattern corresponds to the Epetra
506  // concept of a Graph. Hence, we generate a graph by copying the
507  // sparsity pattern into it, and then build up the matrix from the
508  // graph. This is considerable faster than directly filling elements
509  // into the matrix. Moreover, it consumes less memory, since the
510  // internal reordering is done on ints only, and we can leave the
511  // doubles aside.
512 
513  // for more than one processor, need to specify only row map first and
514  // let the matrix entries decide about the column map (which says which
515  // columns are present in the matrix, not to be confused with the
516  // col_map that tells how the domain dofs of the matrix will be
517  // distributed). for only one processor, we can directly assign the
518  // columns as well. Compare this with bug # 4123 in the Sandia Bugzilla.
519  std_cxx11::shared_ptr<Epetra_CrsGraph> graph;
520  if (input_row_map.Comm().NumProc() > 1)
521  graph.reset (new Epetra_CrsGraph (Copy, input_row_map,
522  &n_entries_per_row[0], true));
523  else
524  graph.reset (new Epetra_CrsGraph (Copy, input_row_map, input_col_map,
525  &n_entries_per_row[0], true));
526 
527  // This functions assumes that the sparsity pattern sits on all
528  // processors (completely). The parallel version uses an Epetra graph
529  // that is already distributed.
530 
531  // now insert the indices
532  std::vector<TrilinosWrappers::types::int_type> row_indices;
533 
534  for (size_type row=first_row; row<last_row; ++row)
535  {
536  const int row_length = sparsity_pattern.row_length(row);
537  if (row_length == 0)
538  continue;
539 
540  row_indices.resize (row_length, -1);
541  {
542  typename SparsityPatternType::iterator p = sparsity_pattern.begin(row);
543  for (size_type col=0; p != sparsity_pattern.end(row); ++p, ++col)
544  row_indices[col] = p->column();
545  }
546  graph->Epetra_CrsGraph::InsertGlobalIndices (row, row_length,
547  &row_indices[0]);
548  }
549 
550  // Eventually, optimize the graph structure (sort indices, make memory
551  // contiguous, etc). note that the documentation of the function indeed
552  // states that we first need to provide the column (domain) map and then
553  // the row (range) map
554  graph->FillComplete(input_col_map, input_row_map);
555  graph->OptimizeStorage();
556 
557  // check whether we got the number of columns right.
558  AssertDimension (sparsity_pattern.n_cols(),
559  static_cast<size_type>(n_global_cols(*graph)));
560  (void)n_global_cols;
561 
562  // And now finally generate the matrix.
563  matrix.reset (new Epetra_FECrsMatrix(Copy, *graph, false));
564  }
565 
566 
567 
568  // for the non-local graph, we need to circumvent the problem that some
569  // processors will not add into the non-local graph at all: We do not want
570  // to insert dummy elements on >5000 processors because that gets very
571  // slow. Thus, we set a flag in Epetra_CrsGraph that sets the correct
572  // flag. Since it is protected, we need to expose this information by
573  // deriving a class from Epetra_CrsGraph for the purpose of creating the
574  // data structure
575  class Epetra_CrsGraphMod : public Epetra_CrsGraph
576  {
577  public:
578  Epetra_CrsGraphMod (const Epetra_Map &row_map,
579  const int *n_entries_per_row)
580  :
581  Epetra_CrsGraph(Copy, row_map, n_entries_per_row, true)
582  {};
583 
584  void SetIndicesAreGlobal()
585  {
586  this->Epetra_CrsGraph::SetIndicesAreGlobal(true);
587  }
588  };
589 
590 
591  // specialization for DynamicSparsityPattern which can provide us with
592  // more information about the non-locally owned rows
593  template <>
594  void
595  reinit_matrix (const Epetra_Map &input_row_map,
596  const Epetra_Map &input_col_map,
597  const DynamicSparsityPattern &sparsity_pattern,
598  const bool exchange_data,
599  std_cxx11::shared_ptr<Epetra_Map> &column_space_map,
600  std_cxx11::shared_ptr<Epetra_FECrsMatrix> &matrix,
601  std_cxx11::shared_ptr<Epetra_CrsMatrix> &nonlocal_matrix,
602  std_cxx11::shared_ptr<Epetra_Export> &nonlocal_matrix_exporter)
603  {
604  matrix.reset();
605  nonlocal_matrix.reset();
606  nonlocal_matrix_exporter.reset();
607 
608  AssertDimension (sparsity_pattern.n_rows(),
609  static_cast<size_type>(n_global_elements(input_row_map)));
610  AssertDimension (sparsity_pattern.n_cols(),
611  static_cast<size_type>(n_global_elements(input_col_map)));
612 
613  column_space_map.reset (new Epetra_Map (input_col_map));
614 
615  IndexSet relevant_rows (sparsity_pattern.row_index_set());
616  // serial case
617  if (relevant_rows.size() == 0)
618  {
619  relevant_rows.set_size(n_global_elements(input_row_map));
620  relevant_rows.add_range(0, n_global_elements(input_row_map));
621  }
622  relevant_rows.compress();
623  Assert(relevant_rows.n_elements() >= static_cast<unsigned int>(input_row_map.NumMyElements()),
624  ExcMessage("Locally relevant rows of sparsity pattern must contain "
625  "all locally owned rows"));
626 
627  // check whether the relevant rows correspond to exactly the same map as
628  // the owned rows. In that case, do not create the nonlocal graph and
629  // fill the columns by demand
630  bool have_ghost_rows = false;
631  {
632  std::vector<::types::global_dof_index> indices;
633  relevant_rows.fill_index_vector(indices);
634  Epetra_Map relevant_map (TrilinosWrappers::types::int_type(-1),
635  TrilinosWrappers::types::int_type(relevant_rows.n_elements()),
636  (indices.empty() ? 0 :
637  reinterpret_cast<TrilinosWrappers::types::int_type *>(&indices[0])),
638  0, input_row_map.Comm());
639  if (relevant_map.SameAs(input_row_map))
640  have_ghost_rows = false;
641  else
642  have_ghost_rows = true;
643  }
644 
645  const unsigned int n_rows = relevant_rows.n_elements();
646  std::vector<TrilinosWrappers::types::int_type> ghost_rows;
647  std::vector<int> n_entries_per_row(input_row_map.NumMyElements());
648  std::vector<int> n_entries_per_ghost_row;
649  for (unsigned int i=0, own=0; i<n_rows; ++i)
650  {
651  const TrilinosWrappers::types::int_type global_row =
652  relevant_rows.nth_index_in_set(i);
653  if (input_row_map.MyGID(global_row))
654  n_entries_per_row[own++] = sparsity_pattern.row_length(global_row);
655  else if (sparsity_pattern.row_length(global_row) > 0)
656  {
657  ghost_rows.push_back(global_row);
658  n_entries_per_ghost_row.push_back(sparsity_pattern.row_length(global_row));
659  }
660  }
661 
662  Epetra_Map off_processor_map(-1, ghost_rows.size(),
663  (ghost_rows.size()>0)?(&ghost_rows[0]):NULL,
664  0, input_row_map.Comm());
665 
666  std_cxx11::shared_ptr<Epetra_CrsGraph> graph;
667  std_cxx11::shared_ptr<Epetra_CrsGraphMod> nonlocal_graph;
668  if (input_row_map.Comm().NumProc() > 1)
669  {
670  graph.reset (new Epetra_CrsGraph (Copy, input_row_map,
671  (n_entries_per_row.size()>0)?(&n_entries_per_row[0]):NULL,
672  exchange_data ? false : true));
673  if (have_ghost_rows == true)
674  nonlocal_graph.reset (new Epetra_CrsGraphMod (off_processor_map,
675  &n_entries_per_ghost_row[0]));
676  }
677  else
678  graph.reset (new Epetra_CrsGraph (Copy, input_row_map, input_col_map,
679  (n_entries_per_row.size()>0)?(&n_entries_per_row[0]):NULL,
680  true));
681 
682  // now insert the indices, select between the right matrix
683  std::vector<TrilinosWrappers::types::int_type> row_indices;
684 
685  for (unsigned int i=0; i<n_rows; ++i)
686  {
687  const TrilinosWrappers::types::int_type global_row =
688  relevant_rows.nth_index_in_set(i);
689  const int row_length = sparsity_pattern.row_length(global_row);
690  if (row_length == 0)
691  continue;
692 
693  row_indices.resize (row_length, -1);
694  for (int col=0; col < row_length; ++col)
695  row_indices[col] = sparsity_pattern.column_number(global_row, col);
696 
697  if (input_row_map.MyGID(global_row))
698  graph->InsertGlobalIndices (global_row, row_length, &row_indices[0]);
699  else
700  {
701  Assert(nonlocal_graph.get() != 0, ExcInternalError());
702  nonlocal_graph->InsertGlobalIndices (global_row, row_length,
703  &row_indices[0]);
704  }
705  }
706 
707  // finalize nonlocal graph and create nonlocal matrix
708  if (nonlocal_graph.get() != 0)
709  {
710  // must make sure the IndicesAreGlobal flag is set on all processors
711  // because some processors might not call InsertGlobalIndices (and
712  // we do not want to insert dummy indices on all processors for
713  // large-scale simulations due to the bad impact on performance)
714  nonlocal_graph->SetIndicesAreGlobal();
715  Assert(nonlocal_graph->IndicesAreGlobal() == true,
716  ExcInternalError());
717  nonlocal_graph->FillComplete(input_col_map, input_row_map);
718  nonlocal_graph->OptimizeStorage();
719 
720  // insert data from nonlocal graph into the final sparsity pattern
721  if (exchange_data)
722  {
723  Epetra_Export exporter(nonlocal_graph->RowMap(), input_row_map);
724  int ierr = graph->Export(*nonlocal_graph, exporter, Add);
725  (void)ierr;
726  Assert (ierr==0, ExcTrilinosError(ierr));
727  }
728 
729  nonlocal_matrix.reset (new Epetra_CrsMatrix(Copy, *nonlocal_graph));
730  }
731 
732  graph->FillComplete(input_col_map, input_row_map);
733  graph->OptimizeStorage();
734 
735  AssertDimension (sparsity_pattern.n_cols(),static_cast<size_type>(
736  n_global_cols(*graph)));
737 
738  matrix.reset (new Epetra_FECrsMatrix(Copy, *graph, false));
739  }
740  }
741 
742 
743 
744  template <typename SparsityPatternType>
745  void
746  SparseMatrix::reinit (const SparsityPatternType &sparsity_pattern)
747  {
748  const Epetra_Map rows (static_cast<TrilinosWrappers::types::int_type>(sparsity_pattern.n_rows()),
749  0,
751  const Epetra_Map columns (static_cast<TrilinosWrappers::types::int_type>(sparsity_pattern.n_cols()),
752  0,
754 
755  reinit_matrix (rows, columns, sparsity_pattern, false,
758  }
759 
760 
761 
762  template <typename SparsityPatternType>
763  void
764  SparseMatrix::reinit (const Epetra_Map &input_map,
765  const SparsityPatternType &sparsity_pattern,
766  const bool exchange_data)
767  {
768  reinit_matrix (input_map, input_map, sparsity_pattern, exchange_data,
771  }
772 
773 
774 
775 
776 
777 
778  template <typename SparsityPatternType>
779  inline
780  void SparseMatrix::reinit (const IndexSet &row_parallel_partitioning,
781  const IndexSet &col_parallel_partitioning,
782  const SparsityPatternType &sparsity_pattern,
783  const MPI_Comm &communicator,
784  const bool exchange_data)
785  {
786  Epetra_Map row_map =
787  row_parallel_partitioning.make_trilinos_map (communicator, false);
788  Epetra_Map col_map =
789  col_parallel_partitioning.make_trilinos_map (communicator, false);
790  reinit_matrix (row_map, col_map, sparsity_pattern, exchange_data,
793 
794  // In the end, the matrix needs to be compressed in order to be really
795  // ready.
796  last_action = Zero;
797  compress(VectorOperation::insert);
798  }
799 
800 
801 
802  template <typename SparsityPatternType>
803  inline
804  void SparseMatrix::reinit (const Epetra_Map &row_map,
805  const Epetra_Map &col_map,
806  const SparsityPatternType &sparsity_pattern,
807  const bool exchange_data)
808  {
809  reinit_matrix (row_map, col_map, sparsity_pattern, exchange_data,
812 
813  // In the end, the matrix needs to be compressed in order to be really
814  // ready.
815  last_action = Zero;
816  compress(VectorOperation::insert);
817  }
818 
819 
820 
821 
822  void
823  SparseMatrix::reinit (const SparsityPattern &sparsity_pattern)
824  {
825  matrix.reset ();
826  nonlocal_matrix_exporter.reset();
827 
828  // reinit with a (parallel) Trilinos sparsity pattern.
829  column_space_map.reset (new Epetra_Map
830  (sparsity_pattern.domain_partitioner()));
831  matrix.reset (new Epetra_FECrsMatrix
832  (Copy, sparsity_pattern.trilinos_sparsity_pattern(), false));
833 
834  if (sparsity_pattern.nonlocal_graph.get() != 0)
835  nonlocal_matrix.reset (new Epetra_CrsMatrix(Copy, *sparsity_pattern.nonlocal_graph));
836  else
837  nonlocal_matrix.reset ();
838 
839  last_action = Zero;
840  compress(VectorOperation::insert);
841  }
842 
843 
844 
845  void
846  SparseMatrix::reinit (const SparseMatrix &sparse_matrix)
847  {
848  if (this == &sparse_matrix)
849  return;
850 
851  column_space_map.reset (new Epetra_Map (sparse_matrix.domain_partitioner()));
852  matrix.reset ();
853  nonlocal_matrix_exporter.reset();
854  matrix.reset (new Epetra_FECrsMatrix
855  (Copy, sparse_matrix.trilinos_sparsity_pattern(), false));
856 
857  if (sparse_matrix.nonlocal_matrix != 0)
858  nonlocal_matrix.reset (new Epetra_CrsMatrix
859  (Copy, sparse_matrix.nonlocal_matrix->Graph()));
860  else
861  nonlocal_matrix.reset();
862 
863  last_action = Zero;
864  compress(VectorOperation::insert);
865  }
866 
867 
868 
869  template <typename number>
870  inline
871  void SparseMatrix::reinit (const IndexSet &row_parallel_partitioning,
872  const IndexSet &col_parallel_partitioning,
873  const ::SparseMatrix<number> &dealii_sparse_matrix,
874  const MPI_Comm &communicator,
875  const double drop_tolerance,
876  const bool copy_values,
877  const ::SparsityPattern *use_this_sparsity)
878  {
879  if (copy_values == false)
880  {
881  // in case we do not copy values, just
882  // call the other function.
883  if (use_this_sparsity == 0)
884  reinit (row_parallel_partitioning, col_parallel_partitioning,
885  dealii_sparse_matrix.get_sparsity_pattern(),
886  communicator, false);
887  else
888  reinit (row_parallel_partitioning, col_parallel_partitioning,
889  *use_this_sparsity, communicator, false);
890  return;
891  }
892 
893  const size_type n_rows = dealii_sparse_matrix.m();
894 
895  AssertDimension (row_parallel_partitioning.size(), n_rows);
896  AssertDimension (col_parallel_partitioning.size(), dealii_sparse_matrix.n());
897 
898  const ::SparsityPattern &sparsity_pattern =
899  (use_this_sparsity!=0)? *use_this_sparsity :
900  dealii_sparse_matrix.get_sparsity_pattern();
901 
902  if (matrix.get() == 0 ||
903  m() != n_rows ||
904  n_nonzero_elements() != sparsity_pattern.n_nonzero_elements())
905  {
906  reinit (row_parallel_partitioning, col_parallel_partitioning,
907  sparsity_pattern, communicator, false);
908  }
909 
910  // fill the values. the same as above: go through all rows of the
911  // matrix, and then all columns. since the sparsity patterns of the
912  // input matrix and the specified sparsity pattern might be different,
913  // need to go through the row for both these sparsity structures
914  // simultaneously in order to really set the correct values.
915  size_type maximum_row_length = matrix->MaxNumEntries();
916  std::vector<size_type> row_indices (maximum_row_length);
917  std::vector<TrilinosScalar> values (maximum_row_length);
918 
919  for (size_type row=0; row<n_rows; ++row)
920  // see if the row is locally stored on this processor
921  if (row_parallel_partitioning.is_element(row) == true)
922  {
923  ::SparsityPattern::iterator select_index =
924  sparsity_pattern.begin(row);
925  typename ::SparseMatrix<number>::const_iterator it =
926  dealii_sparse_matrix.begin(row);
927  size_type col = 0;
928  if (sparsity_pattern.n_rows() == sparsity_pattern.n_cols())
929  {
930  // optimized diagonal
931  AssertDimension(it->column(), row);
932  if (std::fabs(it->value()) > drop_tolerance)
933  {
934  values[col] = it->value();
935  row_indices[col++] = it->column();
936  }
937  ++select_index;
938  ++it;
939  }
940 
941  while (it != dealii_sparse_matrix.end(row) &&
942  select_index != sparsity_pattern.end(row))
943  {
944  while (select_index->column() < it->column() &&
945  select_index != sparsity_pattern.end(row))
946  ++select_index;
947  while (it->column() < select_index->column() &&
948  it != dealii_sparse_matrix.end(row))
949  ++it;
950 
951  if (it == dealii_sparse_matrix.end(row))
952  break;
953  if (std::fabs(it->value()) > drop_tolerance)
954  {
955  values[col] = it->value();
956  row_indices[col++] = it->column();
957  }
958  ++select_index;
959  ++it;
960  }
961  set (row, col, reinterpret_cast<size_type *>(&row_indices[0]),
962  &values[0], false);
963  }
964  compress(VectorOperation::insert);
965  }
966 
967 
968 
969  template <typename number>
970  void
971  SparseMatrix::reinit (const ::SparseMatrix<number> &dealii_sparse_matrix,
972  const double drop_tolerance,
973  const bool copy_values,
974  const ::SparsityPattern *use_this_sparsity)
975  {
976  reinit (complete_index_set(dealii_sparse_matrix.m()),
977  complete_index_set(dealii_sparse_matrix.n()),
978  dealii_sparse_matrix, MPI_COMM_SELF, drop_tolerance,
979  copy_values, use_this_sparsity);
980  }
981 
982 
983 
984  template <typename number>
985  void
986  SparseMatrix::reinit (const Epetra_Map &input_map,
987  const ::SparseMatrix<number> &dealii_sparse_matrix,
988  const double drop_tolerance,
989  const bool copy_values,
990  const ::SparsityPattern *use_this_sparsity)
991  {
992  reinit (IndexSet(input_map), IndexSet(input_map), dealii_sparse_matrix,
993  MPI_COMM_SELF, drop_tolerance, copy_values, use_this_sparsity);
994  }
995 
996 
997 
998  template <typename number>
999  void
1000  SparseMatrix::reinit (const Epetra_Map &input_row_map,
1001  const Epetra_Map &input_col_map,
1002  const ::SparseMatrix<number> &dealii_sparse_matrix,
1003  const double drop_tolerance,
1004  const bool copy_values,
1005  const ::SparsityPattern *use_this_sparsity)
1006  {
1007  reinit (IndexSet(input_row_map), IndexSet(input_col_map),
1008  dealii_sparse_matrix, MPI_COMM_SELF,
1009  drop_tolerance, copy_values, use_this_sparsity);
1010  }
1011 
1012 
1013 
1014  void
1015  SparseMatrix::reinit (const Epetra_CrsMatrix &input_matrix,
1016  const bool copy_values)
1017  {
1018  Assert (input_matrix.Filled()==true,
1019  ExcMessage("Input CrsMatrix has not called FillComplete()!"));
1020 
1021  column_space_map.reset (new Epetra_Map (input_matrix.DomainMap()));
1022 
1023  const Epetra_CrsGraph *graph = &input_matrix.Graph();
1024 
1025  nonlocal_matrix.reset();
1026  nonlocal_matrix_exporter.reset();
1027  matrix.reset ();
1028  matrix.reset (new Epetra_FECrsMatrix(Copy, *graph, false));
1029 
1030  matrix->FillComplete (*column_space_map, input_matrix.RangeMap(), true);
1031 
1032  if (copy_values == true)
1033  {
1034  // point to the first data entry in the two
1035  // matrices and copy the content
1036  const TrilinosScalar *in_values = input_matrix[0];
1037  TrilinosScalar *values = (*matrix)[0];
1038  const size_type my_nonzeros = input_matrix.NumMyNonzeros();
1039  std::memcpy (&values[0], &in_values[0],
1040  my_nonzeros*sizeof (TrilinosScalar));
1041  }
1042 
1043  last_action = Zero;
1044  compress(VectorOperation::insert);
1045  }
1046 
1047 
1048 
1049  void
1050  SparseMatrix::compress (::VectorOperation::values operation)
1051  {
1052 
1053  Epetra_CombineMode mode = last_action;
1054  if (last_action == Zero)
1055  {
1056  if ((operation==::VectorOperation::add) ||
1057  (operation==::VectorOperation::unknown))
1058  mode = Add;
1059  else if (operation==::VectorOperation::insert)
1060  mode = Insert;
1061  }
1062  else
1063  {
1064  Assert(
1065  ((last_action == Add) && (operation!=::VectorOperation::insert))
1066  ||
1067  ((last_action == Insert) && (operation!=::VectorOperation::add)),
1068  ExcMessage("Operation and argument to compress() do not match"));
1069  }
1070 
1071  // flush buffers
1072  int ierr;
1073  if (nonlocal_matrix.get() != 0 && mode == Add)
1074  {
1075  // do only export in case of an add() operation, otherwise the owning
1076  // processor must have set the correct entry
1077  nonlocal_matrix->FillComplete(*column_space_map, matrix->RowMap());
1078  if (nonlocal_matrix_exporter.get() == 0)
1080  (new Epetra_Export(nonlocal_matrix->RowMap(), matrix->RowMap()));
1081  ierr = matrix->Export(*nonlocal_matrix, *nonlocal_matrix_exporter, mode);
1082  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
1083  ierr = matrix->FillComplete(*column_space_map, matrix->RowMap());
1084  nonlocal_matrix->PutScalar(0);
1085  }
1086  else
1087  ierr = matrix->GlobalAssemble (*column_space_map, matrix->RowMap(),
1088  true, mode);
1089 
1090  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
1091 
1092  ierr = matrix->OptimizeStorage ();
1093  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
1094 
1095  last_action = Zero;
1096 
1097  compressed = true;
1098  }
1099 
1100 
1101 
1102  void
1104  {
1105  // When we clear the matrix, reset
1106  // the pointer and generate an
1107  // empty matrix.
1108  column_space_map.reset (new Epetra_Map (0, 0,
1110  matrix.reset (new Epetra_FECrsMatrix(View, *column_space_map, 0));
1111  nonlocal_matrix.reset();
1112  nonlocal_matrix_exporter.reset();
1113 
1114  matrix->FillComplete();
1115 
1116  compressed = true;
1117  }
1118 
1119 
1120 
1121  void
1123  const TrilinosScalar new_diag_value)
1124  {
1125  Assert (matrix->Filled()==true, ExcMatrixNotCompressed());
1126 
1127  // Only do this on the rows owned
1128  // locally on this processor.
1129  int local_row =
1130  matrix->LRID(static_cast<TrilinosWrappers::types::int_type>(row));
1131  if (local_row >= 0)
1132  {
1133  TrilinosScalar *values;
1134  int *col_indices;
1135  int num_entries;
1136  const int ierr = matrix->ExtractMyRowView(local_row, num_entries,
1137  values, col_indices);
1138  (void)ierr;
1139 
1140  Assert (ierr == 0,
1141  ExcTrilinosError(ierr));
1142 
1143  int *diag_find = std::find(col_indices,col_indices+num_entries,local_row);
1144  int diag_index = (int)(diag_find - col_indices);
1145 
1146  for (TrilinosWrappers::types::int_type j=0; j<num_entries; ++j)
1147  if (diag_index != j || new_diag_value == 0)
1148  values[j] = 0.;
1149 
1150  if (diag_find && std::fabs(values[diag_index]) == 0.0 &&
1151  new_diag_value != 0.0)
1152  values[diag_index] = new_diag_value;
1153  }
1154  }
1155 
1156 
1157 
1158  void
1159  SparseMatrix::clear_rows (const std::vector<size_type> &rows,
1160  const TrilinosScalar new_diag_value)
1161  {
1162  for (size_type row=0; row<rows.size(); ++row)
1163  clear_row(rows[row], new_diag_value);
1164  }
1165 
1166 
1167 
1168  TrilinosScalar
1170  const size_type j) const
1171  {
1172  // Extract local indices in
1173  // the matrix.
1174  int trilinos_i = matrix->LRID(static_cast<TrilinosWrappers::types::int_type>(i)),
1175  trilinos_j = 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, ExcAccessToNonLocalElement(i, j, local_range().first,
1187  local_range().second));
1188  }
1189  else
1190  {
1191  // Check whether the matrix has
1192  // already been transformed to local
1193  // indices.
1194  Assert (matrix->Filled(), ExcMatrixNotCompressed());
1195 
1196  // Prepare pointers for extraction
1197  // of a view of the row.
1198  int nnz_present = matrix->NumMyEntries(trilinos_i);
1199  int nnz_extracted;
1200  int *col_indices;
1201  TrilinosScalar *values;
1202 
1203  // Generate the view and make
1204  // sure that we have not generated
1205  // an error.
1206  // TODO Check that col_indices are int and not long long
1207  int ierr = matrix->ExtractMyRowView(trilinos_i, nnz_extracted,
1208  values, col_indices);
1209  (void)ierr;
1210  Assert (ierr==0, ExcTrilinosError(ierr));
1211 
1212  Assert (nnz_present == nnz_extracted,
1213  ExcDimensionMismatch(nnz_present, nnz_extracted));
1214 
1215  // Search the index where we
1216  // look for the value, and then
1217  // finally get it.
1218 
1219  int *el_find = std::find(col_indices, col_indices + nnz_present, trilinos_j);
1220 
1221  int local_col_index = (int)(el_find - col_indices);
1222 
1223  // This is actually the only
1224  // difference to the el(i,j)
1225  // function, which means that
1226  // we throw an exception in
1227  // this case instead of just
1228  // returning zero for an
1229  // element that is not present
1230  // in the sparsity pattern.
1231  if (local_col_index == nnz_present)
1232  {
1233  Assert (false, ExcInvalidIndex (i,j));
1234  }
1235  else
1236  value = values[local_col_index];
1237  }
1238 
1239  return value;
1240  }
1241 
1242 
1243 
1244  TrilinosScalar
1246  const size_type j) const
1247  {
1248  // Extract local indices in
1249  // the matrix.
1250  int trilinos_i = matrix->LRID(static_cast<TrilinosWrappers::types::int_type>(i)),
1251  trilinos_j = matrix->LCID(static_cast<TrilinosWrappers::types::int_type>(j));
1252  TrilinosScalar value = 0.;
1253 
1254  // If the data is not on the
1255  // present processor, we can't
1256  // continue. Just print out zero
1257  // as discussed in the
1258  // documentation of this
1259  // function. if you want error
1260  // checking, use operator().
1261  if ((trilinos_i == -1 ) || (trilinos_j == -1))
1262  return 0.;
1263  else
1264  {
1265  // Check whether the matrix
1266  // already is transformed to
1267  // local indices.
1268  Assert (matrix->Filled(), ExcMatrixNotCompressed());
1269 
1270  // Prepare pointers for extraction
1271  // of a view of the row.
1272  int nnz_present = matrix->NumMyEntries(trilinos_i);
1273  int nnz_extracted;
1274  int *col_indices;
1275  TrilinosScalar *values;
1276 
1277  // Generate the view and make
1278  // sure that we have not generated
1279  // an error.
1280  int ierr = matrix->ExtractMyRowView(trilinos_i, nnz_extracted,
1281  values, col_indices);
1282  (void)ierr;
1283  Assert (ierr==0, ExcTrilinosError(ierr));
1284 
1285  Assert (nnz_present == nnz_extracted,
1286  ExcDimensionMismatch(nnz_present, nnz_extracted));
1287 
1288  // Search the index where we
1289  // look for the value, and then
1290  // finally get it.
1291  int *el_find = std::find(col_indices, col_indices + nnz_present, trilinos_j);
1292 
1293  int local_col_index = (int)(el_find - col_indices);
1294 
1295 
1296  // This is actually the only
1297  // difference to the () function
1298  // querying (i,j), where we throw an
1299  // exception instead of just
1300  // returning zero for an element
1301  // that is not present in the
1302  // sparsity pattern.
1303  if (local_col_index == nnz_present)
1304  value = 0;
1305  else
1306  value = values[local_col_index];
1307  }
1308 
1309  return value;
1310  }
1311 
1312 
1313 
1314  TrilinosScalar
1316  {
1317  Assert (m() == n(), ExcNotQuadratic());
1318 
1319 #ifdef DEBUG
1320  // use operator() in debug mode because
1321  // it checks if this is a valid element
1322  // (in parallel)
1323  return operator()(i,i);
1324 #else
1325  // Trilinos doesn't seem to have a
1326  // more efficient way to access the
1327  // diagonal than by just using the
1328  // standard el(i,j) function.
1329  return el(i,i);
1330 #endif
1331  }
1332 
1333 
1334 
1335  unsigned int
1337  {
1338  Assert (row < m(), ExcInternalError());
1339 
1340  // get a representation of the
1341  // present row
1342  int ncols = -1;
1343  int local_row = matrix->LRID(static_cast<TrilinosWrappers::types::int_type>(row));
1344 
1345  // on the processor who owns this
1346  // row, we'll have a non-negative
1347  // value.
1348  if (local_row >= 0)
1349  {
1350  int ierr = matrix->NumMyRowEntries (local_row, ncols);
1351  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
1352  }
1353 
1354  return ncols;
1355  }
1356 
1357 
1358 
1359  void
1360  SparseMatrix::set (const std::vector<size_type> &row_indices,
1361  const std::vector<size_type> &col_indices,
1362  const FullMatrix<TrilinosScalar> &values,
1363  const bool elide_zero_values)
1364  {
1365  Assert (row_indices.size() == values.m(),
1366  ExcDimensionMismatch(row_indices.size(), values.m()));
1367  Assert (col_indices.size() == values.n(),
1368  ExcDimensionMismatch(col_indices.size(), values.n()));
1369 
1370  for (size_type i=0; i<row_indices.size(); ++i)
1371  set (row_indices[i], col_indices.size(), &col_indices[0], &values(i,0),
1372  elide_zero_values);
1373  }
1374 
1375 
1376 
1377  void
1379  const std::vector<size_type> &col_indices,
1380  const std::vector<TrilinosScalar> &values,
1381  const bool elide_zero_values)
1382  {
1383  Assert (col_indices.size() == values.size(),
1384  ExcDimensionMismatch(col_indices.size(), values.size()));
1385 
1386  set (row, col_indices.size(), &col_indices[0], &values[0],
1387  elide_zero_values);
1388  }
1389 
1390 
1391 
1392  void
1394  const size_type n_cols,
1395  const size_type *col_indices,
1396  const TrilinosScalar *values,
1397  const bool elide_zero_values)
1398  {
1399  AssertIndexRange(row, this->m());
1400 
1401  int ierr;
1402  if (last_action == Add)
1403  {
1404  ierr = matrix->GlobalAssemble (*column_space_map, matrix->RowMap(),
1405  true);
1406 
1407  Assert (ierr == 0, ExcTrilinosError(ierr));
1408  }
1409 
1410  last_action = Insert;
1411 
1412  TrilinosWrappers::types::int_type *col_index_ptr;
1413  TrilinosScalar *col_value_ptr;
1414  TrilinosWrappers::types::int_type n_columns;
1415 
1416  TrilinosScalar short_val_array[100];
1417  TrilinosWrappers::types::int_type short_index_array[100];
1418  std::vector<TrilinosScalar> long_val_array;
1419  std::vector<TrilinosWrappers::types::int_type> long_index_array;
1420 
1421 
1422  // If we don't elide zeros, the pointers are already available... need to
1423  // cast to non-const pointers as that is the format taken by Trilinos (but
1424  // we will not modify const data)
1425  if (elide_zero_values == false)
1426  {
1427  col_index_ptr = (TrilinosWrappers::types::int_type *)col_indices;
1428  col_value_ptr = const_cast<TrilinosScalar *>(values);
1429  n_columns = n_cols;
1430  }
1431  else
1432  {
1433  // Otherwise, extract nonzero values in each row and get the
1434  // respective indices.
1435  if (n_cols > 100)
1436  {
1437  long_val_array.resize(n_cols);
1438  long_index_array.resize(n_cols);
1439  col_index_ptr = &long_index_array[0];
1440  col_value_ptr = &long_val_array[0];
1441  }
1442  else
1443  {
1444  col_index_ptr = &short_index_array[0];
1445  col_value_ptr = &short_val_array[0];
1446  }
1447 
1448  n_columns = 0;
1449  for (size_type j=0; j<n_cols; ++j)
1450  {
1451  const double value = values[j];
1452  AssertIsFinite(value);
1453  if (value != 0)
1454  {
1455  col_index_ptr[n_columns] = col_indices[j];
1456  col_value_ptr[n_columns] = value;
1457  n_columns++;
1458  }
1459  }
1460 
1461  Assert(n_columns <= (TrilinosWrappers::types::int_type)n_cols, ExcInternalError());
1462  }
1463 
1464 
1465  // If the calling matrix owns the row to which we want to insert values,
1466  // we can directly call the Epetra_CrsMatrix input function, which is much
1467  // faster than the Epetra_FECrsMatrix function. We distinguish between two
1468  // cases: the first one is when the matrix is not filled (i.e., it is
1469  // possible to add new elements to the sparsity pattern), and the second
1470  // one is when the pattern is already fixed. In the former case, we add
1471  // the possibility to insert new values, and in the second we just replace
1472  // data.
1473  if (matrix->RowMap().MyGID(static_cast<TrilinosWrappers::types::int_type>(row)) == true)
1474  {
1475  if (matrix->Filled() == false)
1476  {
1477  ierr = matrix->Epetra_CrsMatrix::InsertGlobalValues(
1478  static_cast<TrilinosWrappers::types::int_type>(row),
1479  static_cast<int>(n_columns),const_cast<double *>(col_value_ptr),
1480  col_index_ptr);
1481 
1482  // When inserting elements, we do not want to create exceptions in
1483  // the case when inserting non-local data (since that's what we
1484  // want to do right now).
1485  if (ierr > 0)
1486  ierr = 0;
1487  }
1488  else
1489  ierr = matrix->Epetra_CrsMatrix::ReplaceGlobalValues(row, n_columns,
1490  col_value_ptr,
1491  col_index_ptr);
1492  }
1493  else
1494  {
1495  // When we're at off-processor data, we have to stick with the
1496  // standard Insert/ReplaceGlobalValues function. Nevertheless, the way
1497  // we call it is the fastest one (any other will lead to repeated
1498  // allocation and deallocation of memory in order to call the function
1499  // we already use, which is very inefficient if writing one element at
1500  // a time).
1501  compressed = false;
1502 
1503  if (matrix->Filled() == false)
1504  {
1505  ierr = matrix->InsertGlobalValues (1,
1506  (TrilinosWrappers::types::int_type *)&row,
1507  n_columns, col_index_ptr,
1508  &col_value_ptr,
1509  Epetra_FECrsMatrix::ROW_MAJOR);
1510  if (ierr > 0)
1511  ierr = 0;
1512  }
1513  else
1514  ierr = matrix->ReplaceGlobalValues (1,
1515  (TrilinosWrappers::types::int_type *)&row,
1516  n_columns, col_index_ptr,
1517  &col_value_ptr,
1518  Epetra_FECrsMatrix::ROW_MAJOR);
1519  // use the FECrsMatrix facilities for set even in the case when we
1520  // have explicitly set the off-processor rows because that only works
1521  // properly when adding elements, not when setting them (since we want
1522  // to only touch elements that have been set explicitly, and there is
1523  // no way on the receiving processor to identify them otherwise)
1524  }
1525 
1526  Assert (ierr <= 0, ExcAccessToNonPresentElement(row, col_index_ptr[0]));
1527  AssertThrow (ierr >= 0, ExcTrilinosError(ierr));
1528  }
1529 
1530 
1531 
1532  void
1533  SparseMatrix::add (const std::vector<size_type> &indices,
1534  const FullMatrix<TrilinosScalar> &values,
1535  const bool elide_zero_values)
1536  {
1537  Assert (indices.size() == values.m(),
1538  ExcDimensionMismatch(indices.size(), values.m()));
1539  Assert (values.m() == values.n(), ExcNotQuadratic());
1540 
1541  for (size_type i=0; i<indices.size(); ++i)
1542  add (indices[i], indices.size(), &indices[0], &values(i,0),
1543  elide_zero_values);
1544  }
1545 
1546 
1547 
1548  void
1549  SparseMatrix::add (const std::vector<size_type> &row_indices,
1550  const std::vector<size_type> &col_indices,
1551  const FullMatrix<TrilinosScalar> &values,
1552  const bool elide_zero_values)
1553  {
1554  Assert (row_indices.size() == values.m(),
1555  ExcDimensionMismatch(row_indices.size(), values.m()));
1556  Assert (col_indices.size() == values.n(),
1557  ExcDimensionMismatch(col_indices.size(), values.n()));
1558 
1559  for (size_type i=0; i<row_indices.size(); ++i)
1560  add (row_indices[i], col_indices.size(), &col_indices[0],
1561  &values(i,0), elide_zero_values);
1562  }
1563 
1564 
1565 
1566  void
1568  const std::vector<size_type> &col_indices,
1569  const std::vector<TrilinosScalar> &values,
1570  const bool elide_zero_values)
1571  {
1572  Assert (col_indices.size() == values.size(),
1573  ExcDimensionMismatch(col_indices.size(), values.size()));
1574 
1575  add (row, col_indices.size(), &col_indices[0], &values[0],
1576  elide_zero_values);
1577  }
1578 
1579 
1580 
1581  void
1583  const size_type n_cols,
1584  const size_type *col_indices,
1585  const TrilinosScalar *values,
1586  const bool elide_zero_values,
1587  const bool /*col_indices_are_sorted*/)
1588  {
1589  AssertIndexRange(row, this->m());
1590  int ierr;
1591  if (last_action == Insert)
1592  {
1593  // TODO: this could lead to a dead lock when only one processor
1594  // calls GlobalAssemble.
1595  ierr = matrix->GlobalAssemble(*column_space_map,
1596  matrix->RowMap(), false);
1597 
1598  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
1599  }
1600 
1601  last_action = Add;
1602 
1603  TrilinosWrappers::types::int_type *col_index_ptr;
1604  TrilinosScalar *col_value_ptr;
1605  TrilinosWrappers::types::int_type n_columns;
1606 
1607  double short_val_array[100];
1608  TrilinosWrappers::types::int_type short_index_array[100];
1609  std::vector<TrilinosScalar> long_val_array;
1610  std::vector<TrilinosWrappers::types::int_type> long_index_array;
1611 
1612  // If we don't elide zeros, the pointers are already available... need to
1613  // cast to non-const pointers as that is the format taken by Trilinos (but
1614  // we will not modify const data)
1615  if (elide_zero_values == false)
1616  {
1617  col_index_ptr = (TrilinosWrappers::types::int_type *)col_indices;
1618  col_value_ptr = const_cast<TrilinosScalar *>(values);
1619  n_columns = n_cols;
1620 #ifdef DEBUG
1621  for (size_type j=0; j<n_cols; ++j)
1622  AssertIsFinite(values[j]);
1623 #endif
1624  }
1625  else
1626  {
1627  // Otherwise, extract nonzero values in each row and the corresponding
1628  // index.
1629  if (n_cols > 100)
1630  {
1631  long_val_array.resize(n_cols);
1632  long_index_array.resize(n_cols);
1633  col_index_ptr = &long_index_array[0];
1634  col_value_ptr = &long_val_array[0];
1635  }
1636  else
1637  {
1638  col_index_ptr = &short_index_array[0];
1639  col_value_ptr = &short_val_array[0];
1640  }
1641 
1642  n_columns = 0;
1643  for (size_type j=0; j<n_cols; ++j)
1644  {
1645  const double value = values[j];
1646 
1647  AssertIsFinite(value);
1648  if (value != 0)
1649  {
1650  col_index_ptr[n_columns] = col_indices[j];
1651  col_value_ptr[n_columns] = value;
1652  n_columns++;
1653  }
1654  }
1655 
1656  Assert(n_columns <= (TrilinosWrappers::types::int_type)n_cols, ExcInternalError());
1657 
1658  }
1659 
1660  // If the calling processor owns the row to which we want to add values, we
1661  // can directly call the Epetra_CrsMatrix input function, which is much
1662  // faster than the Epetra_FECrsMatrix function.
1663  if (matrix->RowMap().MyGID(static_cast<TrilinosWrappers::types::int_type>(row)) == true)
1664  {
1665  ierr = matrix->Epetra_CrsMatrix::SumIntoGlobalValues(row, n_columns,
1666  col_value_ptr,
1667  col_index_ptr);
1668  }
1669  else if (nonlocal_matrix.get() != 0)
1670  {
1671  compressed = false;
1672  // this is the case when we have explicitly set the off-processor rows
1673  // and want to create a separate matrix object for them (to retain
1674  // thread-safety)
1675  Assert (nonlocal_matrix->RowMap().LID(static_cast<TrilinosWrappers::types::int_type>(row)) != -1,
1676  ExcMessage("Attempted to write into off-processor matrix row "
1677  "that has not be specified as being writable upon "
1678  "initialization"));
1679  ierr = nonlocal_matrix->SumIntoGlobalValues(row, n_columns,
1680  col_value_ptr,
1681  col_index_ptr);
1682  }
1683  else
1684  {
1685  // When we're at off-processor data, we have to stick with the
1686  // standard SumIntoGlobalValues function. Nevertheless, the way we
1687  // call it is the fastest one (any other will lead to repeated
1688  // allocation and deallocation of memory in order to call the function
1689  // we already use, which is very inefficient if writing one element at
1690  // a time).
1691  compressed = false;
1692 
1693  ierr = matrix->SumIntoGlobalValues (1,
1694  (TrilinosWrappers::types::int_type *)&row, n_columns,
1695  col_index_ptr,
1696  &col_value_ptr,
1697  Epetra_FECrsMatrix::ROW_MAJOR);
1698  }
1699 
1700 #ifdef DEBUG
1701  if (ierr > 0)
1702  {
1703  std::cout << "------------------------------------------"
1704  << std::endl;
1705  std::cout << "Got error " << ierr << " in row " << row
1706  << " of proc " << matrix->RowMap().Comm().MyPID()
1707  << " when trying to add the columns:" << std::endl;
1708  for (TrilinosWrappers::types::int_type i=0; i<n_columns; ++i)
1709  std::cout << col_index_ptr[i] << " ";
1710  std::cout << std::endl << std::endl;
1711  std::cout << "Matrix row "
1712  << (matrix->RowMap().MyGID(static_cast<TrilinosWrappers::types::int_type>(row)) == false ? "(nonlocal part)" : "")
1713  << " has the following indices:" << std::endl;
1714  std::vector<TrilinosWrappers::types::int_type> indices;
1715  const Epetra_CrsGraph *graph =
1716  (nonlocal_matrix.get() != 0 &&
1717  matrix->RowMap().MyGID(static_cast<TrilinosWrappers::types::int_type>(row)) == false) ?
1718  &nonlocal_matrix->Graph() : &matrix->Graph();
1719 
1720  indices.resize(graph->NumGlobalIndices(static_cast<TrilinosWrappers::types::int_type>(row)));
1721  int n_indices = 0;
1722  graph->ExtractGlobalRowCopy(static_cast<TrilinosWrappers::types::int_type>(row),
1723  indices.size(), n_indices, &indices[0]);
1724  AssertDimension(static_cast<unsigned int>(n_indices), indices.size());
1725 
1726  for (TrilinosWrappers::types::int_type i=0; i<n_indices; ++i)
1727  std::cout << indices[i] << " ";
1728  std::cout << std::endl << std::endl;
1729  Assert (ierr <= 0,
1730  ExcAccessToNonPresentElement(row, col_index_ptr[0]));
1731  }
1732 #endif
1733  Assert (ierr >= 0, ExcTrilinosError(ierr));
1734  }
1735 
1736 
1737 
1738  SparseMatrix &
1740  {
1741  Assert (d==0, ExcScalarAssignmentOnlyForZeroValue());
1742  compress (::VectorOperation::unknown); // TODO: why do we do this? Should we not check for is_compressed?
1743 
1744  const int ierr = matrix->PutScalar(d);
1745  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
1746  if (nonlocal_matrix.get() != 0)
1747  nonlocal_matrix->PutScalar(d);
1748 
1749  return *this;
1750  }
1751 
1752 
1753 
1754  void
1755  SparseMatrix::add (const TrilinosScalar factor,
1756  const SparseMatrix &rhs)
1757  {
1758  AssertDimension (rhs.m(), m());
1759  AssertDimension (rhs.n(), n());
1760  AssertDimension (rhs.local_range().first, local_range().first);
1761  AssertDimension (rhs.local_range().second, local_range().second);
1762  Assert(matrix->RowMap().SameAs(rhs.matrix->RowMap()),
1763  ExcMessage("Can only add matrices with same distribution of rows"));
1764  Assert(matrix->Filled() && rhs.matrix->Filled(),
1765  ExcMessage("Addition of matrices only allowed if matrices are "
1766  "filled, i.e., compress() has been called"));
1767 
1768  const std::pair<size_type, size_type>
1769  local_range = rhs.local_range();
1770  const bool same_col_map = matrix->ColMap().SameAs(rhs.matrix->ColMap());
1771 
1772  int ierr;
1773  for (size_type row=local_range.first; row < local_range.second; ++row)
1774  {
1775  const int row_local =
1776  matrix->RowMap().LID(static_cast<TrilinosWrappers::types::int_type>(row));
1777 
1778  // First get a view to the matrix columns of both matrices. Note that
1779  // the data is in local index spaces so we need to be careful not only
1780  // to compare column indices in case they are derived from the same
1781  // map.
1782  int n_entries, rhs_n_entries;
1783  TrilinosScalar *value_ptr, *rhs_value_ptr;
1784  int *index_ptr, *rhs_index_ptr;
1785  ierr = rhs.matrix->ExtractMyRowView (row_local, rhs_n_entries,
1786  rhs_value_ptr, rhs_index_ptr);
1787  (void)ierr;
1788  Assert (ierr == 0, ExcTrilinosError(ierr));
1789 
1790  ierr = matrix->ExtractMyRowView (row_local, n_entries, value_ptr,
1791  index_ptr);
1792  Assert (ierr == 0, ExcTrilinosError(ierr));
1793  bool expensive_checks = (n_entries != rhs_n_entries || !same_col_map);
1794  if (!expensive_checks)
1795  {
1796  // check if the column indices are the same. If yes, can simply
1797  // copy over the data.
1798  expensive_checks = std::memcmp(static_cast<void *>(index_ptr),
1799  static_cast<void *>(rhs_index_ptr),
1800  sizeof(int)*n_entries) != 0;
1801  if (!expensive_checks)
1802  for (int i=0; i<n_entries; ++i)
1803  value_ptr[i] += rhs_value_ptr[i] * factor;
1804  }
1805  // Now to the expensive case where we need to check all column indices
1806  // against each other (transformed into global index space) and where
1807  // we need to make sure that all entries we are about to add into the
1808  // lhs matrix actually exist
1809  if (expensive_checks)
1810  {
1811  for (int i=0; i<rhs_n_entries; ++i)
1812  {
1813  if (rhs_value_ptr[i] == 0.)
1814  continue;
1815  const TrilinosWrappers::types::int_type rhs_global_col =
1816  global_column_index(*rhs.matrix, rhs_index_ptr[i]);
1817  int local_col = matrix->ColMap().LID(rhs_global_col);
1818  int *local_index = Utilities::lower_bound(index_ptr,
1819  index_ptr+n_entries,
1820  local_col);
1821  Assert(local_index != index_ptr + n_entries &&
1822  *local_index == local_col,
1823  ExcMessage("Adding the entries from the other matrix "
1824  "failed, because the sparsity pattern "
1825  "of that matrix includes more elements than the "
1826  "calling matrix, which is not allowed."));
1827  value_ptr[local_index-index_ptr] += factor * rhs_value_ptr[i];
1828  }
1829  }
1830  }
1831  }
1832 
1833 
1834 
1835  void
1837  {
1838  // This only flips a flag that tells
1839  // Trilinos that any vmult operation
1840  // should be done with the
1841  // transpose. However, the matrix
1842  // structure is not reset.
1843  int ierr;
1844 
1845  if (!matrix->UseTranspose())
1846  {
1847  ierr = matrix->SetUseTranspose (true);
1848  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
1849  }
1850  else
1851  {
1852  ierr = matrix->SetUseTranspose (false);
1853  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
1854  }
1855  }
1856 
1857 
1858 
1859  SparseMatrix &
1860  SparseMatrix::operator *= (const TrilinosScalar a)
1861  {
1862  const int ierr = matrix->Scale (a);
1863  Assert (ierr == 0, ExcTrilinosError(ierr));
1864  (void)ierr; // removes -Wunused-variable in optimized mode
1865 
1866  return *this;
1867  }
1868 
1869 
1870 
1871  SparseMatrix &
1872  SparseMatrix::operator /= (const TrilinosScalar a)
1873  {
1874  Assert (a !=0, ExcDivideByZero());
1875 
1876  const TrilinosScalar factor = 1./a;
1877 
1878  const int ierr = matrix->Scale (factor);
1879  Assert (ierr == 0, ExcTrilinosError(ierr));
1880  (void)ierr; // removes -Wunused-variable in optimized mode
1881 
1882  return *this;
1883  }
1884 
1885 
1886 
1887  TrilinosScalar
1889  {
1890  Assert (matrix->Filled(), ExcMatrixNotCompressed());
1891  return matrix->NormOne();
1892  }
1893 
1894 
1895 
1896  TrilinosScalar
1898  {
1899  Assert (matrix->Filled(), ExcMatrixNotCompressed());
1900  return matrix->NormInf();
1901  }
1902 
1903 
1904 
1905  TrilinosScalar
1907  {
1908  Assert (matrix->Filled(), ExcMatrixNotCompressed());
1909  return matrix->NormFrobenius();
1910  }
1911 
1912 
1913 
1914  namespace internal
1915  {
1916  namespace SparseMatrix
1917  {
1918  template <typename VectorType>
1919  inline
1920  void check_vector_map_equality(const Epetra_CrsMatrix &,
1921  const VectorType &,
1922  const VectorType &)
1923  {
1924  }
1925 
1926  inline
1927  void check_vector_map_equality(const Epetra_CrsMatrix &m,
1929  const TrilinosWrappers::MPI::Vector &out)
1930  {
1931  Assert (in.vector_partitioner().SameAs(m.DomainMap()) == true,
1932  ExcMessage ("Column map of matrix does not fit with vector map!"));
1933  Assert (out.vector_partitioner().SameAs(m.RangeMap()) == true,
1934  ExcMessage ("Row map of matrix does not fit with vector map!"));
1935  (void)m;
1936  (void)in;
1937  (void)out;
1938  }
1939  }
1940  }
1941 
1942 
1943  template <typename VectorType>
1944  void
1945  SparseMatrix::vmult (VectorType &dst,
1946  const VectorType &src) const
1947  {
1948  Assert (&src != &dst, ExcSourceEqualsDestination());
1949  Assert (matrix->Filled(), ExcMatrixNotCompressed());
1950  (void)src;
1951  (void)dst;
1952 
1953  internal::SparseMatrix::check_vector_map_equality(*matrix, src, dst);
1954  const size_type dst_local_size = dst.end() - dst.begin();
1955  AssertDimension (dst_local_size, static_cast<size_type>(matrix->RangeMap().NumMyPoints()));
1956  const size_type src_local_size = src.end() - src.begin();
1957  AssertDimension (src_local_size, static_cast<size_type>(matrix->DomainMap().NumMyPoints()));
1958 
1959  Epetra_MultiVector tril_dst (View, matrix->RangeMap(), dst.begin(),
1960  dst_local_size, 1);
1961  Epetra_MultiVector tril_src (View, matrix->DomainMap(),
1962  const_cast<TrilinosScalar *>(src.begin()),
1963  src_local_size, 1);
1964 
1965  const int ierr = matrix->Multiply (false, tril_src, tril_dst);
1966  Assert (ierr == 0, ExcTrilinosError(ierr));
1967  (void)ierr; // removes -Wunused-variable in optimized mode
1968  }
1969 
1970 
1971 
1972  template <typename VectorType>
1973  void
1974  SparseMatrix::Tvmult (VectorType &dst,
1975  const VectorType &src) const
1976  {
1977  Assert (&src != &dst, ExcSourceEqualsDestination());
1978  Assert (matrix->Filled(), ExcMatrixNotCompressed());
1979 
1980  internal::SparseMatrix::check_vector_map_equality(*matrix, dst, src);
1981  const size_type dst_local_size = dst.end() - dst.begin();
1982  AssertDimension (dst_local_size, static_cast<size_type>(matrix->DomainMap().NumMyPoints()));
1983  const size_type src_local_size = src.end() - src.begin();
1984  AssertDimension (src_local_size, static_cast<size_type>(matrix->RangeMap().NumMyPoints()));
1985 
1986  Epetra_MultiVector tril_dst (View, matrix->DomainMap(), dst.begin(),
1987  dst_local_size, 1);
1988  Epetra_MultiVector tril_src (View, matrix->RangeMap(),
1989  const_cast<double *>(src.begin()),
1990  src_local_size, 1);
1991 
1992  const int ierr = matrix->Multiply (true, tril_src, tril_dst);
1993  Assert (ierr == 0, ExcTrilinosError(ierr));
1994  (void)ierr; // removes -Wunused-variable in optimized mode
1995  }
1996 
1997 
1998 
1999  template <typename VectorType>
2000  void
2001  SparseMatrix::vmult_add (VectorType &dst,
2002  const VectorType &src) const
2003  {
2004  Assert (&src != &dst, ExcSourceEqualsDestination());
2005 
2006  // Reinit a temporary vector with fast argument set, which does not
2007  // overwrite the content (to save time). However, the
2008  // TrilinosWrappers::Vector classes do not support this, so create a
2009  // deal.II local vector that has this fast setting. It will be accepted in
2010  // vmult because it only checks the local size.
2011  ::Vector<TrilinosScalar> temp_vector;
2012  temp_vector.reinit(dst.end()-dst.begin(), true);
2013  ::VectorView<TrilinosScalar> src_view(src.end()-src.begin(), src.begin());
2014  ::VectorView<TrilinosScalar> dst_view(dst.end()-dst.begin(), dst.begin());
2015  vmult (temp_vector, static_cast<const ::Vector<TrilinosScalar>&>(src_view));
2016  if (dst_view.size() > 0)
2017  dst_view += temp_vector;
2018  }
2019 
2020 
2021 
2022  template <typename VectorType>
2023  void
2024  SparseMatrix::Tvmult_add (VectorType &dst,
2025  const VectorType &src) const
2026  {
2027  Assert (&src != &dst, ExcSourceEqualsDestination());
2028 
2029  // Reinit a temporary vector with fast argument set, which does not
2030  // overwrite the content (to save time). However, the
2031  // TrilinosWrappers::Vector classes do not support this, so create a
2032  // deal.II local vector that has this fast setting. It will be accepted in
2033  // vmult because it only checks the local size.
2034  ::Vector<TrilinosScalar> temp_vector;
2035  temp_vector.reinit(dst.end()-dst.begin(), true);
2036  ::VectorView<TrilinosScalar> src_view(src.end()-src.begin(), src.begin());
2037  ::VectorView<TrilinosScalar> dst_view(dst.end()-dst.begin(), dst.begin());
2038  Tvmult (temp_vector, static_cast<const ::Vector<TrilinosScalar>&>(src_view));
2039  if (dst_view.size() > 0)
2040  dst_view += temp_vector;
2041  }
2042 
2043 
2044 
2045  TrilinosScalar
2047  {
2048  Assert (matrix->RowMap().SameAs(matrix->DomainMap()),
2049  ExcNotQuadratic());
2050 
2051  VectorBase temp_vector;
2052  temp_vector.reinit(v, true);
2053 
2054  vmult (temp_vector, v);
2055  return temp_vector*v;
2056  }
2057 
2058 
2059 
2060  TrilinosScalar
2062  const VectorBase &v) const
2063  {
2064  Assert (matrix->RowMap().SameAs(matrix->DomainMap()),
2065  ExcNotQuadratic());
2066 
2067  VectorBase temp_vector;
2068  temp_vector.reinit(v, true);
2069 
2070  vmult (temp_vector, v);
2071  return u*temp_vector;
2072  }
2073 
2074 
2075 
2076  TrilinosScalar
2078  const VectorBase &x,
2079  const VectorBase &b) const
2080  {
2081  vmult (dst, x);
2082  dst -= b;
2083  dst *= -1.;
2084 
2085  return dst.l2_norm();
2086  }
2087 
2088 
2089 
2090  namespace internals
2091  {
2092  typedef ::types::global_dof_index size_type;
2093 
2094  void perform_mmult (const SparseMatrix &inputleft,
2095  const SparseMatrix &inputright,
2096  SparseMatrix &result,
2097  const VectorBase &V,
2098  const bool transpose_left)
2099  {
2100 #ifdef DEAL_II_WITH_64BIT_INDICES
2101  Assert(false,ExcNotImplemented())
2102 #endif
2103  const bool use_vector = (V.size() == inputright.m() ? true : false);
2104  if (transpose_left == false)
2105  {
2106  Assert (inputleft.n() == inputright.m(),
2107  ExcDimensionMismatch(inputleft.n(), inputright.m()));
2108  Assert (inputleft.domain_partitioner().SameAs(inputright.range_partitioner()),
2109  ExcMessage ("Parallel partitioning of A and B does not fit."));
2110  }
2111  else
2112  {
2113  Assert (inputleft.m() == inputright.m(),
2114  ExcDimensionMismatch(inputleft.m(), inputright.m()));
2115  Assert (inputleft.range_partitioner().SameAs(inputright.range_partitioner()),
2116  ExcMessage ("Parallel partitioning of A and B does not fit."));
2117  }
2118 
2119  result.clear();
2120 
2121  // create a suitable operator B: in case
2122  // we do not use a vector, all we need to
2123  // do is to set the pointer. Otherwise,
2124  // we insert the data from B, but
2125  // multiply each row with the respective
2126  // vector element.
2127  Teuchos::RCP<Epetra_CrsMatrix> mod_B;
2128  if (use_vector == false)
2129  {
2130  mod_B = Teuchos::rcp(const_cast<Epetra_CrsMatrix *>
2131  (&inputright.trilinos_matrix()),
2132  false);
2133  }
2134  else
2135  {
2136  mod_B = Teuchos::rcp(new Epetra_CrsMatrix
2137  (Copy, inputright.trilinos_sparsity_pattern()),
2138  true);
2139  mod_B->FillComplete(inputright.domain_partitioner(),
2140  inputright.range_partitioner());
2141  Assert (inputright.local_range() == V.local_range(),
2142  ExcMessage ("Parallel distribution of matrix B and vector V "
2143  "does not match."));
2144 
2145  const int local_N = inputright.local_size();
2146  for (int i=0; i<local_N; ++i)
2147  {
2148  int N_entries = -1;
2149  double *new_data, *B_data;
2150  mod_B->ExtractMyRowView (i, N_entries, new_data);
2151  inputright.trilinos_matrix().ExtractMyRowView (i, N_entries, B_data);
2152  double value = V.trilinos_vector()[0][i];
2153  for (TrilinosWrappers::types::int_type j=0; j<N_entries; ++j)
2154  new_data[j] = value * B_data[j];
2155  }
2156  }
2157 
2158  // use ML built-in method for performing
2159  // the matrix-matrix product.
2160  // create ML operators on top of the
2161  // Epetra matrices. if we use a
2162  // transposed matrix, let ML know it
2163  ML_Comm *comm;
2164  ML_Comm_Create(&comm);
2165 #ifdef ML_MPI
2166  const Epetra_MpiComm *epcomm = dynamic_cast<const Epetra_MpiComm *>(&(inputleft.trilinos_matrix().Comm()));
2167  // Get the MPI communicator, as it may not be MPI_COMM_W0RLD, and update the ML comm object
2168  if (epcomm) ML_Comm_Set_UsrComm(comm,epcomm->Comm());
2169 #endif
2170  ML_Operator *A_ = ML_Operator_Create(comm);
2171  ML_Operator *B_ = ML_Operator_Create(comm);
2172  ML_Operator *C_ = ML_Operator_Create(comm);
2173  SparseMatrix transposed_mat;
2174 
2175  if (transpose_left == false)
2176  ML_Operator_WrapEpetraCrsMatrix
2177  (const_cast<Epetra_CrsMatrix *>(&inputleft.trilinos_matrix()),A_,
2178  false);
2179  else
2180  {
2181  // create transposed matrix
2182  SparsityPattern sparsity_transposed (inputleft.domain_partitioner(),
2183  inputleft.range_partitioner());
2184  Assert (inputleft.domain_partitioner().LinearMap() == true,
2185  ExcMessage("Matrix must be partitioned contiguously between procs."));
2186  for (unsigned int i=0; i<inputleft.local_size(); ++i)
2187  {
2188  int num_entries, * indices;
2189  inputleft.trilinos_sparsity_pattern().ExtractMyRowView(i, num_entries,
2190  indices);
2191  Assert (num_entries >= 0, ExcInternalError());
2192 #ifndef DEAL_II_WITH_64BIT_INDICES
2193  const size_type GID = inputleft.trilinos_matrix().RowMap().GID(i);
2194  for (TrilinosWrappers::types::int_type j=0; j<num_entries; ++j)
2195  sparsity_transposed.add (inputleft.trilinos_matrix().ColMap().GID(indices[j]),
2196  GID);
2197 #else
2198  const size_type GID = inputleft.trilinos_matrix().RowMap().GID64(i);
2199  for (TrilinosWrappers::types::int_type j=0; j<num_entries; ++j)
2200  sparsity_transposed.add (inputleft.trilinos_matrix().ColMap().GID64(indices[j]),
2201  GID);
2202 #endif
2203  }
2204 
2205  sparsity_transposed.compress();
2206  transposed_mat.reinit (sparsity_transposed);
2207  for (unsigned int i=0; i<inputleft.local_size(); ++i)
2208  {
2209  int num_entries, * indices;
2210  double *values;
2211  inputleft.trilinos_matrix().ExtractMyRowView(i, num_entries,
2212  values, indices);
2213  Assert (num_entries >= 0, ExcInternalError());
2214 #ifndef DEAL_II_WITH_64BIT_INDICES
2215  const size_type GID = inputleft.trilinos_matrix().RowMap().GID(i);
2216  for (TrilinosWrappers::types::int_type j=0; j<num_entries; ++j)
2217  transposed_mat.set (inputleft.trilinos_matrix().ColMap().GID(indices[j]),
2218  GID, values[j]);
2219 #else
2220  const size_type GID = inputleft.trilinos_matrix().RowMap().GID64(i);
2221  for (TrilinosWrappers::types::int_type j=0; j<num_entries; ++j)
2222  transposed_mat.set (inputleft.trilinos_matrix().ColMap().GID64(indices[j]),
2223  GID, values[j]);
2224 #endif
2225  }
2226  transposed_mat.compress(VectorOperation::insert);
2227  ML_Operator_WrapEpetraCrsMatrix
2228  (const_cast<Epetra_CrsMatrix *>(&transposed_mat.trilinos_matrix()),
2229  A_,false);
2230  }
2231  ML_Operator_WrapEpetraCrsMatrix(mod_B.get(),B_,false);
2232 
2233  // We implement the multiplication by
2234  // hand in a similar way as is done in
2235  // ml/src/Operator/ml_rap.c for a triple
2236  // matrix product. This means that the
2237  // code is very similar to the one found
2238  // in ml/src/Operator/ml_rap.c
2239 
2240  // import data if necessary
2241  ML_Operator *Btmp, *Ctmp, *Ctmp2, *tptr;
2242  ML_CommInfoOP *getrow_comm;
2243  int max_per_proc;
2244  TrilinosWrappers::types::int_type N_input_vector = B_->invec_leng;
2245  getrow_comm = B_->getrow->pre_comm;
2246  if ( getrow_comm != NULL)
2247  for (TrilinosWrappers::types::int_type i = 0; i < getrow_comm->N_neighbors; i++)
2248  for (TrilinosWrappers::types::int_type j = 0; j < getrow_comm->neighbors[i].N_send; j++)
2249  AssertThrow (getrow_comm->neighbors[i].send_list[j] < N_input_vector,
2250  ExcInternalError());
2251 
2252  ML_create_unique_col_id(N_input_vector, &(B_->getrow->loc_glob_map),
2253  getrow_comm, &max_per_proc, B_->comm);
2254  B_->getrow->use_loc_glob_map = ML_YES;
2255  if (A_->getrow->pre_comm != NULL)
2256  ML_exchange_rows( B_, &Btmp, A_->getrow->pre_comm);
2257  else Btmp = B_;
2258 
2259  // perform matrix-matrix product
2260  ML_matmat_mult(A_, Btmp , &Ctmp);
2261 
2262  // release temporary structures we needed
2263  // for multiplication
2264  ML_free(B_->getrow->loc_glob_map);
2265  B_->getrow->loc_glob_map = NULL;
2266  B_->getrow->use_loc_glob_map = ML_NO;
2267  if (A_->getrow->pre_comm != NULL)
2268  {
2269  tptr = Btmp;
2270  while ( (tptr!= NULL) && (tptr->sub_matrix != B_))
2271  tptr = tptr->sub_matrix;
2272  if (tptr != NULL) tptr->sub_matrix = NULL;
2273  ML_RECUR_CSR_MSRdata_Destroy(Btmp);
2274  ML_Operator_Destroy(&Btmp);
2275  }
2276 
2277  // make correct data structures
2278  if (A_->getrow->post_comm != NULL)
2279  ML_exchange_rows(Ctmp, &Ctmp2, A_->getrow->post_comm);
2280  else
2281  Ctmp2 = Ctmp;
2282 
2283  ML_back_to_csrlocal(Ctmp2, C_, max_per_proc);
2284 
2285  ML_RECUR_CSR_MSRdata_Destroy (Ctmp);
2286  ML_Operator_Destroy (&Ctmp);
2287 
2288  if (A_->getrow->post_comm != NULL)
2289  {
2290  ML_RECUR_CSR_MSRdata_Destroy(Ctmp2);
2291  ML_Operator_Destroy (&Ctmp2);
2292  }
2293 
2294  // create an Epetra matrix from the ML
2295  // matrix that we got as a result.
2296  Epetra_CrsMatrix *C_mat;
2297  ML_Operator2EpetraCrsMatrix(C_, C_mat);
2298  C_mat->FillComplete();
2299  C_mat->OptimizeStorage();
2300  result.reinit (*C_mat);
2301 
2302  // destroy allocated memory
2303  delete C_mat;
2304  ML_Operator_Destroy (&A_);
2305  ML_Operator_Destroy (&B_);
2306  ML_Operator_Destroy (&C_);
2307  ML_Comm_Destroy (&comm);
2308  }
2309  }
2310 
2311 
2312  void
2314  const SparseMatrix &B,
2315  const VectorBase &V) const
2316  {
2317 #ifdef DEAL_II_WITH_64BIT_INDICES
2318  Assert(false,ExcNotImplemented())
2319 #endif
2320  internals::perform_mmult (*this, B, C, V, false);
2321  }
2322 
2323 
2324 
2325  void
2327  const SparseMatrix &B,
2328  const VectorBase &V) const
2329  {
2330 #ifdef DEAL_II_WITH_64BIT_INDICES
2331  Assert(false,ExcNotImplemented())
2332 #endif
2333  internals::perform_mmult (*this, B, C, V, true);
2334  }
2335 
2336 
2337 
2338  void
2340  {
2341  Assert (false, ExcNotImplemented());
2342  }
2343 
2344 
2345 
2346  // As of now, no particularly neat
2347  // ouput is generated in case of
2348  // multiple processors.
2349  void
2350  SparseMatrix::print (std::ostream &out,
2351  const bool print_detailed_trilinos_information) const
2352  {
2353  if (print_detailed_trilinos_information == true)
2354  out << *matrix;
2355  else
2356  {
2357  double *values;
2358  int *indices;
2359  int num_entries;
2360 
2361  for (int i=0; i<matrix->NumMyRows(); ++i)
2362  {
2363  matrix->ExtractMyRowView (i, num_entries, values, indices);
2364  for (TrilinosWrappers::types::int_type j=0; j<num_entries; ++j)
2365  out << "(" << global_row_index(*matrix,i) << ","
2366  << global_column_index(*matrix,indices[j]) << ") "
2367  << values[j] << std::endl;
2368  }
2369  }
2370 
2371  AssertThrow (out, ExcIO());
2372  }
2373 
2374 
2375 
2378  {
2379  size_type static_memory = sizeof(this) + sizeof (*matrix)
2380  + sizeof(*matrix->Graph().DataPtr());
2381  return ((sizeof(TrilinosScalar)+sizeof(TrilinosWrappers::types::int_type))*
2382  matrix->NumMyNonzeros() + sizeof(int)*local_size() + static_memory);
2383  }
2384 
2385 
2386 
2387  const Epetra_Map &
2389  {
2390  return matrix->DomainMap();
2391  }
2392 
2393 
2394 
2395  const Epetra_Map &
2397  {
2398  return matrix->RangeMap();
2399  }
2400 
2401 
2402 
2403  const Epetra_Map &
2405  {
2406  return matrix->RowMap();
2407  }
2408 
2409 
2410 
2411  const Epetra_Map &
2413  {
2414  return matrix->ColMap();
2415  }
2416 
2417 
2418 
2420  {
2421 
2422 #ifdef DEAL_II_WITH_MPI
2423 
2424  const Epetra_MpiComm *mpi_comm
2425  = dynamic_cast<const Epetra_MpiComm *>(&matrix->RangeMap().Comm());
2426  return mpi_comm->Comm();
2427 #else
2428 
2429  return MPI_COMM_SELF;
2430 
2431 #endif
2432 
2433  }
2434 }
2435 
2436 
2437 
2438 // explicit instantiations
2439 #include "trilinos_sparse_matrix.inst"
2440 
2441 
2442 // TODO: put these instantiations into generic file
2443 namespace TrilinosWrappers
2444 {
2445  template void
2446  SparseMatrix::reinit (const ::SparsityPattern &);
2447 
2448  template void
2450 
2451  template void
2452  SparseMatrix::reinit (const Epetra_Map &,
2453  const ::SparsityPattern &,
2454  const bool);
2455  template void
2456  SparseMatrix::reinit (const Epetra_Map &,
2457  const DynamicSparsityPattern &,
2458  const bool);
2459  template void
2460  SparseMatrix::reinit (const Epetra_Map &,
2461  const Epetra_Map &,
2462  const ::SparsityPattern &,
2463  const bool);
2464  template void
2465  SparseMatrix::reinit (const Epetra_Map &,
2466  const Epetra_Map &,
2467  const DynamicSparsityPattern &,
2468  const bool);
2469  template void
2470  SparseMatrix::reinit (const IndexSet &,
2471  const IndexSet &,
2472  const ::SparsityPattern &,
2473  const MPI_Comm &,
2474  const bool);
2475  template void
2476  SparseMatrix::reinit (const IndexSet &,
2477  const IndexSet &,
2478  const DynamicSparsityPattern &,
2479  const MPI_Comm &,
2480  const bool);
2481 
2482  template void
2484  const VectorBase &) const;
2485  template void
2487  const Vector &) const;
2488  template void
2490  const MPI::Vector &) const;
2491  template void
2493  const ::Vector<double> &) const;
2494  template void
2496  const ::parallel::distributed::Vector<double> &) const;
2497  template void
2499  const VectorBase &) const;
2500  template void
2502  const Vector &) const;
2503  template void
2505  const MPI::Vector &) const;
2506  template void
2508  const ::Vector<double> &) const;
2509  template void
2511  const ::parallel::distributed::Vector<double> &) const;
2512  template void
2514  const VectorBase &) const;
2515  template void
2517  const Vector &) const;
2518  template void
2520  const MPI::Vector &) const;
2521  template void
2523  const ::Vector<double> &) const;
2524  template void
2526  const ::parallel::distributed::Vector<double> &) const;
2527  template void
2529  const VectorBase &) const;
2530  template void
2532  const Vector &) const;
2533  template void
2535  const MPI::Vector &) const;
2536  template void
2538  const ::Vector<double> &) const;
2539  template void
2541  const ::parallel::distributed::Vector<double> &) const;
2542 }
2543 
2544 DEAL_II_NAMESPACE_CLOSE
2545 
2546 #endif // DEAL_II_WITH_TRILINOS
std_cxx11::shared_ptr< std::vector< TrilinosScalar > > value_cache
TrilinosScalar diag_element(const size_type i) const
Iterator lower_bound(Iterator first, Iterator last, const T &val)
Definition: utilities.h:604
TrilinosScalar el(const size_type i, const size_type j) const
size_type m() const
std_cxx11::shared_ptr< std::vector< size_type > > colnum_cache
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1052
void clear_row(const size_type row, const TrilinosScalar new_diag_value=0)
size_type column_number(const size_type row, const size_type index) const
void vmult_add(VectorType &dst, const VectorType &src) const
const Epetra_FECrsGraph & trilinos_sparsity_pattern() const
::ExceptionBase & ExcMessage(std::string arg1)
::types::global_dof_index size_type
#define AssertIndexRange(index, range)
Definition: exceptions.h:1081
TrilinosScalar residual(VectorBase &dst, const VectorBase &x, const VectorBase &b) const
std::pair< size_type, size_type > local_range() const
std_cxx11::shared_ptr< Epetra_CrsGraph > nonlocal_graph
TrilinosScalar operator()(const size_type i, const size_type j) const
void reinit(const VectorBase &v, const bool omit_zeroing_entries=false)
#define AssertThrow(cond, exc)
Definition: exceptions.h:358
const Epetra_Comm & comm_self()
Definition: utilities.cc:733
void clear_rows(const std::vector< size_type > &rows, const TrilinosScalar new_diag_value=0)
void Tvmult(VectorType &dst, const VectorType &src) const
SparseMatrix & operator=(const double d)
std_cxx11::shared_ptr< Epetra_FECrsMatrix > matrix
Epetra_Map make_trilinos_map(const MPI_Comm &communicator=MPI_COMM_WORLD, const bool overlapping=false) const
Definition: index_set.cc:484
SparseMatrix & operator*=(const TrilinosScalar factor)
size_type size() const
Definition: index_set.h:1247
const Epetra_Map & col_partitioner() const DEAL_II_DEPRECATED
std_cxx11::shared_ptr< Epetra_CrsMatrix > nonlocal_matrix
size_type n() const
Definition: types.h:30
void add(const size_type i, const size_type j, const TrilinosScalar value)
#define Assert(cond, exc)
Definition: exceptions.h:294
std_cxx11::shared_ptr< Epetra_Map > column_space_map
void Tmmult(SparseMatrix &C, const SparseMatrix &B, const VectorBase &V=VectorBase()) const
const Epetra_Map & vector_partitioner() const
const IndexSet & row_index_set() const
void print(std::ostream &out, const bool write_extended_trilinos_info=false) const
void mmult(SparseMatrix &C, const SparseMatrix &B, const VectorBase &V=VectorBase()) const
const Epetra_CrsMatrix & trilinos_matrix() const
void vmult(VectorType &dst, const VectorType &src) const
size_type n_nonzero_elements() const
const Epetra_Map & domain_partitioner() const DEAL_II_DEPRECATED
unsigned int row_length(const size_type row) const
void set_size(const size_type size)
Definition: index_set.h:1234
void Tvmult_add(VectorType &dst, const VectorType &src) const
Definition: mpi.h:55
const Epetra_Map & row_partitioner() const DEAL_II_DEPRECATED
void reinit(const SparsityPatternType &sparsity_pattern)
TrilinosScalar matrix_norm_square(const VectorBase &v) const
size_type row_length(const size_type row) const
void copy_from(const SparseMatrix &source)
SparseMatrix & operator/=(const TrilinosScalar factor)
std::pair< size_type, size_type > local_range() const
const Epetra_MultiVector & trilinos_vector() const
TrilinosScalar matrix_scalar_product(const VectorBase &u, const VectorBase &v) const
const Epetra_CrsGraph & trilinos_sparsity_pattern() const
const Epetra_Map & domain_partitioner() const DEAL_II_DEPRECATED
bool is_element(const size_type index) const
Definition: index_set.h:1317
void reinit(const size_type n, const bool omit_zeroing_entries=false)
void compress(::VectorOperation::values operation)
void reinit(const size_type m, const size_type n, const size_type n_entries_per_row=0)
void set(const size_type i, const size_type j, const TrilinosScalar value)
#define AssertIsFinite(number)
Definition: exceptions.h:1096
const Epetra_Map & range_partitioner() const DEAL_II_DEPRECATED
unsigned int local_size() const
real_type l2_norm() const
std_cxx11::shared_ptr< Epetra_Export > nonlocal_matrix_exporter