]> https://gitweb.dealii.org/ - dealii.git/commitdiff
make documentation available, document and remove exceptions
authorGuido Kanschat <dr.guido.kanschat@gmail.com>
Fri, 7 Aug 2009 16:51:00 +0000 (16:51 +0000)
committerGuido Kanschat <dr.guido.kanschat@gmail.com>
Fri, 7 Aug 2009 16:51:00 +0000 (16:51 +0000)
git-svn-id: https://svn.dealii.org/trunk@19197 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/base/include/base/exceptions.h
deal.II/lac/include/lac/sparse_matrix.h
deal.II/lac/include/lac/sparsity_pattern.h
deal.II/lac/source/sparsity_pattern.cc

index e03628133f9eeecf4c96cb24e280c35aca2a22a9..95b1c3a4f96776afc1a2e1d64358170dd56322e2 100644 (file)
@@ -800,6 +800,23 @@ namespace StandardExceptions
                  << " by " << arg2
                  << " has remainder different from zero");
                  
+                                  /**
+                                   * This exception is thrown if the
+                                   * iterator you access has
+                                   * corrupted data. It might for
+                                   * instance be, that the container
+                                   * it refers does not have an entry
+                                   * at the point the iterator
+                                   * refers.
+                                   *
+                                   * Typically, this will be an
+                                   * internal error of deal.II,
+                                   * because the increment and
+                                   * decrement operators should never
+                                   * yield an invalid iterator.
+                                   */
+  DeclException0 (ExcInvalidIterator);
+  
                                   /**
                                    * This exception is thrown if the
                                    * iterator you incremented or
index 3f85476cc6df6a200977c6a8c8b403fe1d0bec15..277a2137088aa935f244aec3c4c94fa2c1d16ac3 100644 (file)
@@ -61,12 +61,39 @@ namespace SparseMatrixIterators
                                    * matrix.
                                    *
                                    * The general template is not
-                                   * implemented, only the specializations
-                                   * for the two possible values of the
-                                   * second template argument.
+                                   * implemented, only the
+                                   * specializations for the two
+                                   * possible values of the second
+                                   * template argument. Therefore,
+                                   * the interface listed here only
+                                   * serves as a template provided
+                                   * since doxygen does not link the
+                                   * specializations.
                                    */
   template <typename number, bool Constness>
-  class Accessor;
+  class Accessor : public SparsityPatternIterators::Accessor
+  {
+    public:
+                                      /**
+                                       * Value of this matrix entry.
+                                       */
+      number value() const;
+   
+                                      /**
+                                       * Value of this matrix entry.
+                                       */
+      number& value();
+   
+                                      /**
+                                       * Return a reference to the matrix
+                                       * into which this accessor
+                                       * points. Note that in the present
+                                       * case, this is a constant
+                                       * reference.
+                                       */
+      const SparseMatrix<number>& get_matrix () const;
+  };
+  
   
   
                                   /**
index b4c51d34dd37ede458a1c8171c9d8fdf98f9a4e2..27e8647348c2e486793969ebae5929b4db4e216d 100644 (file)
@@ -301,10 +301,6 @@ namespace SparsityPatternIterators
                                       /** @addtogroup Exceptions
                                        * @{ */
        
-                                      /**
-                                       * Exception
-                                       */
-      DeclException0 (ExcInvalidIterator);
                                       //@}        
     protected:
                                       /**
@@ -1435,20 +1431,8 @@ class SparsityPattern : public Subscriptor
                                     /** @addtogroup Exceptions
                                      * @{ */
                                     /**
-                                     * Exception
-                                     */
-    DeclException1 (ExcInvalidNumber,
-                   int,
-                   << "The provided number is invalid here: " << arg1);
-                                    /**
-                                     * Exception
-                                     */
-    DeclException2 (ExcInvalidIndex,
-                   int, int,
-                   << "The given index " << arg1
-                   << " should be less than " << arg2 << ".");
-                                    /**
-                                     * Exception
+                                     * You tried to add an element to
+                                     * a row, but there was no space left.
                                      */
     DeclException2 (ExcNotEnoughSpace,
                    int, int,
@@ -1457,20 +1441,23 @@ class SparsityPattern : public Subscriptor
                    << "(Maximum number of entries for this row: "
                    << arg2 << "; maybe the matrix is already compressed?)");
                                     /**
-                                     * Exception
+                                     * The operation is only allowed
+                                     * after the SparsityPattern has
+                                     * been set up and compress() was
+                                     * called.
                                      */
     DeclException0 (ExcNotCompressed);
                                     /**
-                                     * Exception
+                                     * This operation changes the
+                                     * structure of the
+                                     * SparsityPattern and is not
+                                     * possible after compress() has
+                                     * been called.
                                      */
     DeclException0 (ExcMatrixIsCompressed);
                                     /**
                                      * Exception
                                      */
