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_sparsity_pattern.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/mpi.h>
22 # include <deal.II/base/utilities.h>
23 
26 
27 # include <Epetra_Export.h>
28 
30 
31 namespace TrilinosWrappers
32 {
33  namespace SparsityPatternIterators
34  {
35  void
37  {
38  // if we are asked to visit the past-the-end line, then simply
39  // release all our caches and go on with life
40  if (this->a_row == sparsity_pattern->n_rows())
41  {
42  colnum_cache.reset();
43  return;
44  }
45 
46  // otherwise first flush Trilinos caches if necessary
49 
50  colnum_cache = std::make_shared<std::vector<size_type>>(
52 
53  if (colnum_cache->size() > 0)
54  {
55  // get a representation of the present row
56  int ncols;
57  const int ierr = sparsity_pattern->graph->ExtractGlobalRowCopy(
58  this->a_row,
59  colnum_cache->size(),
60  ncols,
61  reinterpret_cast<TrilinosWrappers::types::int_type *>(
62  const_cast<size_type *>(colnum_cache->data())));
63  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
64  AssertThrow(static_cast<std::vector<size_type>::size_type>(ncols) ==
65  colnum_cache->size(),
67  }
68  }
69  } // namespace SparsityPatternIterators
70 
71 
72  // The constructor is actually the
73  // only point where we have to check
74  // whether we build a serial or a
75  // parallel Trilinos matrix.
76  // Actually, it does not even matter
77  // how many threads there are, but
78  // only if we use an MPI compiler or
79  // a standard compiler. So, even one
80  // thread on a configuration with
81  // MPI will still get a parallel
82  // interface.
84  {
86  std_cxx14::make_unique<Epetra_Map>(TrilinosWrappers::types::int_type(0),
89  graph = std_cxx14::make_unique<Epetra_FECrsGraph>(View,
92  0);
93  graph->FillComplete();
94  }
95 
96 
97 
99  const size_type n,
100  const size_type n_entries_per_row)
101  {
102  reinit(m, n, n_entries_per_row);
103  }
104 
105 
106 
108  const size_type m,
109  const size_type n,
110  const std::vector<size_type> &n_entries_per_row)
111  {
112  reinit(m, n, n_entries_per_row);
113  }
114 
115 
116 
118  : Subscriptor(std::move(other))
119  , column_space_map(std::move(other.column_space_map))
120  , graph(std::move(other.graph))
121  , nonlocal_graph(std::move(other.nonlocal_graph))
122  {}
123 
124 
125 
126  // Copy function only works if the
127  // sparsity pattern is empty.
129  : Subscriptor()
130  , column_space_map(new Epetra_Map(TrilinosWrappers::types::int_type(0),
132  Utilities::Trilinos::comm_self()))
133  , graph(
134  new Epetra_FECrsGraph(View, *column_space_map, *column_space_map, 0))
135  {
136  (void)input_sparsity;
137  Assert(input_sparsity.n_rows() == 0,
138  ExcMessage(
139  "Copy constructor only works for empty sparsity patterns."));
140  }
141 
142 
143 
144  SparsityPattern::SparsityPattern(const IndexSet &parallel_partitioning,
145  const MPI_Comm &communicator,
146  const size_type n_entries_per_row)
147  {
148  reinit(parallel_partitioning,
149  parallel_partitioning,
150  communicator,
151  n_entries_per_row);
152  }
153 
154 
155 
157  const IndexSet & parallel_partitioning,
158  const MPI_Comm & communicator,
159  const std::vector<size_type> &n_entries_per_row)
160  {
161  reinit(parallel_partitioning,
162  parallel_partitioning,
163  communicator,
164  n_entries_per_row);
165  }
166 
167 
168 
169  SparsityPattern::SparsityPattern(const IndexSet &row_parallel_partitioning,
170  const IndexSet &col_parallel_partitioning,
171  const MPI_Comm &communicator,
172  const size_type n_entries_per_row)
173  {
174  reinit(row_parallel_partitioning,
175  col_parallel_partitioning,
176  communicator,
177  n_entries_per_row);
178  }
179 
180 
181 
183  const IndexSet & row_parallel_partitioning,
184  const IndexSet & col_parallel_partitioning,
185  const MPI_Comm & communicator,
186  const std::vector<size_type> &n_entries_per_row)
187  {
188  reinit(row_parallel_partitioning,
189  col_parallel_partitioning,
190  communicator,
191  n_entries_per_row);
192  }
193 
194 
195 
196  SparsityPattern::SparsityPattern(const IndexSet &row_parallel_partitioning,
197  const IndexSet &col_parallel_partitioning,
198  const IndexSet &writable_rows,
199  const MPI_Comm &communicator,
200  const size_type n_max_entries_per_row)
201  {
202  reinit(row_parallel_partitioning,
203  col_parallel_partitioning,
204  writable_rows,
205  communicator,
206  n_max_entries_per_row);
207  }
208 
209 
210 
211  void
213  const size_type n,
214  const size_type n_entries_per_row)
215  {
218  MPI_COMM_SELF,
219  n_entries_per_row);
220  }
221 
222 
223 
224  void
226  const size_type n,
227  const std::vector<size_type> &n_entries_per_row)
228  {
231  MPI_COMM_SELF,
232  n_entries_per_row);
233  }
234 
235 
236 
237  namespace
238  {
240 
241  void
242  reinit_sp(const Epetra_Map & row_map,
243  const Epetra_Map & col_map,
244  const size_type n_entries_per_row,
245  std::unique_ptr<Epetra_Map> & column_space_map,
246  std::unique_ptr<Epetra_FECrsGraph> &graph,
247  std::unique_ptr<Epetra_CrsGraph> & nonlocal_graph)
248  {
249  Assert(row_map.IsOneToOne(),
250  ExcMessage("Row map must be 1-to-1, i.e., no overlap between "
251  "the maps of different processors."));
252  Assert(col_map.IsOneToOne(),
253  ExcMessage("Column map must be 1-to-1, i.e., no overlap between "
254  "the maps of different processors."));
255 
256  nonlocal_graph.reset();
257  graph.reset();
258  column_space_map = std_cxx14::make_unique<Epetra_Map>(col_map);
259 
260  // for more than one processor, need to specify only row map first and
261  // let the matrix entries decide about the column map (which says which
262  // columns are present in the matrix, not to be confused with the
263  // col_map that tells how the domain dofs of the matrix will be
264  // distributed). for only one processor, we can directly assign the
265  // columns as well. If we use a recent Trilinos version, we can also
266  // require building a non-local graph which gives us thread-safe
267  // initialization.
268  if (row_map.Comm().NumProc() > 1)
269  graph = std_cxx14::make_unique<Epetra_FECrsGraph>(
270  Copy, row_map, n_entries_per_row, false
271  // TODO: Check which new Trilinos version supports this...
272  // Remember to change tests/trilinos/assemble_matrix_parallel_07, too.
273  //#if DEAL_II_TRILINOS_VERSION_GTE(11,14,0)
274  // , true
275  //#endif
276  );
277  else
278  graph = std_cxx14::make_unique<Epetra_FECrsGraph>(
279  Copy, row_map, col_map, n_entries_per_row, false);
280  }
281 
282 
283 
284  void
285  reinit_sp(const Epetra_Map & row_map,
286  const Epetra_Map & col_map,
287  const std::vector<size_type> & n_entries_per_row,
288  std::unique_ptr<Epetra_Map> & column_space_map,
289  std::unique_ptr<Epetra_FECrsGraph> &graph,
290  std::unique_ptr<Epetra_CrsGraph> & nonlocal_graph)
291  {
292  Assert(row_map.IsOneToOne(),
293  ExcMessage("Row map must be 1-to-1, i.e., no overlap between "
294  "the maps of different processors."));
295  Assert(col_map.IsOneToOne(),
296  ExcMessage("Column map must be 1-to-1, i.e., no overlap between "
297  "the maps of different processors."));
298 
299  // release memory before reallocation
300  nonlocal_graph.reset();
301  graph.reset();
302  AssertDimension(n_entries_per_row.size(),
304 
305  column_space_map = std_cxx14::make_unique<Epetra_Map>(col_map);
306  std::vector<int> local_entries_per_row(
309  for (unsigned int i = 0; i < local_entries_per_row.size(); ++i)
310  local_entries_per_row[i] =
311  n_entries_per_row[TrilinosWrappers::min_my_gid(row_map) + i];
312 
313  if (row_map.Comm().NumProc() > 1)
314  graph = std_cxx14::make_unique<Epetra_FECrsGraph>(
315  Copy, row_map, local_entries_per_row.data(), false
316  // TODO: Check which new Trilinos version supports this...
317  // Remember to change tests/trilinos/assemble_matrix_parallel_07, too.
318  //#if DEAL_II_TRILINOS_VERSION_GTE(11,14,0)
319  // , true
320  //#endif
321  );
322  else
323  graph = std_cxx14::make_unique<Epetra_FECrsGraph>(
324  Copy, row_map, col_map, local_entries_per_row.data(), false);
325  }
326 
327 
328 
329  template <typename SparsityPatternType>
330  void
331  reinit_sp(const Epetra_Map & row_map,
332  const Epetra_Map & col_map,
333  const SparsityPatternType & sp,
334  const bool exchange_data,
335  std::unique_ptr<Epetra_Map> & column_space_map,
336  std::unique_ptr<Epetra_FECrsGraph> &graph,
337  std::unique_ptr<Epetra_CrsGraph> & nonlocal_graph)
338  {
339  nonlocal_graph.reset();
340  graph.reset();
341 
342  AssertDimension(sp.n_rows(),
344  AssertDimension(sp.n_cols(),
346 
347  column_space_map = std_cxx14::make_unique<Epetra_Map>(col_map);
348 
349  Assert(row_map.LinearMap() == true,
350  ExcMessage(
351  "This function only works if the row map is contiguous."));
352 
353  const size_type first_row = TrilinosWrappers::min_my_gid(row_map),
354  last_row = TrilinosWrappers::max_my_gid(row_map) + 1;
355  std::vector<int> n_entries_per_row(last_row - first_row);
356 
357  // Trilinos wants the row length as an int this is hopefully never going
358  // to be a problem.
359  for (size_type row = first_row; row < last_row; ++row)
360  n_entries_per_row[row - first_row] =
361  static_cast<int>(sp.row_length(row));
362 
363  if (row_map.Comm().NumProc() > 1)
364  graph = std_cxx14::make_unique<Epetra_FECrsGraph>(
365  Copy, row_map, n_entries_per_row.data(), false);
366  else
367  graph = std_cxx14::make_unique<Epetra_FECrsGraph>(
368  Copy, row_map, col_map, n_entries_per_row.data(), false);
369 
370  AssertDimension(sp.n_rows(), n_global_rows(*graph));
371 
372  std::vector<TrilinosWrappers::types::int_type> row_indices;
373 
374  // Include possibility to exchange data since DynamicSparsityPattern is
375  // able to do so
376  if (exchange_data == false)
377  for (size_type row = first_row; row < last_row; ++row)
378  {
379  const TrilinosWrappers::types::int_type row_length =
380  sp.row_length(row);
381  if (row_length == 0)
382  continue;
383 
384  row_indices.resize(row_length, -1);
385  {
386  typename SparsityPatternType::iterator p = sp.begin(row);
387  // avoid incrementing p over the end of the current row because
388  // it is slow for DynamicSparsityPattern in parallel
389  for (int col = 0; col < row_length;)
390  {
391  row_indices[col++] = p->column();
392  if (col < row_length)
393  ++p;
394  }
395  }
396  graph->Epetra_CrsGraph::InsertGlobalIndices(row,
397  row_length,
398  row_indices.data());
399  }
400  else
401  for (size_type row = 0; row < sp.n_rows(); ++row)
402  {
403  const TrilinosWrappers::types::int_type row_length =
404  sp.row_length(row);
405  if (row_length == 0)
406  continue;
407 
408  row_indices.resize(row_length, -1);
409  {
410  typename SparsityPatternType::iterator p = sp.begin(row);
411  // avoid incrementing p over the end of the current row because
412  // it is slow for DynamicSparsityPattern in parallel
413  for (int col = 0; col < row_length;)
414  {
415  row_indices[col++] = p->column();
416  if (col < row_length)
417  ++p;
418  }
419  }
420  const TrilinosWrappers::types::int_type trilinos_row = row;
421  graph->InsertGlobalIndices(1,
422  &trilinos_row,
423  row_length,
424  row_indices.data());
425  }
426 
427  // TODO A dynamic_cast fails here, this is suspicious.
428  const auto &range_map =
429  static_cast<const Epetra_Map &>(graph->RangeMap()); // NOLINT
430  int ierr = graph->GlobalAssemble(*column_space_map, range_map, true);
431  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
432 
433  ierr = graph->OptimizeStorage();
434  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
435  }
436  } // namespace
437 
438 
439 
440  void
441  SparsityPattern::reinit(const IndexSet &parallel_partitioning,
442  const MPI_Comm &communicator,
443  const size_type n_entries_per_row)
444  {
445  Epetra_Map map =
446  parallel_partitioning.make_trilinos_map(communicator, false);
447  reinit_sp(
448  map, map, n_entries_per_row, column_space_map, graph, nonlocal_graph);
449  }
450 
451 
452 
453  void
454  SparsityPattern::reinit(const IndexSet & parallel_partitioning,
455  const MPI_Comm & communicator,
456  const std::vector<size_type> &n_entries_per_row)
457  {
458  Epetra_Map map =
459  parallel_partitioning.make_trilinos_map(communicator, false);
460  reinit_sp(
461  map, map, n_entries_per_row, column_space_map, graph, nonlocal_graph);
462  }
463 
464 
465 
466  void
467  SparsityPattern::reinit(const IndexSet &row_parallel_partitioning,
468  const IndexSet &col_parallel_partitioning,
469  const MPI_Comm &communicator,
470  const size_type n_entries_per_row)
471  {
472  Epetra_Map row_map =
473  row_parallel_partitioning.make_trilinos_map(communicator, false);
474  Epetra_Map col_map =
475  col_parallel_partitioning.make_trilinos_map(communicator, false);
476  reinit_sp(row_map,
477  col_map,
478  n_entries_per_row,
480  graph,
482  }
483 
484 
485 
486  void
487  SparsityPattern::reinit(const IndexSet &row_parallel_partitioning,
488  const IndexSet &col_parallel_partitioning,
489  const MPI_Comm &communicator,
490  const std::vector<size_type> &n_entries_per_row)
491  {
492  Epetra_Map row_map =
493  row_parallel_partitioning.make_trilinos_map(communicator, false);
494  Epetra_Map col_map =
495  col_parallel_partitioning.make_trilinos_map(communicator, false);
496  reinit_sp(row_map,
497  col_map,
498  n_entries_per_row,
500  graph,
502  }
503 
504 
505 
506  void
507  SparsityPattern::reinit(const IndexSet &row_parallel_partitioning,
508  const IndexSet &col_parallel_partitioning,
509  const IndexSet &writable_rows,
510  const MPI_Comm &communicator,
511  const size_type n_entries_per_row)
512  {
513  Epetra_Map row_map =
514  row_parallel_partitioning.make_trilinos_map(communicator, false);
515  Epetra_Map col_map =
516  col_parallel_partitioning.make_trilinos_map(communicator, false);
517  reinit_sp(row_map,
518  col_map,
519  n_entries_per_row,
521  graph,
523 
524  IndexSet nonlocal_partitioner = writable_rows;
525  AssertDimension(nonlocal_partitioner.size(),
526  row_parallel_partitioning.size());
527 # ifdef DEBUG
528  {
529  IndexSet tmp = writable_rows & row_parallel_partitioning;
530  Assert(tmp == row_parallel_partitioning,
531  ExcMessage(
532  "The set of writable rows passed to this method does not "
533  "contain the locally owned rows, which is not allowed."));
534  }
535 # endif
536  nonlocal_partitioner.subtract_set(row_parallel_partitioning);
537  if (Utilities::MPI::n_mpi_processes(communicator) > 1)
538  {
539  Epetra_Map nonlocal_map =
540  nonlocal_partitioner.make_trilinos_map(communicator, true);
542  std_cxx14::make_unique<Epetra_CrsGraph>(Copy, nonlocal_map, 0);
543  }
544  else
545  Assert(nonlocal_partitioner.n_elements() == 0, ExcInternalError());
546  }
547 
548 
549 
550  template <typename SparsityPatternType>
551  void
553  const IndexSet & row_parallel_partitioning,
554  const IndexSet & col_parallel_partitioning,
555  const SparsityPatternType &nontrilinos_sparsity_pattern,
556  const MPI_Comm & communicator,
557  const bool exchange_data)
558  {
559  Epetra_Map row_map =
560  row_parallel_partitioning.make_trilinos_map(communicator, false);
561  Epetra_Map col_map =
562  col_parallel_partitioning.make_trilinos_map(communicator, false);
563  reinit_sp(row_map,
564  col_map,
565  nontrilinos_sparsity_pattern,
566  exchange_data,
568  graph,
570  }
571 
572 
573 
574  template <typename SparsityPatternType>
575  void
577  const IndexSet & parallel_partitioning,
578  const SparsityPatternType &nontrilinos_sparsity_pattern,
579  const MPI_Comm & communicator,
580  const bool exchange_data)
581  {
582  Epetra_Map map =
583  parallel_partitioning.make_trilinos_map(communicator, false);
584  reinit_sp(map,
585  map,
586  nontrilinos_sparsity_pattern,
587  exchange_data,
589  graph,
591  }
592 
593 
594 
597  {
598  Assert(false, ExcNotImplemented());
599  return *this;
600  }
601 
602 
603 
604  void
606  {
607  column_space_map = std_cxx14::make_unique<Epetra_Map>(*sp.column_space_map);
608  graph = std_cxx14::make_unique<Epetra_FECrsGraph>(*sp.graph);
609 
610  if (sp.nonlocal_graph.get() != nullptr)
612  std_cxx14::make_unique<Epetra_CrsGraph>(*sp.nonlocal_graph);
613  else
614  nonlocal_graph.reset();
615  }
616 
617 
618 
619  template <typename SparsityPatternType>
620  void
621  SparsityPattern::copy_from(const SparsityPatternType &sp)
622  {
623  const Epetra_Map rows(TrilinosWrappers::types::int_type(sp.n_rows()),
624  0,
626  const Epetra_Map columns(TrilinosWrappers::types::int_type(sp.n_cols()),
627  0,
629 
630  reinit_sp(
631  rows, columns, sp, false, column_space_map, graph, nonlocal_graph);
632  }
633 
634 
635 
636  void
638  {
639  // When we clear the matrix, reset
640  // the pointer and generate an
641  // empty sparsity pattern.
643  std_cxx14::make_unique<Epetra_Map>(TrilinosWrappers::types::int_type(0),
646  graph = std_cxx14::make_unique<Epetra_FECrsGraph>(View,
649  0);
650  graph->FillComplete();
651 
652  nonlocal_graph.reset();
653  }
654 
655 
656 
657  void
659  {
660  int ierr;
662  if (nonlocal_graph.get() != nullptr)
663  {
664  if (nonlocal_graph->IndicesAreGlobal() == false &&
665  nonlocal_graph->RowMap().NumMyElements() > 0 &&
667  {
668  // Insert dummy element at (row, column) that corresponds to row 0
669  // in local index counting.
673 
674  // in case we have a square sparsity pattern, add the entry on the
675  // diagonal
678  column = row;
679  // if not, take a column index that we have ourselves since we
680  // know for sure it is there (and it will not create spurious
681  // messages to many ranks like putting index 0 on many processors)
682  else if (column_space_map->NumMyElements() > 0)
684  ierr = nonlocal_graph->InsertGlobalIndices(row, 1, &column);
685  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
686  }
687  Assert(nonlocal_graph->RowMap().NumMyElements() == 0 ||
689  nonlocal_graph->IndicesAreGlobal() == true,
690  ExcInternalError());
691 
692  ierr =
693  nonlocal_graph->FillComplete(*column_space_map, graph->RangeMap());
694  AssertThrow(ierr >= 0, ExcTrilinosError(ierr));
695  ierr = nonlocal_graph->OptimizeStorage();
696  AssertThrow(ierr >= 0, ExcTrilinosError(ierr));
697  Epetra_Export exporter(nonlocal_graph->RowMap(), graph->RowMap());
698  ierr = graph->Export(*nonlocal_graph, exporter, Add);
699  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
700  ierr = graph->FillComplete(*column_space_map, graph->RangeMap());
701  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
702  }
703  else
704  {
705  // TODO A dynamic_cast fails here, this is suspicious.
706  const auto &range_map =
707  static_cast<const Epetra_Map &>(graph->RangeMap()); // NOLINT
708  ierr = graph->GlobalAssemble(*column_space_map, range_map, true);
709  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
710  }
711 
712  ierr = graph->OptimizeStorage();
713  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
714  }
715 
716 
717 
718  bool
720  {
721  return graph->RowMap().LID(
722  static_cast<TrilinosWrappers::types::int_type>(i)) != -1;
723  }
724 
725 
726 
727  bool
729  {
730  if (!row_is_stored_locally(i))
731  {
732  return false;
733  }
734  else
735  {
736  // Extract local indices in
737  // the matrix.
738  int trilinos_i =
739  graph->LRID(static_cast<TrilinosWrappers::types::int_type>(i)),
740  trilinos_j =
741  graph->LCID(static_cast<TrilinosWrappers::types::int_type>(j));
742 
743  // Check whether the matrix
744  // already is transformed to
745  // local indices.
746  if (graph->Filled() == false)
747  {
748  int nnz_present = graph->NumGlobalIndices(i);
749  int nnz_extracted;
751 
752  // Generate the view and make
753  // sure that we have not generated
754  // an error.
755  // TODO: trilinos_i is the local row index -> it is an int but
756  // ExtractGlobalRowView requires trilinos_i to be the global row
757  // index and thus it should be a long long int
758  int ierr = graph->ExtractGlobalRowView(trilinos_i,
759  nnz_extracted,
760  col_indices);
761  (void)ierr;
762  Assert(ierr == 0, ExcTrilinosError(ierr));
763  Assert(nnz_present == nnz_extracted,
764  ExcDimensionMismatch(nnz_present, nnz_extracted));
765 
766  // Search the index
767  const std::ptrdiff_t local_col_index =
768  std::find(col_indices, col_indices + nnz_present, trilinos_j) -
769  col_indices;
770 
771  if (local_col_index == nnz_present)
772  return false;
773  }
774  else
775  {
776  // Prepare pointers for extraction
777  // of a view of the row.
778  int nnz_present = graph->NumGlobalIndices(i);
779  int nnz_extracted;
780  int *col_indices;
781 
782  // Generate the view and make
783  // sure that we have not generated
784  // an error.
785  int ierr =
786  graph->ExtractMyRowView(trilinos_i, nnz_extracted, col_indices);
787  (void)ierr;
788  Assert(ierr == 0, ExcTrilinosError(ierr));
789 
790  Assert(nnz_present == nnz_extracted,
791  ExcDimensionMismatch(nnz_present, nnz_extracted));
792 
793  // Search the index
794  const std::ptrdiff_t local_col_index =
795  std::find(col_indices, col_indices + nnz_present, trilinos_j) -
796  col_indices;
797 
798  if (local_col_index == nnz_present)
799  return false;
800  }
801  }
802 
803  return true;
804  }
805 
806 
807 
810  {
811  size_type local_b = 0;
813  for (int i = 0; i < static_cast<int>(local_size()); ++i)
814  {
815  int *indices;
816  int num_entries;
817  graph->ExtractMyRowView(i, num_entries, indices);
818  for (unsigned int j = 0; j < static_cast<unsigned int>(num_entries);
819  ++j)
820  {
821  if (static_cast<size_type>(std::abs(i - indices[j])) > local_b)
822  local_b = std::abs(i - indices[j]);
823  }
824  }
825  graph->Comm().MaxAll(reinterpret_cast<TrilinosWrappers::types::int_type *>(
826  &local_b),
827  &global_b,
828  1);
829  return static_cast<size_type>(global_b);
830  }
831 
832 
833 
836  {
838  return n_rows;
839  }
840 
841 
842 
845  {
847  if (graph->Filled() == true)
849  else
851 
852  return n_cols;
853  }
854 
855 
856 
857  unsigned int
859  {
860  int n_rows = graph->NumMyRows();
861 
862  return n_rows;
863  }
864 
865 
866 
867  std::pair<SparsityPattern::size_type, SparsityPattern::size_type>
869  {
872  end = TrilinosWrappers::max_my_gid(graph->RowMap()) + 1;
873 
874  return std::make_pair(begin, end);
875  }
876 
877 
878 
881  {
883 
884  return static_cast<size_type>(nnz);
885  }
886 
887 
888 
889  unsigned int
891  {
892  int nnz = graph->MaxNumIndices();
893 
894  return static_cast<unsigned int>(nnz);
895  }
896 
897 
898 
901  {
902  Assert(row < n_rows(), ExcInternalError());
903 
904  // Get a representation of the where the present row is located on
905  // the current processor
907  graph->LRID(static_cast<TrilinosWrappers::types::int_type>(row));
908 
909  // On the processor who owns this row, we'll have a non-negative
910  // value for `local_row` and can ask for the length of the row.
911  if (local_row >= 0)
912  return graph->NumMyIndices(local_row);
913  else
914  return static_cast<size_type>(-1);
915  }
916 
917 
918 
919  const Epetra_Map &
921  {
922  // TODO A dynamic_cast fails here, this is suspicious.
923  const auto &domain_map =
924  static_cast<const Epetra_Map &>(graph->DomainMap()); // NOLINT
925  return domain_map;
926  }
927 
928 
929 
930  const Epetra_Map &
932  {
933  // TODO A dynamic_cast fails here, this is suspicious.
934  const auto &range_map =
935  static_cast<const Epetra_Map &>(graph->RangeMap()); // NOLINT
936  return range_map;
937  }
938 
939 
940 
941  MPI_Comm
943  {
944 # ifdef DEAL_II_WITH_MPI
945 
946  const Epetra_MpiComm *mpi_comm =
947  dynamic_cast<const Epetra_MpiComm *>(&graph->RangeMap().Comm());
948  Assert(mpi_comm != nullptr, ExcInternalError());
949  return mpi_comm->Comm();
950 # else
951 
952  return MPI_COMM_SELF;
953 
954 # endif
955  }
956 
957 
958 
959  void
961  {
962  Assert(false, ExcNotImplemented());
963  }
964 
965 
966 
967  // As of now, no particularly neat
968  // output is generated in case of
969  // multiple processors.
970  void
971  SparsityPattern::print(std::ostream &out,
972  const bool write_extended_trilinos_info) const
973  {
974  if (write_extended_trilinos_info)
975  out << *graph;
976  else
977  {
978  int *indices;
979  int num_entries;
980 
981  for (int i = 0; i < graph->NumMyRows(); ++i)
982  {
983  graph->ExtractMyRowView(i, num_entries, indices);
984  for (int j = 0; j < num_entries; ++j)
985  out << "(" << TrilinosWrappers::global_index(graph->RowMap(), i)
986  << ","
987  << TrilinosWrappers::global_index(graph->ColMap(), indices[j])
988  << ") " << std::endl;
989  }
990  }
991 
992  AssertThrow(out, ExcIO());
993  }
994 
995 
996 
997  void
998  SparsityPattern::print_gnuplot(std::ostream &out) const
999  {
1000  Assert(graph->Filled() == true, ExcInternalError());
1001  for (::types::global_dof_index row = 0; row < local_size(); ++row)
1002  {
1003  int *indices;
1004  int num_entries;
1005  graph->ExtractMyRowView(row, num_entries, indices);
1006 
1007  Assert(num_entries >= 0, ExcInternalError());
1008  // avoid sign comparison warning
1009  const ::types::global_dof_index num_entries_ = num_entries;
1010  for (::types::global_dof_index j = 0; j < num_entries_; ++j)
1011  // while matrix entries are usually
1012  // written (i,j), with i vertical and
1013  // j horizontal, gnuplot output is
1014  // x-y, that is we have to exchange
1015  // the order of output
1016  out << static_cast<int>(
1017  TrilinosWrappers::global_index(graph->ColMap(), indices[j]))
1018  << " "
1019  << -static_cast<int>(
1020  TrilinosWrappers::global_index(graph->RowMap(), row))
1021  << std::endl;
1022  }
1023 
1024  AssertThrow(out, ExcIO());
1025  }
1026 
1027  // TODO: Implement!
1028  std::size_t
1030  {
1031  Assert(false, ExcNotImplemented());
1032  return 0;
1033  }
1034 
1035 
1036 # ifndef DOXYGEN
1037  // explicit instantiations
1038  //
1039  template void
1040  SparsityPattern::copy_from(const ::SparsityPattern &);
1041  template void
1042  SparsityPattern::copy_from(const ::DynamicSparsityPattern &);
1043 
1044  template void
1046  const ::SparsityPattern &,
1047  const MPI_Comm &,
1048  bool);
1049  template void
1051  const ::DynamicSparsityPattern &,
1052  const MPI_Comm &,
1053  bool);
1054 
1055 
1056  template void
1058  const IndexSet &,
1059  const ::SparsityPattern &,
1060  const MPI_Comm &,
1061  bool);
1062  template void
1064  const IndexSet &,
1065  const ::DynamicSparsityPattern &,
1066  const MPI_Comm &,
1067  bool);
1068 # endif
1069 
1070 } // namespace TrilinosWrappers
1071 
1073 
1074 #endif // DEAL_II_WITH_TRILINOS
IndexSet
Definition: index_set.h:74
IndexSet::subtract_set
void subtract_set(const IndexSet &other)
Definition: index_set.cc:258
dynamic_sparsity_pattern.h
TrilinosWrappers::SparsityPattern::n_nonzero_elements
size_type n_nonzero_elements() const
Definition: trilinos_sparsity_pattern.cc:880
StandardExceptions::ExcIO
static ::ExceptionBase & ExcIO()
TrilinosWrappers::SparsityPattern::range_partitioner
const Epetra_Map & range_partitioner() const
Definition: trilinos_sparsity_pattern.cc:931
trilinos_sparsity_pattern.h
TrilinosWrappers::SparsityPattern::row_length
size_type row_length(const size_type row) const
Definition: trilinos_sparsity_pattern.cc:900
TrilinosWrappers::n_global_entries
TrilinosWrappers::types::int_type n_global_entries(const Epetra_CrsGraph &graph)
Definition: trilinos_index_access.h:140
StandardExceptions::ExcNotImplemented
static ::ExceptionBase & ExcNotImplemented()
TrilinosWrappers::SparsityPattern::row_is_stored_locally
bool row_is_stored_locally(const size_type i) const
Definition: trilinos_sparsity_pattern.cc:719
IndexSet::size
size_type size() const
Definition: index_set.h:1635
TrilinosWrappers::SparsityPatternIterators::Accessor::sparsity_pattern
SparsityPattern * sparsity_pattern
Definition: trilinos_sparsity_pattern.h:133
utilities.h
TrilinosWrappers::SparsityPattern::domain_partitioner
const Epetra_Map & domain_partitioner() const
Definition: trilinos_sparsity_pattern.cc:920
TrilinosWrappers::SparsityPattern::local_size
unsigned int local_size() const
Definition: trilinos_sparsity_pattern.cc:858
TrilinosWrappers::SparsityPattern::column_space_map
std::unique_ptr< Epetra_Map > column_space_map
Definition: trilinos_sparsity_pattern.h:995
TrilinosWrappers
Definition: types.h:161
TrilinosWrappers::max_my_gid
TrilinosWrappers::types::int_type max_my_gid(const Epetra_BlockMap &map)
Definition: trilinos_index_access.h:68
TrilinosWrappers::SparsityPattern::exists
bool exists(const size_type i, const size_type j) const
Definition: trilinos_sparsity_pattern.cc:728
trilinos_index_access.h
MPI_Comm
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::SparsityPattern::size_type
::types::global_dof_index size_type
Definition: trilinos_sparsity_pattern.h:285
TrilinosWrappers::n_global_cols
TrilinosWrappers::types::int_type n_global_cols(const Epetra_CrsGraph &graph)
Definition: trilinos_index_access.h:126
TrilinosWrappers::SparsityPatternIterators::Accessor::colnum_cache
std::shared_ptr< const std::vector< size_type > > colnum_cache
Definition: trilinos_sparsity_pattern.h:157
TrilinosWrappers::SparsityPattern::ExcTrilinosError
static ::ExceptionBase & ExcTrilinosError(int arg1)
TrilinosWrappers::SparsityPattern
Definition: trilinos_sparsity_pattern.h:279
TrilinosWrappers::SparsityPattern::SparsityPattern
SparsityPattern()
Definition: trilinos_sparsity_pattern.cc:83
Subscriptor
Definition: subscriptor.h:62
IndexSet::n_elements
size_type n_elements() const
Definition: index_set.h:1833
TrilinosWrappers::SparsityPattern::write_ascii
void write_ascii()
Definition: trilinos_sparsity_pattern.cc:960
TrilinosWrappers::n_global_rows
TrilinosWrappers::types::int_type n_global_rows(const Epetra_CrsGraph &graph)
Definition: trilinos_index_access.h:112
StandardExceptions::ExcMessage
static ::ExceptionBase & ExcMessage(std::string arg1)
types::global_dof_index
unsigned int global_dof_index
Definition: types.h:76
TrilinosWrappers::SparsityPattern::memory_consumption
std::size_t memory_consumption() const
Definition: trilinos_sparsity_pattern.cc:1029
mpi.h
DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:358
TrilinosWrappers::SparsityPattern::begin
const_iterator begin() const
TrilinosWrappers::SparsityPattern::n_rows
size_type n_rows() const
Definition: trilinos_sparsity_pattern.cc:835
TrilinosWrappers::SparsityPattern::print
void print(std::ostream &out, const bool write_extended_trilinos_info=false) const
Definition: trilinos_sparsity_pattern.cc:971
TrilinosWrappers::SparsityPattern::print_gnuplot
void print_gnuplot(std::ostream &out) const
Definition: trilinos_sparsity_pattern.cc:998
SparsityPatternIterators
Definition: sparsity_pattern.h:111
TrilinosWrappers::SparsityPattern::bandwidth
size_type bandwidth() const
Definition: trilinos_sparsity_pattern.cc:809
TrilinosWrappers::SparsityPattern::clear
void clear()
Definition: trilinos_sparsity_pattern.cc:637
AssertDimension
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1579
sparsity_pattern.h
TrilinosWrappers::SparsityPatternIterators::Accessor::visit_present_row
void visit_present_row()
Definition: trilinos_sparsity_pattern.cc:36
TrilinosWrappers::SparsityPattern::get_mpi_communicator
MPI_Comm get_mpi_communicator() const
Definition: trilinos_sparsity_pattern.cc:942
Utilities::Trilinos::comm_self
const Epetra_Comm & comm_self()
Definition: utilities.cc:1114
types
Definition: types.h:31
complete_index_set
IndexSet complete_index_set(const IndexSet::size_type N)
Definition: index_set.h:1014
TrilinosWrappers::types::int_type
int int_type
Definition: types.h:174
TrilinosWrappers::SparsityPattern::compress
void compress()
Definition: trilinos_sparsity_pattern.cc:658
TrilinosWrappers::SparsityPattern::graph
std::unique_ptr< Epetra_FECrsGraph > graph
Definition: trilinos_sparsity_pattern.h:1002
unsigned int
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
LinearAlgebra::CUDAWrappers::kernel::size_type
types::global_dof_index size_type
Definition: cuda_kernels.h:45
LACExceptions::ExcTrilinosError
static ::ExceptionBase & ExcTrilinosError(int arg1)
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::SparsityPattern::max_entries_per_row
unsigned int max_entries_per_row() const
Definition: trilinos_sparsity_pattern.cc:890
TrilinosWrappers::SparsityPattern::is_compressed
bool is_compressed() const
Utilities::MPI::n_mpi_processes
unsigned int n_mpi_processes(const MPI_Comm &mpi_communicator)
Definition: mpi.cc:117
TrilinosWrappers::SparsityPattern::copy_from
void copy_from(const SparsityPattern &input_sparsity_pattern)
Definition: trilinos_sparsity_pattern.cc:605
TrilinosWrappers::SparsityPattern::end
const_iterator end() const
StandardExceptions::ExcDimensionMismatch
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
TrilinosWrappers::SparsityPatternIterators::Accessor::a_row
size_type a_row
Definition: trilinos_sparsity_pattern.h:138
TrilinosWrappers::SparsityPattern::operator=
SparsityPattern & operator=(const SparsityPattern &input_sparsity_pattern)
Definition: trilinos_sparsity_pattern.cc:596
TrilinosWrappers::SparsityPattern::reinit
void reinit(const size_type m, const size_type n, const size_type n_entries_per_row=0)
Definition: trilinos_sparsity_pattern.cc:212
DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:359
TrilinosWrappers::SparsityPattern::local_range
std::pair< size_type, size_type > local_range() const
Definition: trilinos_sparsity_pattern.cc:868
TrilinosWrappers::global_index
TrilinosWrappers::types::int_type global_index(const Epetra_BlockMap &map, const ::types::global_dof_index i)
Definition: trilinos_index_access.h:82
AssertThrow
#define AssertThrow(cond, exc)
Definition: exceptions.h:1531
TrilinosWrappers::SparsityPattern::n_cols
size_type n_cols() const
Definition: trilinos_sparsity_pattern.cc:844
FEValuesViews::View
typename ::internal::FEValuesViews::ViewType< dim, spacedim, Extractor >::type View
Definition: fe_values.h:1980
Utilities
Definition: cuda.h:31
TrilinosWrappers::n_global_elements
TrilinosWrappers::types::int_type n_global_elements(const Epetra_BlockMap &map)
Definition: trilinos_index_access.h:40
int