]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Fix a few things that did not work when instantiating
authorWolfgang Bangerth <bangerth@math.tamu.edu>
Tue, 22 Jun 2004 16:10:35 +0000 (16:10 +0000)
committerWolfgang Bangerth <bangerth@math.tamu.edu>
Tue, 22 Jun 2004 16:10:35 +0000 (16:10 +0000)
vmult(Vector<double>,Vector<float>) and similar functions.

git-svn-id: https://svn.dealii.org/trunk@9452 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/lac/include/lac/sparse_matrix.templates.h

index d7a3648b93b839849965a5074d5a3e74de8fa0a8..b732c9dc7aa373ea806f327cdd264697b88c17af 100644 (file)
@@ -14,6 +14,7 @@
 #define __deal2__sparse_matrix_templates_h
 
 
+#include <base/template_constraints.h>
 #include <lac/sparse_matrix.h>
 #include <lac/vector.h>
 #include <lac/full_matrix.h>
@@ -79,7 +80,8 @@ SparseMatrix<number>::operator = (const SparseMatrix<number> &m)
 
 
 template <typename number>
-SparseMatrix<number>::SparseMatrix (const SparsityPattern &c) :
+SparseMatrix<number>::SparseMatrix (const SparsityPattern &c)
+                :
                cols(0),
                val(0),
                max_len(0)
@@ -295,7 +297,7 @@ SparseMatrix<number>::vmult (OutVector& dst,
   Assert(m() == dst.size(), ExcDimensionMismatch(m(),dst.size()));
   Assert(n() == src.size(), ExcDimensionMismatch(n(),src.size()));
 
-  Assert (&src != &dst, ExcSourceEqualsDestination());
+  Assert (!PointerComparison::equal(&src, &dst), ExcSourceEqualsDestination());
   
   const unsigned int n_rows = m();
 
@@ -328,7 +330,7 @@ SparseMatrix<number>::vmult (OutVector& dst,
         const unsigned int) const;
       const mem_fun_p comp
        = (&SparseMatrix<number>::
-           template threaded_vmult<InVector,OutVector>);
+           template threaded_vmult<OutVector,InVector>);
       Threads::ThreadGroup<> threads;
       for (unsigned int i=0; i<n_threads; ++i)
        threads += Threads::spawn (*this, comp) (dst, src,
@@ -394,7 +396,7 @@ SparseMatrix<number>::Tvmult (OutVector& dst,
   Assert(n() == dst.size(), ExcDimensionMismatch(n(),dst.size()));
   Assert(m() == src.size(), ExcDimensionMismatch(m(),src.size()));
 
-  Assert (&src != &dst, ExcSourceEqualsDestination());
+  Assert (!PointerComparison::equal(&src, &dst), ExcSourceEqualsDestination());
 
   dst = 0;
 
@@ -420,7 +422,7 @@ SparseMatrix<number>::vmult_add (OutVector& dst,
   Assert(m() == dst.size(), ExcDimensionMismatch(m(),dst.size()));
   Assert(n() == src.size(), ExcDimensionMismatch(n(),src.size()));
 
-  Assert (&src != &dst, ExcSourceEqualsDestination());
+  Assert (!PointerComparison::equal(&src, &dst), ExcSourceEqualsDestination());
 
   const unsigned int  n_rows     = m();
   const number       *val_ptr    = &val[cols->rowstart[0]];
@@ -448,7 +450,7 @@ SparseMatrix<number>::Tvmult_add (OutVector& dst,
   Assert(n() == dst.size(), ExcDimensionMismatch(n(),dst.size()));
   Assert(m() == src.size(), ExcDimensionMismatch(m(),src.size()));
 
-  Assert (&src != &dst, ExcSourceEqualsDestination());
+  Assert (!PointerComparison::equal(&src, &dst), ExcSourceEqualsDestination());
 
   for (unsigned int i=0;i<m();i++)
     for (unsigned int j=cols->rowstart[i]; j<cols->rowstart[i+1] ;j++)

In the beginning the Universe was created. This has made a lot of people very angry and has been widely regarded as a bad move.

Douglas Adams


Typeset in Trocchi and Trocchi Bold Sans Serif.