// pattern as the input matrix
if (&this->get_sparsity_pattern() == &matrix.get_sparsity_pattern())
{
- const somenumber *input_ptr = matrix.val;
- number *this_ptr = this->val;
+ const somenumber *input_ptr = matrix.val.get();
+ number *this_ptr = this->val.get();
const number *const end_ptr = this_ptr + this->n_nonzero_elements();
if (types_are_equal<somenumber, number>::value == true)
std::memcpy (this_ptr, input_ptr, this->n_nonzero_elements()*sizeof(number));
const std::size_t *const ia = sparsity.rowstart;
const size_type *const ja = sparsity.colnums;
- number *luval = this->SparseMatrix<number>::val;
+ number *luval = this->SparseMatrix<number>::val.get();
const size_type N = this->m();
size_type jrow = 0;
const size_type *const first_after_diagonal = this->prebuilt_lower_bound[row];
somenumber dst_row = dst(row);
- const number *luval = this->SparseMatrix<number>::val +
+ const number *luval = this->SparseMatrix<number>::val.get() +
(rowstart - column_numbers);
for (const size_type *col=rowstart; col!=first_after_diagonal; ++col, ++luval)
dst_row -= *luval * dst(*col);
const size_type *const first_after_diagonal = this->prebuilt_lower_bound[row];
somenumber dst_row = dst(row);
- const number *luval = this->SparseMatrix<number>::val +
+ const number *luval = this->SparseMatrix<number>::val.get() +
(first_after_diagonal - column_numbers);
for (const size_type *col=first_after_diagonal; col!=rowend; ++col, ++luval)
dst_row -= *luval * dst(*col);
const size_type *const first_after_diagonal = this->prebuilt_lower_bound[row];
const somenumber dst_row = dst (row);
- const number *luval = this->SparseMatrix<number>::val +
+ const number *luval = this->SparseMatrix<number>::val.get() +
(first_after_diagonal - column_numbers);
for (const size_type *col=first_after_diagonal; col!=rowend; ++col, ++luval)
tmp(*col) += *luval * dst_row;
const size_type *const first_after_diagonal = this->prebuilt_lower_bound[row];
const somenumber dst_row = dst (row);
- const number *luval = this->SparseMatrix<number>::val +
+ const number *luval = this->SparseMatrix<number>::val.get() +
(rowstart - column_numbers);
for (const size_type *col=rowstart; col!=first_after_diagonal; ++col, ++luval)
tmp(*col) += *luval * dst_row;
#include <deal.II/lac/exceptions.h>
#include <deal.II/lac/vector.h>
+#include <memory>
+
+
DEAL_II_NAMESPACE_OPEN
template <typename number> class Vector;
SmartPointer<const SparsityPattern,SparseMatrix<number> > cols;
/**
- * Array of values for all the nonzero entries. The position within the
- * matrix, i.e. the row and column number for a given entry can only be
- * deduced using the sparsity pattern. The same holds for the more common
- * operation of finding an entry by its coordinates.
+ * Array of values for all the nonzero entries. The position of an
+ * entry within the matrix, i.e., the row and column number for a
+ * given value in this array can only be deduced using the sparsity
+ * pattern. The same holds for the more common operation of finding
+ * an entry by its coordinates.
*/
- number *val;
+ std::unique_ptr<number[]> val;
/**
* Allocated size of #val. This can be larger than the actually used part if
:
Subscriptor(std::move(m)),
cols(m.cols),
- val(m.val),
+ val(std::move(m.val)),
max_len(m.max_len)
{
m.cols = nullptr;
SparseMatrix<number>::operator = (SparseMatrix<number> &&m)
{
cols = m.cols;
- val = m.val;
+ val = std::move(m.val);
max_len = m.max_len;
m.cols = nullptr;
SparseMatrix<number>::~SparseMatrix ()
{
cols = nullptr;
-
- if (val != nullptr)
- delete[] val;
}
std::bind(&internal::SparseMatrix::template
zero_subrange<number>,
std::placeholders::_1, std::placeholders::_2,
- val),
+ val.get()),
grain_size);
else if (matrix_size > 0)
std::memset (&val[0], 0, matrix_size*sizeof(number));
if (cols->empty())
{
- if (val != nullptr)
- delete[] val;
- val = nullptr;
+ val.reset ();
max_len = 0;
return;
}
const std::size_t N = cols->n_nonzero_elements();
if (N > max_len || max_len == 0)
{
- if (val != nullptr)
- delete[] val;
- val = new number[N];
+ val.reset (new number[N]);
max_len = N;
}
SparseMatrix<number>::clear ()
{
cols = nullptr;
- if (val) delete[] val;
- val = nullptr;
+ val.reset ();
max_len = 0;
}
std::bind (&internal::SparseMatrix::vmult_on_subrange
<number,InVector,OutVector>,
std::placeholders::_1, std::placeholders::_2,
- val,
+ val.get(),
cols->rowstart,
cols->colnums,
std::cref(src),
std::bind (&internal::SparseMatrix::vmult_on_subrange
<number,InVector,OutVector>,
std::placeholders::_1, std::placeholders::_2,
- val,
+ val.get(),
cols->rowstart,
cols->colnums,
std::cref(src),
(std::bind (&internal::SparseMatrix::matrix_norm_sqr_on_subrange
<number,Vector<somenumber> >,
std::placeholders::_1, std::placeholders::_2,
- val, cols->rowstart, cols->colnums,
+ val.get(), cols->rowstart, cols->colnums,
std::cref(v)),
0, m(),
internal::SparseMatrix::minimum_parallel_grain_size);
(std::bind (&internal::SparseMatrix::matrix_scalar_product_on_subrange
<number,Vector<somenumber> >,
std::placeholders::_1, std::placeholders::_2,
- val, cols->rowstart, cols->colnums,
+ val.get(), cols->rowstart, cols->colnums,
std::cref(u),
std::cref(v)),
0, m(),
(std::bind (&internal::SparseMatrix::residual_sqr_on_subrange
<number,Vector<somenumber>,Vector<somenumber> >,
std::placeholders::_1, std::placeholders::_2,
- val, cols->rowstart, cols->colnums,
+ val.get(), cols->rowstart, cols->colnums,
std::cref(u),
std::cref(b),
std::ref(dst)),
AssertThrow (c == '[', ExcIO());
// reallocate space
- delete[] val;
- val = new number[max_len];
+ val.reset (new number[max_len]);
// then read data
in.read (reinterpret_cast<char *>(&val[0]),