<ol>
+ <li> <p> New: There are now assignment operators from <code
+ class="class">BlockVector</code> to <code
+ class="class">Vector</code> and back.
+ <br>
+ (WB 2006/03/30)
+ </p>
+
<li> <p> Improved: <code class="class">BlockSparsityPattern</code> can be
- initialized directly using the vector generated by
- <code class="class">DoFTools</code>::<code
- class="member">compute_row_length_vector</code>.
- <br>
- (GK 2006/03/30)
- </p>
+ initialized directly using the vector generated by
+ <code class="class">DoFTools</code>::<code
+ class="member">compute_row_length_vector</code>.
+ <br>
+ (GK 2006/03/30)
+ </p>
<li> <p>
Improved: All matrix (and some vector) classes now check whether
BlockVector &
operator= (const BlockVector<Number2> &V);
- /**
+ /**
+ * Copy a regular vector into a
+ * block vector.
+ */
+ BlockVector &
+ operator= (const Vector<Number> &V);
+
+ /**
* Reinitialize the BlockVector to
* contain <tt>num_blocks</tt> blocks of
* size <tt>block_size</tt> each.
+template <typename Number>
+inline
+BlockVector<Number> &
+BlockVector<Number>::operator = (const Vector<Number> &v)
+{
+ BaseClass::operator = (v);
+ return *this;
+}
+
+
+
template <typename Number>
template <typename Number2>
inline
BlockVectorBase &
operator= (const BlockVectorBase<VectorType2> &V);
+ /**
+ * Copy operator from non-block
+ * vectors to block vectors.
+ */
+ BlockVectorBase &
+ operator = (const VectorType &v);
+
/**
* Check for equality of two block vector
* types. This operation is only allowed
}
+
template <class VectorType>
BlockVectorBase<VectorType>&
BlockVectorBase<VectorType>::operator = (const value_type s)
ExcDimensionMismatch(n_blocks(), v.n_blocks()));
for (unsigned int i=0;i<n_blocks();++i)
- {
- components[i] = v.components[i];
- }
+ components[i] = v.components[i];
+
return *this;
}
+template <class VectorType>
+BlockVectorBase<VectorType>&
+BlockVectorBase<VectorType>::operator = (const VectorType &v)
+{
+ Assert (size() == v.size(),
+ ExcDimensionMismatch(size(), v.size()));
+
+ unsigned int index_v = 0;
+ for (unsigned int b=0;b<n_blocks();++b)
+ for (unsigned int i=0; i<block(b).size(); ++i, ++index_v)
+ block(b)(i) = v(index_v);
+
+ return *this;
+}
+
+
+
template <class VectorType>
template <class VectorType2>
inline
template<typename number> class LAPACKFullMatrix;
+template <typename> class BlockVector;
+
+
+
/*! @addtogroup Vectors
*@{
*/
template <typename Number2>
Vector<Number> & operator= (const Vector<Number2> &v);
+ /**
+ * Copy operator for assigning a
+ * block vector to a regular
+ * vector.
+ */
+ Vector<Number> & operator= (const BlockVector<Number> &v);
+
#ifdef DEAL_II_USE_PETSC
/**
* Another copy operator: copy the values
#include <lac/vector.h>
+#include <lac/block_vector.h>
#ifdef DEAL_II_USE_PETSC
# include <lac/petsc_vector.h>
template <typename Number>
-Vector<Number>&
+Vector<Number> &
Vector<Number>::operator = (const Vector<Number>& v)
{
if (v.vec_size != vec_size)
template <typename Number>
template <typename Number2>
-Vector<Number>&
+Vector<Number> &
Vector<Number>::operator = (const Vector<Number2>& v)
{
if (v.size() != vec_size)
+template <typename Number>
+Vector<Number> &
+Vector<Number>::operator = (const BlockVector<Number>& v)
+{
+ if (v.size() != vec_size)
+ reinit (v.size(), true);
+
+ unsigned int this_index = 0;
+ for (unsigned int b=0; b<v.n_blocks(); ++b)
+ for (unsigned int i=0; i<v.block(b).size(); ++i, ++this_index)
+ val[this_index] = v.block(b)(i);
+
+ return *this;
+}
+
+
+
#ifdef DEAL_II_USE_PETSC
template <typename Number>
full_matrix householder tridiagonal_matrix \
sparsity_pattern sparse_matrices matrices \
block_sparsity_pattern_01 \
- block_vector block_vector_iterator block_matrices \
+ block_vector block_vector_* block_matrices \
matrix_lib matrix_out \
solver eigen \
- sparse_ilu sparse_ilu_t lapack lapack_fill block_minres
+ sparse_ilu sparse_ilu_t lapack lapack_fill block_minres \
+ amg_*
# from above list of regular expressions, generate the real set of
--- /dev/null
+//--------------------------------------------------------------------
+// $Id$
+// Version: $Name$
+//
+// Copyright (C) 2006 by the deal.II authors
+//
+// This file is subject to QPL and may not be distributed
+// without copyright and license information. Please refer
+// to the file deal.II/doc/license.html for the text and
+// further information on this license.
+//
+//--------------------------------------------------------------------
+
+// check assignment between block vectors and regular vectors
+
+#include "../tests.h"
+#include <base/logstream.h>
+#include <lac/block_vector.h>
+#include <fstream>
+#include <iostream>
+#include <vector>
+#include <algorithm>
+#include <numeric>
+#include <utility>
+#include <cmath>
+
+template <typename Vector1, typename Vector2>
+bool operator == (const Vector1 &v1,
+ const Vector2 &v2)
+{
+ if (v1.size() != v2.size())
+ return false;
+ for (unsigned int i=0; i<v1.size(); ++i)
+ if (v1(i) != v2(i))
+ return false;
+ return true;
+}
+
+
+void test ()
+{
+ std::vector<unsigned int> ivector(4);
+ ivector[0] = 2;
+ ivector[1] = 4;
+ ivector[2] = 3;
+ ivector[3] = 5;
+
+ BlockVector<double> v1(ivector);
+ Vector<double> v2(v1.size());
+
+ for (unsigned int i=0; i<v1.size(); ++i)
+ v1(i) = 1+i*i;
+
+ v2 = v1;
+ Assert (v1==v2, ExcInternalError());
+
+ BlockVector<double> v3 (ivector);
+ v3 = v2;
+ Assert (v3==v2, ExcInternalError());
+ Assert (v3==v1, ExcInternalError());
+
+ deallog << "OK" << std::endl;
+}
+
+
+
+
+
+int main ()
+{
+ std::ofstream logfile("block_vector_vector_assign/output");
+ logfile.setf(std::ios::fixed);
+ logfile.precision(3);
+ deallog.attach(logfile);
+ deallog.depth_console(0);
+ deallog.threshold_double(1.e-10);
+
+ try
+ {
+ test ();
+ }
+ catch (std::exception &e)
+ {
+ std::cerr << std::endl << std::endl
+ << "----------------------------------------------------"
+ << std::endl;
+ std::cerr << "Exception on processing: " << e.what() << std::endl
+ << "Aborting!" << std::endl
+ << "----------------------------------------------------"
+ << std::endl;
+ // abort
+ return 2;
+ }
+ catch (...)
+ {
+ std::cerr << std::endl << std::endl
+ << "----------------------------------------------------"
+ << std::endl;
+ std::cerr << "Unknown exception!" << std::endl
+ << "Aborting!" << std::endl
+ << "----------------------------------------------------"
+ << std::endl;
+ // abort
+ return 3;
+ };
+
+
+ return 0;
+}
+
--- /dev/null
+
+DEAL::OK