-    DeclException0 (ExcEmptyObject);
-                                    /**
-                                     * Exception
-                                     */
     DeclException0 (ExcInvalidConstructorCall);
                                     /**
                                      * This exception is thrown if
@@ -1496,13 +1483,6 @@ class SparsityPattern : public Subscriptor
                     int,
                     << "The number of partitions you gave is " << arg1
                     << ", but must be greater than zero.");
-                                     /**
-                                      * Exception
-                                      */
-    DeclException2 (ExcInvalidArraySize,
-                    int, int,
-                    << "The array has size " << arg1 << " but should have size "
-                    << arg2);
                                     //@}
   private:
                                     /**
@@ -1687,6 +1667,17 @@ namespace SparsityPatternIterators
   {}
 
 
+  inline
+  bool
+  Accessor::is_valid_entry () const
+  {
+    return (sparsity_pattern
+           ->get_column_numbers()[sparsity_pattern
+                                  ->get_rowstart_indices()[a_row]+a_index]
+           != SparsityPattern::invalid_entry);
+  }
+
+
   inline
   unsigned int
   Accessor::row() const
@@ -1720,17 +1711,6 @@ namespace SparsityPatternIterators
 
 
 
-  inline
-  bool
-  Accessor::is_valid_entry () const
-  {
-    return (sparsity_pattern
-           ->get_column_numbers()[sparsity_pattern
-                                  ->get_rowstart_indices()[a_row]+a_index]
-           != SparsityPattern::invalid_entry);
-  }
-
-
 
   inline
   bool
@@ -2094,7 +2074,7 @@ SparsityPattern::copy_from (const unsigned int    n_rows,
        {
          const unsigned int col
            = internal::SparsityPatternTools::get_column_index_from_iterator(*j);
-         Assert (col < n_cols, ExcInvalidIndex(col,n_cols));
+         Assert (col < n_cols, ExcIndexRange(col,0,n_cols));
          
          if ((col!=row) || !is_square)
            *cols++ = col;
index 4ea77db4a2e82bd328753b68a4e9c91bbced877e..09cdff51120ba33bda2894c71032169599fa38dd 100644 (file)
@@ -2,7 +2,7 @@
 //    $Id$
 //    Version: $Name$
 //
-//    Copyright (C) 2000, 2001, 2002, 2003, 2004, 2005, 2006, 2007, 2008 by the deal.II authors
+//    Copyright (C) 2000, 2001, 2002, 2003, 2004, 2005, 2006, 2007, 2008, 2009 by the deal.II authors
 //
 //    This file is subject to QPL and may not be  distributed
 //    without copyright and license information. Please refer
@@ -279,7 +279,7 @@ SparsityPattern::reinit (
   const VectorSlice<const std::vector<unsigned int> >&row_lengths,
   const bool optimize_diag)
 {
-  Assert (row_lengths.size() == m, ExcInvalidNumber (m));
+  AssertDimension (row_lengths.size(), m);
          
   rows = m;
   cols = n;
@@ -589,7 +589,7 @@ SparsityPattern::copy_from (const CompressedSparsityPattern &csp,
       for (unsigned int j=0; j<row_length; ++j)
         {
           const unsigned int col = csp.column_number(row,j);
-         Assert (col < csp.n_cols(), ExcInvalidIndex(col,csp.n_cols()));
+         Assert (col < csp.n_cols(), ExcIndexRange(col,0,csp.n_cols()));
          
          if ((col!=row) || !is_square)
            *cols++ = col;
@@ -645,8 +645,6 @@ SparsityPattern::copy_from (const CompressedSetSparsityPattern &csp,
       for (; col_num != csp.row_end (row); ++col_num)
        {
          const unsigned int col = *col_num;
-         //Assert (col < csp.n_cols(), ExcInvalidIndex(col,csp.n_cols()));
-         
          if ((col!=row) || !is_square)
            *cols++ = col;
        }
@@ -698,7 +696,7 @@ SparsityPattern::copy_from (const CompressedSimpleSparsityPattern &csp,
       for (unsigned int j=0; j<row_length; ++j)
        {
          const unsigned int col = csp.column_number(row,j);
-          Assert (col < csp.n_cols(), ExcInvalidIndex(col,csp.n_cols()));
+          Assert (col < csp.n_cols(), ExcIndexRange(col,0,csp.n_cols()));
 
           if ((col!=row) || !is_square)
             *cols++ = col;
@@ -816,8 +814,8 @@ SparsityPattern::operator () (const unsigned int i,
                              const unsigned int j) const
 {
   Assert ((rowstart!=0) && (colnums!=0), ExcEmptyObject());  
-  Assert (i<rows, ExcInvalidIndex(i,rows));
-  Assert (j<cols, ExcInvalidIndex(j,cols));
+  Assert (i<rows, ExcIndexRange(i,0,rows));
+  Assert (j<cols, ExcIndexRange(j,0,cols));
   Assert (compressed, ExcNotCompressed());
 
                                   // let's see whether there is
@@ -864,8 +862,8 @@ SparsityPattern::add (const unsigned int i,
                      const unsigned int j)
 {
   Assert ((rowstart!=0) && (colnums!=0), ExcEmptyObject());  
-  Assert (i<rows, ExcInvalidIndex(i,rows));
-  Assert (j<cols, ExcInvalidIndex(j,cols));
+  Assert (i<rows, ExcIndexRange(i,0,rows));
+  Assert (j<cols, ExcIndexRange(j,0,cols));
   Assert (compressed==false, ExcMatrixIsCompressed());
 
   for (unsigned int k=rowstart[i]; k<rowstart[i+1]; k++)

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.