]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Make instantiation scheme a little simpler.
authorWolfgang Bangerth <bangerth@math.tamu.edu>
Mon, 12 Jul 2004 13:59:10 +0000 (13:59 +0000)
committerWolfgang Bangerth <bangerth@math.tamu.edu>
Mon, 12 Jul 2004 13:59:10 +0000 (13:59 +0000)
git-svn-id: https://svn.dealii.org/trunk@9507 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/lac/source/filtered_matrix.cc
deal.II/lac/source/filtered_matrix.in.h [deleted file]

index 82eeb008bec3c10a6ecd733cce9c7aff4e85b470..3bd0d036f6211d15b3cc5ca1dae508bd5a1d89cc 100644 (file)
@@ -2,7 +2,7 @@
 //    $Id$
 //    Version: $Name$
 //
-//    Copyright (C) 2001, 2002, 2003 by the deal.II authors
+//    Copyright (C) 2001, 2002, 2003, 2004 by the deal.II authors
 //
 //    This file is subject to QPL and may not be  distributed
 //    without copyright and license information. Please refer
 
 #include <lac/filtered_matrix.templates.h>
 
-//TODO: Maybe split in several files?
 
-#define TYPEMAT double
-#define TYPEVEC double
-#include "filtered_matrix.in.h"
+// Following, we define some partial specializations. since C++ does not allow
+// this properly, we have to work around it using some macro ugliness
 
-template class FilteredMatrix<SparseMatrix<TYPEMAT>,Vector<TYPEVEC> >;
-template class FilteredMatrix<BlockSparseMatrix<TYPEMAT>,BlockVector<TYPEVEC> >;
+#define DEFINE_SPECIALIZATIONS(TYPE) \
+template <>   \
+void   \
+FilteredMatrix<SparseMatrix<TYPE>,Vector<TYPE> >::   \
+get_column_entries (const unsigned int           index,   \
+                   std::vector<IndexValuePair> &column_entries,   \
+                   const bool                   matrix_is_symmetric) const   \
+{   \
+                                  /* depending on whether the matrix */  \
+                                  /* can be assumed symmetric or not,*/   \
+                                  /* either use a fast or a slow */  \
+                                  /* algorithm */  \
+  if (matrix_is_symmetric == true)   \
+                                    /* ok, matrix is symmetric. we */  \
+                                    /* may determine the matrix */  \
+                                    /* entries in this column by */  \
+                                    /* looking at the matrix entries */  \
+                                    /* in this row which is */  \
+                                    /* significantly faster since we */  \
+                                    /* can traverse them linearly and */  \
+                                    /* do not have to check each row */  \
+                                    /* for the possible existence of */  \
+                                    /* a matrix entry */  \
+    {   \
+      const unsigned int *   \
+       col_nums   = &(matrix->get_sparsity_pattern().get_column_numbers()   \
+                      [matrix->get_sparsity_pattern().get_rowstart_indices()[index]]);   \
+      const unsigned int   \
+       row_length = matrix->get_sparsity_pattern().row_length(index);   \
+   \
+      for (unsigned int i=0; i<row_length; ++i)   \
+       {   \
+         const unsigned int c = *(col_nums+i);   \
+   \
+                                          /* if not diagonal entry,*/   \
+                                          /* add to list */  \
+         if (c != index)   \
+           column_entries.push_back (std::make_pair(c, (*matrix)(c,index)));   \
+       };   \
+    }   \
+  else   \
+    {   \
+                                      /* otherwise check each row for */  \
+                                      /* occurrence of an entry in */  \
+                                      /* this column */  \
+      for (unsigned int row=0; row<n(); ++row)   \
+       if (row != index)   \
+         {   \
+           const unsigned int   \
+             global_index = matrix->get_sparsity_pattern()(row,index);   \
+           if (global_index != SparsityPattern::invalid_entry)   \
+             column_entries.push_back (std::make_pair(row,   \
+                                                      (*matrix)(row,index)));   \
+         };   \
+    };   \
+}   \
+   \
+   \
+   \
+template <>   \
+void   \
+FilteredMatrix<BlockSparseMatrix<TYPE>,BlockVector<TYPE> >::   \
+get_column_entries (const unsigned int           /*index*/,   \
+                   std::vector<IndexValuePair> &/*column_entries*/,   \
+                   const bool                   /*matrix_is_symmetric*/) const   \
+{   \
+                                  /* presently not implemented, but */  \
+                                  /* should be fairly simple to do */  \
+  Assert (false, ExcNotImplemented());   \
+}   \
+   \
+   \
+   \
+   \
+template <>   \
+void   \
+FilteredMatrix<SparseMatrix<TYPE>,Vector<TYPE> >::   \
+allocate_tmp_vector ()    \
+{   \
+  Threads::ThreadMutex::ScopedLock lock (tmp_mutex);   \
+  tmp_vector.reinit (matrix->n(), true);   \
+}   \
+   \
+   \
+   \
+template <>   \
+void   \
+FilteredMatrix<BlockSparseMatrix<TYPE>,BlockVector<TYPE> >::   \
+allocate_tmp_vector ()    \
+{   \
+  std::vector<unsigned int> block_sizes (matrix->n_block_rows());   \
+  for (unsigned int i=0; i<block_sizes.size(); ++i)   \
+    block_sizes[i] = matrix->block(i,i).n();   \
+     \
+  Threads::ThreadMutex::ScopedLock lock (tmp_mutex);   \
+  tmp_vector.reinit (block_sizes, true);   \
+}
 
