17 #include <deal.II/lac/chunk_sparsity_pattern.h> 18 #include <deal.II/lac/dynamic_sparsity_pattern.h> 19 #include <deal.II/lac/full_matrix.h> 22 DEAL_II_NAMESPACE_OPEN
35 chunk_size (s.chunk_size),
36 sparsity_pattern(s.sparsity_pattern)
39 ExcMessage(
"This constructor can only be called if the provided argument " 40 "is the sparsity pattern for an empty matrix. This constructor can " 41 "not be used to copy-construct a non-empty sparsity pattern."));
63 const std::vector<size_type> &row_lengths,
84 const std::vector<size_type > &row_lengths,
103 ExcMessage(
"This operator can only be called if the provided argument " 104 "is the sparsity pattern for an empty matrix. This operator can " 105 "not be used to copy a non-empty sparsity pattern."));
108 ExcMessage(
"This operator can only be called if the current object is " 128 const std::vector<size_type> row_lengths (m, max_per_row);
138 const VectorSlice<
const std::vector<size_type> > &row_lengths,
162 std::vector<unsigned int> chunk_row_lengths (m_chunks, 0);
164 chunk_row_lengths[i/
chunk_size] += row_lengths[i];
169 if (m != n && m_chunks == n_chunks)
170 for (
unsigned int i=0; i<m_chunks; ++i)
171 ++chunk_row_lengths[i];
188 template <
typename SparsityPatternType>
212 for (
size_type row = 0; row<dsp.n_rows(); ++row)
219 for (
typename SparsityPatternType::iterator col_num = dsp.begin(row);
220 col_num != dsp.end(row); ++col_num)
221 temporary_sp.
add (reduced_row, col_num->column()/
chunk_size);
230 template <
typename number>
240 std::vector<size_type> entries_per_row (matrix.m(), 0);
241 for (
size_type row=0; row<matrix.m(); ++row)
243 for (
size_type col=0; col<matrix.n(); ++col)
244 if (matrix(row,col) != 0)
245 ++entries_per_row[row];
247 if ((matrix.m() == matrix.n())
249 (matrix(row,row) == 0))
250 ++entries_per_row[row];
253 reinit (matrix.m(), matrix.n(),
258 for (
size_type row=0; row<matrix.m(); ++row)
259 for (
size_type col=0; col<matrix.n(); ++col)
260 if (matrix(row,col) != 0)
273 const std::vector<size_type> &row_lengths,
287 template <
typename SparsityPatternType>
288 void copy_sparsity (
const SparsityPatternType &src,
304 template <
typename Sparsity>
307 (
const unsigned int m,
308 const unsigned int n,
309 const Sparsity &sparsity_pattern_for_chunks,
310 const unsigned int chunk_size_in,
313 Assert (m > (sparsity_pattern_for_chunks.n_rows()-1) * chunk_size_in &&
314 m <= sparsity_pattern_for_chunks.n_rows() * chunk_size_in,
315 ExcMessage(
"Number of rows m is not compatible with chunk size " 316 "and number of rows in sparsity pattern for the chunks."));
317 Assert (n > (sparsity_pattern_for_chunks.n_cols()-1) * chunk_size_in &&
318 n <= sparsity_pattern_for_chunks.n_cols() * chunk_size_in,
319 ExcMessage(
"Number of columns m is not compatible with chunk size " 320 "and number of columns in sparsity pattern for the chunks."));
322 internal::copy_sparsity(sparsity_pattern_for_chunks, sparsity_pattern);
323 chunk_size = chunk_size_in;
398 for ( ; p !=
end; ++p)
499 out <<
']' << std::endl;
618 return (
sizeof(*
this) +
629 void ChunkSparsityPattern::create_from<SparsityPattern>
636 void ChunkSparsityPattern::create_from<DynamicSparsityPattern>
649 DEAL_II_NAMESPACE_CLOSE
void block_write(std::ostream &out) const
void block_write(std::ostream &out) const
void copy_from(const size_type n_rows, const size_type n_cols, const ForwardIterator begin, const ForwardIterator end, const size_type chunk_size)
std::size_t memory_consumption() const
static ::ExceptionBase & ExcEmptyObject()
std::size_t n_nonzero_elements() const
static ::ExceptionBase & ExcInvalidNumber(int arg1)
static ::ExceptionBase & ExcInvalidIndex(int arg1, int arg2)
static ::ExceptionBase & ExcIO()
static const size_type invalid_entry
void add(const size_type i, const size_type j)
std::size_t memory_consumption() const
void block_read(std::istream &in)
bool exists(const size_type i, const size_type j) const
bool exists(const size_type i, const size_type j) const
#define AssertThrow(cond, exc)
void add(const size_type i, const size_type j)
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
SparsityPattern sparsity_pattern
size_type bandwidth() const
static ::ExceptionBase & ExcMessage(std::string arg1)
void add(const size_type i, const size_type j)
#define Assert(cond, exc)
size_type row_length(const size_type row) const
size_type max_entries_per_row() const
size_type n_nonzero_elements() const
void copy_from(const size_type n_rows, const size_type n_cols, const ForwardIterator begin, const ForwardIterator end)
ChunkSparsityPattern & operator=(const ChunkSparsityPattern &)
void reinit(const size_type m, const size_type n, const unsigned int max_per_row)
types::global_dof_index size_type
static ::ExceptionBase & ExcNotQuadratic()
bool stores_only_added_elements() const
void block_read(std::istream &in)
size_type bandwidth() const
bool stores_only_added_elements() const
void print(std::ostream &out) const
void create_from(const unsigned int m, const unsigned int n, const Sparsity &sparsity_pattern_for_chunks, const unsigned int chunk_size, const bool optimize_diagonal=true)
void print_gnuplot(std::ostream &out) const
void reinit(const size_type m, const size_type n, const size_type max_per_row, const size_type chunk_size)
unsigned int row_length(const size_type row) const
size_type max_entries_per_row() const