Reference documentation for deal.II version 8.5.1
sparsity_tools.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2008 - 2016 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 #ifndef dealii__sparsity_tools_h
17 #define dealii__sparsity_tools_h
18 
19 
20 #include <deal.II/base/config.h>
21 #include <deal.II/base/exceptions.h>
22 #include <deal.II/lac/block_sparsity_pattern.h>
23 #include <deal.II/lac/dynamic_sparsity_pattern.h>
24 #include <deal.II/lac/sparsity_pattern.h>
25 
26 #include <vector>
27 
28 #ifdef DEAL_II_WITH_MPI
29 #include <mpi.h>
30 #include <deal.II/base/index_set.h>
31 #endif
32 
33 DEAL_II_NAMESPACE_OPEN
34 
35 
46 namespace SparsityTools
47 {
80  void partition (const SparsityPattern &sparsity_pattern,
81  const unsigned int n_partitions,
82  std::vector<unsigned int> &partition_indices);
83 
133  void
135  std::vector<DynamicSparsityPattern::size_type> &new_indices,
136  const std::vector<DynamicSparsityPattern::size_type> &starting_indices = std::vector<DynamicSparsityPattern::size_type>());
137 
143  void
144  reorder_Cuthill_McKee (const SparsityPattern &sparsity,
145  std::vector<SparsityPattern::size_type> &new_indices,
146  const std::vector<SparsityPattern::size_type> &starting_indices = std::vector<SparsityPattern::size_type>()) DEAL_II_DEPRECATED;
147 
169  void
171  std::vector<DynamicSparsityPattern::size_type> &new_indices);
172 
173 #ifdef DEAL_II_WITH_MPI
174 
200  const std::vector<DynamicSparsityPattern::size_type> &rows_per_cpu,
201  const MPI_Comm &mpi_comm,
202  const IndexSet &myrange);
203 
218  const std::vector<IndexSet> &owned_set_per_cpu,
219  const MPI_Comm &mpi_comm,
220  const IndexSet &myrange);
221 
222 #endif
223 
224 
229  "The function you called requires METIS, but you did not "
230  "configure deal.II with METIS.");
231 
236  int,
237  << "The number of partitions you gave is " << arg1
238  << ", but must be greater than zero.");
239 
244  int,
245  << " An error with error number " << arg1
246  << " occurred while calling a METIS function");
247 
252  int, int,
253  << "The array has size " << arg1 << " but should have size "
254  << arg2);
255 }
256 
261 DEAL_II_NAMESPACE_CLOSE
262 
263 #endif
#define DeclException2(Exception2, type1, type2, outsequence)
Definition: exceptions.h:576
static ::ExceptionBase & ExcMETISNotInstalled()
void reorder_Cuthill_McKee(const DynamicSparsityPattern &sparsity, std::vector< DynamicSparsityPattern::size_type > &new_indices, const std::vector< DynamicSparsityPattern::size_type > &starting_indices=std::vector< DynamicSparsityPattern::size_type >())
void partition(const SparsityPattern &sparsity_pattern, const unsigned int n_partitions, std::vector< unsigned int > &partition_indices)
static ::ExceptionBase & ExcInvalidNumberOfPartitions(int arg1)
#define DeclException1(Exception1, type1, outsequence)
Definition: exceptions.h:564
void distribute_sparsity_pattern(DynamicSparsityPattern &dsp, const std::vector< DynamicSparsityPattern::size_type > &rows_per_cpu, const MPI_Comm &mpi_comm, const IndexSet &myrange)
static ::ExceptionBase & ExcInvalidArraySize(int arg1, int arg2)
#define DeclExceptionMsg(Exception, defaulttext)
Definition: exceptions.h:553
void reorder_hierarchical(const DynamicSparsityPattern &sparsity, std::vector< DynamicSparsityPattern::size_type > &new_indices)
static ::ExceptionBase & ExcMETISError(int arg1)