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