-#undef TYPEVEC
-#undef TYPEMAT
 
+// now instantiate
+DEFINE_SPECIALIZATIONS(double)
+template class FilteredMatrix<SparseMatrix<double>,Vector<double> >;
+template class FilteredMatrix<BlockSparseMatrix<double>,BlockVector<double> >;
 
-
-#define TYPEMAT float
-#define TYPEVEC float
-#include "filtered_matrix.in.h"
-
-template class FilteredMatrix<SparseMatrix<TYPEMAT>,Vector<TYPEVEC> >;
-template class FilteredMatrix<BlockSparseMatrix<TYPEMAT>,BlockVector<TYPEVEC> >;
-
-#undef TYPEVEC
-#undef TYPEMAT
+DEFINE_SPECIALIZATIONS(float)
+template class FilteredMatrix<SparseMatrix<float>,Vector<float> >;
+template class FilteredMatrix<BlockSparseMatrix<float>,BlockVector<float> >;
 
diff --git a/deal.II/lac/source/filtered_matrix.in.h b/deal.II/lac/source/filtered_matrix.in.h
deleted file mode 100644 (file)
index 70d10e2..0000000
+++ /dev/null
@@ -1,117 +0,0 @@
-//---------------------------------------------------------------------------
-//    $Id$
-//    Version: $Name$
-//
-//    Copyright (C) 2001, 2002, 2003 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.
-//
-//---------------------------------------------------------------------------
-
-
-// File included by driver file fitered_matrix.?.cc
-// There, TYPEMAT and TYPEVEC must be set to the number types
-
-
-//TODO: Check if this cannot be implemented with iterators
-//      and applied to SparseMatrixEZ
-
-template <>
-void
-FilteredMatrix<SparseMatrix<TYPEMAT>,Vector<TYPEVEC> >::
-get_column_entries (const unsigned int           index,
-                   std::vector<IndexValuePair> &column_entries,
-                   const bool                   matrix_is_symmetric) const
-{
-                                  // depending on whether the matrix
-                                  // can be assumed symmetric or not,
-                                  // either use a fast or a slow
-                                  // algorithm
-  if (matrix_is_symmetric == true)
-                                    // ok, matrix is symmetric. we
-                                    // may determine the matrix
-                                    // entries in this column by
-                                    // looking at the matrix entries
-                                    // in this row which is
-                                    // significantly faster since we
-                                    // can traverse them linearly and
-                                    // do not have to check each row
-                                    // for the possible existence of
-                                    // a matrix entry
-    {
-      const unsigned int *
-       col_nums   = &(matrix->get_sparsity_pattern().get_column_numbers()
-                      [matrix->get_sparsity_pattern().get_rowstart_indices()[index]]);
-      const unsigned int
-       row_length = matrix->get_sparsity_pattern().row_length(index);
-
-      for (unsigned int i=0; i<row_length; ++i)
-       {
-         const unsigned int c = *(col_nums+i);
-
-                                          // if not diagonal entry,
-                                          // add to list
-         if (c != index)
-           column_entries.push_back (std::make_pair(c, (*matrix)(c,index)));
-       };
-    }
-  else
-    {
-                                      // otherwise check each row for
-                                      // occurrence of an entry in
-                                      // this column
-      for (unsigned int row=0; row<n(); ++row)
-       if (row != index)
-         {
-           const unsigned int
-             global_index = matrix->get_sparsity_pattern()(row,index);
-           if (global_index != SparsityPattern::invalid_entry)
-             column_entries.push_back (std::make_pair(row,
-                                                      (*matrix)(row,index)));
-         };
-    };
-}
-
-
-
-template <>
-void
-FilteredMatrix<BlockSparseMatrix<TYPEMAT>,BlockVector<TYPEVEC> >::
-get_column_entries (const unsigned int           /*index*/,
-                   std::vector<IndexValuePair> &/*column_entries*/,
-                   const bool                   /*matrix_is_symmetric*/) const
-{
-                                  // presently not implemented, but
-                                  // should be fairly simple to do
-  Assert (false, ExcNotImplemented());
-}
-
-
-
-
-template <>
-void
-FilteredMatrix<SparseMatrix<TYPEMAT>,Vector<TYPEVEC> >::
-allocate_tmp_vector () 
-{
-  Threads::ThreadMutex::ScopedLock lock (tmp_mutex);
-  tmp_vector.reinit (matrix->n(), true);
-}
-
-
-
-template <>
-void
-FilteredMatrix<BlockSparseMatrix<TYPEMAT>,BlockVector<TYPEVEC> >::
-allocate_tmp_vector () 
-{
-  std::vector<unsigned int> block_sizes (matrix->n_block_rows());
-  for (unsigned int i=0; i<block_sizes.size(); ++i)
-    block_sizes[i] = matrix->block(i,i).n();
-  
-  Threads::ThreadMutex::ScopedLock lock (tmp_mutex);
-  tmp_vector.reinit (block_sizes, true);
-}

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.